Loading...
Searching...
No Matches
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 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
30#include "HashSet.H"
32#include "adjointSensitivity.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
40 defineTypeNameAndDebug(ShapeSensitivitiesBase, 0);
41}
42
43
44// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
45
47{
48 // Allocate distance solver if needed
50 {
51 eikonalSolver_.reset
52 (
54 (
55 mesh_,
56 this->dict(),
59 )
60 );
61 }
62}
63
64
66(
67 bool (objective::*hasFunction)() const
68)
69{
70 bool hasMult(false);
71 const PtrList<objective>& objectives =
72 adjointSolver_.getObjectiveManager().getObjectiveFunctions();
73 for (const objective& func : objectives)
74 {
75 hasMult = hasMult || (func.*hasFunction)();
76 }
77 return hasMult;
78}
79
80
82{
83 // Wall face sensitivity projected to normal
84 if (wallFaceSensNormalPtr_)
85 {
86 constructAndWriteSensitivityField<scalar>
87 (
88 wallFaceSensNormalPtr_,
89 "faceSensNormal" + suffix_
90 );
91 }
92
93 if (writeAllSurfaceFiles_)
94 {
95 // Wall face sensitivity vectors
96 if (wallFaceSensVecPtr_)
97 {
98 constructAndWriteSensitivityField<vector>
99 (
100 wallFaceSensVecPtr_,
101 "faceSensVec" + suffix_
102 );
103 }
104
105 // Normal sens as vectors
106 if (wallFaceSensNormalVecPtr_)
107 {
108 constructAndWriteSensitivityField<vector>
109 (
110 wallFaceSensNormalVecPtr_,
111 "faceSensNormalVec" + suffix_
112 );
113 }
114 }
115}
116
117
119{
120 // Wall point sensitivity projected to normal
121 if (wallPointSensNormalPtr_)
122 {
123 constructAndWriteSensitivtyPointField<scalar>
124 (
125 wallPointSensNormalPtr_,
126 "pointSensNormal" + suffix_
127 );
128 }
129
130 // Write point-based sensitivities, if present
131 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
132
133 if (writeAllSurfaceFiles_)
134 {
135 // Wall point sensitivity vectors
136 if (wallPointSensVecPtr_)
137 {
138 constructAndWriteSensitivtyPointField<vector>
139 (
140 wallPointSensVecPtr_,
141 "pointSensVec" + suffix_
142 );
143 }
144
145 // Normal point as vectors
146 if (wallPointSensNormalVecPtr_)
147 {
148 constructAndWriteSensitivtyPointField<vector>
149 (
150 wallPointSensNormalVecPtr_,
151 "pointSensNormalVec" + suffix_
152 );
153 }
154 }
155}
156
157
159{
160 // Face-based boundary sens
161 if (wallFaceSensVecPtr_)
162 {
163 wallFaceSensVecPtr_() = vector::zero;
164 }
165 if (wallFaceSensNormalVecPtr_)
166 {
167 wallFaceSensNormalVecPtr_() = vector::zero;
168 }
169 if (wallFaceSensNormalPtr_)
170 {
171 wallFaceSensNormalPtr_() = scalar(0);
172 }
173
174 // Point-based boundary sens
175 if (wallPointSensVecPtr_)
176 {
177 for (vectorField& patchSens : wallPointSensVecPtr_())
178 {
179 patchSens = vector::zero;
180 }
181 }
182 if (wallPointSensNormalVecPtr_)
183 {
184 for (vectorField& patchSens : wallPointSensNormalVecPtr_())
185 {
186 patchSens = vector::zero;
187 }
188 }
189 if (wallPointSensNormalPtr_)
190 {
191 for (scalarField& patchSens : wallPointSensNormalPtr_())
193 patchSens = scalar(0);
194 }
195 }
196}
197
198
200{
201 gradDxDbMult_.reset
202 (
204 (
206 (
207 "gradDxDbMult",
208 mesh_.time().timeName(),
209 mesh_,
212 ),
213 mesh_,
215 )
216 );
217 if (hasMultiplier(&objective::hasDivDxDbMult))
218 {
219 divDxDbMult_.reset(new scalarField(mesh_.nCells(), Zero));
220 }
221 if (hasMultiplier(&objective::hasdSdbMult))
222 {
223 dSfdbMult_.reset(createZeroBoundaryPtr<vector>(mesh_));
224 }
225 if (hasMultiplier(&objective::hasdndbMult))
226 {
227 dnfdbMult_.reset(createZeroBoundaryPtr<vector>(mesh_));
228 }
229 if (hasMultiplier(&objective::hasdxdbDirectMult))
230 {
232 }
233 bcDxDbMult_.reset(createZeroBoundaryPtr<vector>(mesh_));
234 optionsDxDbMult_.reset(new vectorField(mesh_.nCells(), Zero));
235}
236
237
239{
240 gradDxDbMult_() = Zero;
241 if (divDxDbMult_)
242 {
243 divDxDbMult_() = Zero;
244 }
245 if (eikonalSolver_)
246 {
247 eikonalSolver_->reset();
248 }
249 if (dxdbMult_)
250 {
251 dxdbMult_() = Zero;
252 }
253 if (dSfdbMult_)
254 {
255 dSfdbMult_() = Zero;
256 }
257 if (dnfdbMult_)
258 {
259 dnfdbMult_() = Zero;
260 }
261 if (dxdbDirectMult_)
262 {
263 dxdbDirectMult_() = Zero;
264 }
265 if (pointDxDbDirectMult_)
266 {
267 for (vectorField& field : pointDxDbDirectMult_())
268 {
269 field = Zero;
270 }
271 }
274}
275
276
277// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
278
279Foam::ShapeSensitivitiesBase::ShapeSensitivitiesBase
280(
281 const fvMesh& mesh,
282 const dictionary& dict,
284)
285:
287 sensitivityPatchIDs_
288 (
289 mesh.boundaryMesh().patchSet
290 (
291 dict.optionalSubDict(mesh.name()).
292 get<wordRes>("patches", keyType::REGEX_RECURSIVE)
293 )
294 ),
295 writeAllSurfaceFiles_
296 (
297 dict.getOrDefault<bool>("writeAllSurfaceFiles", false)
298 ),
299 wallFaceSensVecPtr_(nullptr),
300 wallFaceSensNormalPtr_(nullptr),
301 wallFaceSensNormalVecPtr_(nullptr),
302
303 wallPointSensVecPtr_(nullptr),
304 wallPointSensNormalPtr_(nullptr),
305 wallPointSensNormalVecPtr_(nullptr)
306{
309}
310
311
312// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
313
315{
317 {
318 sensitivityPatchIDs_ =
319 mesh_.boundaryMesh().patchSet
320 (
321 dict_.optionalSubDict(mesh_.name()).
322 get<wordRes>("patches", keyType::REGEX_RECURSIVE)
323 );
324 writeAllSurfaceFiles_ =
325 dict_.getOrDefault<bool>("writeAllSurfaceFiles", false);
326
327 if (includeDistance_)
328 {
329 if (eikonalSolver_)
330 {
331 eikonalSolver_().readDict(dict);
332 }
333 else
334 {
335 allocateEikonalSolver();
336 }
337 }
338
339 return true;
341
342 return false;
343}
344
345
348{
350}
351
352
354{
355 // Accumulate multiplier of grad(dxdb)
356 adjointSolver_.accumulateGradDxDbMultiplier(gradDxDbMult_(), dt);
357
358 // Accumulate multiplier of div(dxdb)
359 adjointSolver_.accumulateDivDxDbMultiplier(divDxDbMult_, dt);
360
361 // Terms from fvOptions - missing contributions from turbulence models
362 adjointSolver_.accumulateOptionsDxDbMultiplier(optionsDxDbMult_(), dt);
363
364 // Accumulate source for the adjoint to the eikonal equation
365 if (eikonalSolver_)
366 {
367 eikonalSolver_->accumulateIntegrand(dt);
368 }
369
370 // Accumulate direct sensitivities
371 adjointSolver_.accumulateGeometryVariationsMultipliers
372 (
373 dSfdbMult_,
374 dnfdbMult_,
375 dxdbDirectMult_,
376 pointDxDbDirectMult_,
377 geometryVariationIntegrationPatches(),
378 dt
379 );
381 // Accumulate sensitivities due to boundary conditions
382 adjointSolver_.accumulateBCSensitivityIntegrand
384}
385
386
394
395void Foam::ShapeSensitivitiesBase::write(const word& baseName)
396{
400}
401
402
405{
406 if (wallFaceSensVecPtr_)
407 {
408 return
409 constructVolSensitivtyField<vector>
410 (
411 wallFaceSensVecPtr_,
412 "faceSensVec" + suffix_
413 );
414 }
415 else
416 {
418 << " no faceSensVec boundary field. Returning zero" << endl;
419
420 return
421 tmp<volVectorField>
422 (
424 (
425 mesh_,
426 "faceSensVec" + suffix_,
427 dimless
428 ).ptr()
429 );
430 }
431}
432
433
436{
437 if (wallFaceSensNormalPtr_)
438 {
439 return
440 constructVolSensitivtyField<scalar>
441 (
442 wallFaceSensNormalPtr_,
443 "faceSensNormal" + suffix_
444 );
445 }
446 else
447 {
449 << " no wallFaceSensNormal boundary field. Returning zero" << endl;
450
451 return
453 (
455 (
456 mesh_,
457 "faceSensNormal" + suffix_, dimless
458 ).ptr()
459 );
460 }
461}
462
463
466{
467 if (wallFaceSensNormalVecPtr_)
468 {
469 return
470 constructVolSensitivtyField<vector>
471 (
472 wallFaceSensNormalVecPtr_,
473 "faceSensNormalVec" + suffix_
474 );
475 }
476 else
477 {
479 << " no wallFaceSensNormalVec boundary field. Returning zero"
480 << endl;
481
482 return
484 (
486 (
487 mesh_,
488 "faceSensNormalVec" + suffix_,
489 dimless
490 ).ptr()
491 );
492 }
493}
494
495
498{
499 tmp<volVectorField> tWallFaceSensVec = getWallFaceSensVec();
501
502 return (volPointInter.interpolate(tWallFaceSensVec));
503}
504
505
508{
509 tmp<volScalarField> tWallFaceSensNormal = getWallFaceSensNormal();
511
512 return (volPointInter.interpolate(tWallFaceSensNormal));
513}
514
515
518{
519 tmp<volVectorField> tWallFaceSensNormalVec = getWallFaceSensNormalVec();
521
522 return (volPointInter.interpolate(tWallFaceSensNormalVec));
523}
524
525
531
532
538
539
542{
543 return wallFaceSensNormalVecPtr_();
544}
545
546
547// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Base class supporting Shape sensitivity derivatives.
tmp< volVectorField > getWallFaceSensVec()
Get wall face sensitivity vectors field.
autoPtr< pointBoundaryVectorField > wallPointSensNormalVecPtr_
Normal sens as vectors.
void clearSensitivities()
Zero sensitivity fields and their constituents.
bool hasMultiplier(bool(objective::*hasFunction)() const)
Check if any of the available objective has a certain multiplier, provided through a function object.
bool writeAllSurfaceFiles_
Whether to write all surface sensitivity fields.
tmp< volVectorField > getWallFaceSensNormalVec()
Get wall face normal sens as vectors field.
autoPtr< boundaryVectorField > wallFaceSensVecPtr_
Wall face sens w.r.t. (x,y.z).
labelHashSet sensitivityPatchIDs_
Patches on which to compute shape sensitivities.
virtual const boundaryScalarField & getWallFaceSensNormalBoundary() const
Get wall face sensitivity projected to normal field.
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.
autoPtr< pointBoundaryVectorField > wallPointSensVecPtr_
Wall point sens w.r.t. (x,y.z).
tmp< pointVectorField > getWallPointSensVec()
Get wall point sensitivity vectors field.
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
tmp< pointScalarField > getWallPointSensNormal()
Get wall point sensitivity projected to normal field.
void allocateEikonalSolver()
Allocate the adjoint eikonal solver.
void clearSurfaceFields()
Clear surface/point fields.
tmp< GeometricField< Type, fvPatchField, volMesh > > constructVolSensitivtyField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
void allocateMultipliers()
Allocate multiplier fields.
tmp< pointVectorField > getWallPointSensNormalVec()
Get wall point sens as vectors field.
autoPtr< pointBoundaryScalarField > wallPointSensNormalPtr_
Wall point sens projected to normal.
tmp< volScalarField > getWallFaceSensNormal()
Get wall face sensitivity projected to normal field.
void constructAndWriteSensitivityField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
autoPtr< boundaryScalarField > wallFaceSensNormalPtr_
Wall face sens projected to normal.
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
void constructAndWriteSensitivtyPointField(const autoPtr< List< Field< Type > > > &sensFieldPtr, const word &name) const
Constructs pointField based on boundaryField and writes it.
autoPtr< boundaryVectorField > wallFaceSensNormalVecPtr_
Normal sens as vectors.
virtual const labelHashSet & geometryVariationIntegrationPatches() const
Return set of patches on which to compute direct sensitivities.
void clearMultipliers()
Clear multipliers.
static const Form zero
Solver of the adjoint to the eikonal PDE.
Abstract base class for adjoint-based sensitivities.
autoPtr< scalarField > divDxDbMult_
Multiplier of div(dx/db).
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
bool includeDistance_
Include distance variation in sensitivity computations.
autoPtr< boundaryVectorField > dxdbDirectMult_
Multiplier of dCf/db, found in the objective function.
autoPtr< pointBoundaryVectorField > pointDxDbDirectMult_
Multiplier of dx/db computed at points, found in the objective function.
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
autoPtr< boundaryVectorField > bcDxDbMult_
Multiplier of dx/db, coming from boundary conditions that depend on the geometry, like rotatingWallVe...
adjointSolver & adjointSolver_
Reference to the underlaying adjoint solver.
autoPtr< vectorField > optionsDxDbMult_
dx/db multiplier coming from fvOptions
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
autoPtr< boundaryVectorField > dxdbMult_
Multiplier of face dx/db.
autoPtr< boundaryVectorField > dnfdbMult_
Multiplier of dnf/db.
autoPtr< boundaryVectorField > dSfdbMult_
Multiplier of dSf/db.
autoPtr< volTensorField > gradDxDbMult_
Multiplier of grad(dx/b).
word suffix_
Append this word to files related to the sensitivities.
Base class for adjoint solvers.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for handling keywords in dictionaries.
Definition keyType.H:69
@ REGEX_RECURSIVE
Definition keyType.H:88
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
bool hasdSdbMult() const noexcept
Definition objectiveI.H:181
bool hasdndbMult() const noexcept
Definition objectiveI.H:187
bool hasDivDxDbMult() const noexcept
Definition objectiveI.H:211
bool hasdxdbDirectMult() const noexcept
Definition objectiveI.H:199
const fvMesh & mesh_
Definition sensitivity.H:68
dictionary dict_
Definition sensitivity.H:69
const fvMesh & mesh() const
Return reference to mesh.
const dictionary & dict() const
Return the construction dictionary.
A class for managing temporary objects.
Definition tmp.H:75
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:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
rDeltaTY field()
dynamicFvMesh & mesh
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
volScalarField::Boundary boundaryScalarField
volFields
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
autoPtr< GeometricField< Type, fvPatchField, volMesh > > createZeroFieldPtr(const fvMesh &mesh, const word &name, const dimensionSet dims, bool printAllocation=false)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > createZeroBoundaryPtr(const fvMesh &mesh, bool printAllocation=false)
volVectorField::Boundary boundaryVectorField
dimensioned< tensor > dimensionedTensor
Dimensioned tensor obtained from generic dimensioned type.
dictionary dict