shapeSensitivitiesBase.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2007-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 FOSS GP
10 Copyright (C) 2019-2020 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
38}
39
40
41// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42
44{
45 // Wall face sensitivity projected to normal
47 {
48 constructAndWriteSensitivityField<scalar>
49 (
51 "faceSensNormal" + surfaceFieldSuffix_
52 );
53 }
54
56 {
57 // Wall face sensitivity vectors
59 {
60 constructAndWriteSensitivityField<vector>
61 (
63 "faceSensVec" + surfaceFieldSuffix_
64 );
65 }
66
67 // Normal sens as vectors
69 {
70 constructAndWriteSensitivityField<vector>
71 (
73 "faceSensNormalVec" + surfaceFieldSuffix_
74 );
75 }
76 }
77}
78
79
81{
82 // Wall point sensitivity projected to normal
83 if (wallPointSensNormalPtr_)
84 {
85 constructAndWriteSensitivtyPointField<scalar>
86 (
87 wallPointSensNormalPtr_,
88 "pointSensNormal" + surfaceFieldSuffix_
89 );
90 }
91
92 // Write point-based sensitivities, if present
93 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
94
95 if (writeAllSurfaceFiles_)
96 {
97 // Wall point sensitivity vectors
98 if (wallPointSensVecPtr_)
99 {
100 constructAndWriteSensitivtyPointField<vector>
101 (
102 wallPointSensVecPtr_,
103 "pointSensVec" + surfaceFieldSuffix_
104 );
105 }
106
107 // Normal point as vectors
108 if (wallPointSensNormalVecPtr_)
109 {
110 constructAndWriteSensitivtyPointField<vector>
111 (
112 wallPointSensNormalVecPtr_,
113 "pointSensNormalVec" + surfaceFieldSuffix_
114 );
115 }
116 }
117}
118
119
120
121
122// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123
125(
126 const fvMesh& mesh,
127 const dictionary& dict
128)
129:
130 meshShape_(mesh),
131 surfaceFieldSuffix_(word::null),
132 writeAllSurfaceFiles_
133 (
134 dict.getOrDefault<bool>
135 (
136 "writeAllSurfaceFiles",
137 false
138 )
139 ),
140 sensitivityPatchIDs_
141 (
142 mesh.boundaryMesh().patchSet
143 (
144 dict.get<wordRes>("patches", keyType::REGEX_RECURSIVE)
145 )
146 ),
147 wallFaceSensVecPtr_(nullptr),
148 wallFaceSensNormalPtr_(nullptr),
149 wallFaceSensNormalVecPtr_(nullptr),
150
151 wallPointSensVecPtr_(nullptr),
152 wallPointSensNormalPtr_(nullptr),
153 wallPointSensNormalVecPtr_(nullptr)
154{}
155
156
157// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158
161{
162 return sensitivityPatchIDs_;
163}
164
165
167(
168 const labelHashSet& sensPatchIDs
169)
170{
171 sensitivityPatchIDs_ = sensPatchIDs;
172}
173
174
176{
177 // Face-based boundary sens
178 if (wallFaceSensVecPtr_)
179 {
180 wallFaceSensVecPtr_() = vector::zero;
181 }
182 if (wallFaceSensNormalVecPtr_)
183 {
184 wallFaceSensNormalVecPtr_() = vector::zero;
185 }
186 if (wallFaceSensNormalPtr_)
187 {
188 wallFaceSensNormalPtr_() = scalar(0);
189 }
190
191 // Point-based boundary sens
192 if (wallPointSensVecPtr_)
193 {
194 for (vectorField& patchSens : wallPointSensVecPtr_())
195 {
196 patchSens = vector::zero;
197 }
198 }
199 if (wallPointSensNormalVecPtr_)
200 {
201 for (vectorField& patchSens : wallPointSensNormalVecPtr_())
202 {
203 patchSens = vector::zero;
204 }
205 }
206 if (wallPointSensNormalPtr_)
207 {
208 for (scalarField& patchSens : wallPointSensNormalPtr_())
209 {
210 patchSens = scalar(0);
211 }
212 }
213}
214
215
217{
218 writeFaceBasedSens();
219 writePointBasedSens();
220}
221
222
224{
225 surfaceFieldSuffix_ = suffix;
226}
227
228
231{
232 if (wallFaceSensVecPtr_)
233 {
234 return
235 constructVolSensitivtyField<vector>
236 (
237 wallFaceSensVecPtr_,
238 "faceSensVec" + surfaceFieldSuffix_
239 );
240 }
241 else
242 {
244 << " no faceSensVec boundary field. Returning zero" << endl;
245
246 return
248 (
249 createZeroFieldPtr<vector>
250 (
251 meshShape_,
252 "faceSensVec" + surfaceFieldSuffix_,
253 dimless
254 ).ptr()
255 );
256 }
257}
258
259
262{
263 if (wallFaceSensNormalPtr_)
264 {
265 return
266 constructVolSensitivtyField<scalar>
267 (
268 wallFaceSensNormalPtr_,
269 "faceSensNormal" + surfaceFieldSuffix_
270 );
271 }
272 else
273 {
275 << " no wallFaceSensNormal boundary field. Returning zero" << endl;
276
277 return
279 (
280 createZeroFieldPtr<scalar>
281 (
282 meshShape_,
283 "faceSensNormal" + surfaceFieldSuffix_, dimless
284 ).ptr()
285 );
286 }
287}
288
289
292{
293 if (wallFaceSensNormalVecPtr_)
294 {
295 return
296 constructVolSensitivtyField<vector>
297 (
298 wallFaceSensNormalVecPtr_,
299 "faceSensNormalVec" + surfaceFieldSuffix_
300 );
301 }
302 else
303 {
305 << " no wallFaceSensNormalVec boundary field. Returning zero"
306 << endl;
307
308 return
310 (
311 createZeroFieldPtr<vector>
312 (
313 meshShape_,
314 "faceSensNormalVec" + surfaceFieldSuffix_,
315 dimless
316 ).ptr()
317 );
318 }
319}
320
321
324{
325 tmp<volVectorField> tWallFaceSensVec = getWallFaceSensVec();
326 volPointInterpolation volPointInter(meshShape_);
327
328 return (volPointInter.interpolate(tWallFaceSensVec));
329}
330
331
334{
335 tmp<volScalarField> tWallFaceSensNormal = getWallFaceSensNormal();
336 volPointInterpolation volPointInter(meshShape_);
337
338 return (volPointInter.interpolate(tWallFaceSensNormal));
339}
340
341
344{
345 tmp<volVectorField> tWallFaceSensNormalVec = getWallFaceSensNormalVec();
346 volPointInterpolation volPointInter(meshShape_);
347
348 return (volPointInter.interpolate(tWallFaceSensNormalVec));
349}
350
351
354{
355 return wallFaceSensVecPtr_();
356}
357
358
361{
362 return wallFaceSensNormalPtr_();
363}
364
365
368{
369 return wallFaceSensNormalVecPtr_();
370}
371
372
373// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for handling keywords in dictionaries.
Definition: keyType.H:71
tmp< volVectorField > getWallFaceSensVec()
Get wall face sensitivity vectors field.
void clearSensitivities()
Zero sensitivity fields and their constituents.
tmp< volVectorField > getWallFaceSensNormalVec()
Get wall face normal sens as vectors field.
autoPtr< boundaryVectorField > wallFaceSensVecPtr_
Wall face sens w.r.t. (x,y.z)
virtual const boundaryScalarField & getWallFaceSensNormalBoundary() const
Get wall face sensitivity projected to normal field.
void setSuffix(const word &suffix)
Set suffix.
virtual const boundaryVectorField & getWallFaceSensNormalVecBoundary() const
Get wall face normal sens as vectors field.
void writeFaceBasedSens() const
Write face-based sensitivities, if present.
void writePointBasedSens() const
Write point-based sensitivities, if present.
virtual const boundaryVectorField & getWallFaceSensVecBoundary() const
Get wall face sensitivity vectors field.
tmp< pointVectorField > getWallPointSensVec()
Get wall point sensitivity vectors field.
tmp< pointScalarField > getWallPointSensNormal()
Get wall point sensitivity projected to normal field.
tmp< pointVectorField > getWallPointSensNormalVec()
Get wall point sens as vectors field.
tmp< volScalarField > getWallFaceSensNormal()
Get wall face sensitivity projected to normal field.
void write()
Write sensitivity fields.
autoPtr< boundaryScalarField > wallFaceSensNormalPtr_
Wall face sens projected to normal.
void setSensitivityPatchIDs(const labelHashSet &sensPatchIDs)
Overwrite sensitivityPatchIDs.
const labelHashSet & sensitivityPatchIDs() const
Get patch IDs on which sensitivities are computed.
autoPtr< boundaryVectorField > wallFaceSensNormalVecPtr_
Normal sens as vectors.
A class for managing temporary objects.
Definition: tmp.H:65
Interpolate from cell centres to points (vertices) using inverse distance weighting.
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate volField using inverse distance weighting.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
bool
Definition: EEqn.H:20
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
const dimensionSet dimless
Dimensionless.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dictionary dict
A non-counting (dummy) refCount.
Definition: refCount.H:59