Loading...
Searching...
No Matches
ShapeSensitivitiesBase.H
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 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
28Class
29 Foam::ShapeSensitivitiesBase
30
31Description
32 Base class supporting Shape sensitivity derivatives.
33
34 Reference:
35 \verbatim
36 For the FI formulation see
37 Kavvadias, I., Papoutsis-Kiachagias, E., & Giannakoglou, K. (2015).
38 On the proper treatment of grid sensitivities in continuous adjoint
39 methods for shape optimization.
40 Journal of Computational Physics, 301, 1–18.
41 http://doi.org/10.1016/j.jcp.2015.08.012
42
43 The ESI formulation is derived in a slightly different way than the
44 one described in this paper, to provide a common mathematical
45 formulation for both low- and high-Re meshes and to produce numerically
46 identical results as the FI formulation. In brief, the boundary-bound
47 part of the sensitivities is the patchInternalField of the tensor
48 multiplying grad(dxdb) in the FI formulation.
49
50 \endverbatim
51
52SourceFiles
53 ShapeSensitivitiesBase.C
54
55\*---------------------------------------------------------------------------*/
56
57#ifndef ShapeSensitivitiesBase_H
58#define ShapeSensitivitiesBase_H
59
60#include "adjointSensitivity.H"
61#include "objective.H"
63#include "pointMesh.H"
64#include "pointPatchField.H"
65#include "pointPatchFieldsFwd.H"
67#include "boundaryFieldsFwd.H"
68#include "createZeroField.H"
69
70
71// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72
73namespace Foam
74{
76/*---------------------------------------------------------------------------*\
77 Class ShapeSensitivitiesBase Declaration
78\*---------------------------------------------------------------------------*/
79
80class ShapeSensitivitiesBase
81 :
82 public adjointSensitivity
83{
84
85protected:
86
87 // Protected data
88
89 //- Patches on which to compute shape sensitivities
91
92 //- Whether to write all surface sensitivity fields
94
95 // autoPtrs for fields holding sensitivities.
96 // Not all of them are required for each case
97
98 // Boundary sensitivities at faces. Shape opt & flow control
99 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
100
101 //- Wall face sens w.r.t. (x,y.z)
104 //- Wall face sens projected to normal
106
107 //- Normal sens as vectors
109
110 // Boundary sensitivities at points. Shape opt
111 //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
112
113 //- Wall point sens w.r.t. (x,y.z)
115
116 //- Wall point sens projected to normal
118
119 //- Normal sens as vectors
122
123 // Protected Member Functions
124
125 //- Allocate the adjoint eikonal solver
127
128 //- Check if any of the available objective has a certain multiplier,
129 //- provided through a function object
130 bool hasMultiplier(bool (objective::*hasFunction)() const);
132 //- Constructs volField based on boundaryField and writes it
133 template<class Type>
135 (
136 const autoPtr
137 <
139 >& sensFieldPtr,
140 const word& name
141 ) const;
142
143 //- Constructs pointField based on boundaryField and writes it
144 template<class Type>
146 (
147 const autoPtr<List<Field<Type>>>& sensFieldPtr,
148 const word& name
149 ) const;
150
151 //- Constructs volField based on boundaryField and writes it
152 template<class Type>
155 (
156 const autoPtr
157 <
159 >& sensFieldPtr,
160 const word& name
161 ) const;
162
163 //- Write face-based sensitivities, if present
164 void writeFaceBasedSens() const;
165
166 //- Write point-based sensitivities, if present
167 void writePointBasedSens() const;
168
169 //- Clear surface/point fields
170 void clearSurfaceFields();
171
172 //- Allocate multiplier fields
173 void allocateMultipliers();
174
175 //- Clear multipliers
176 void clearMultipliers();
177
178
179private:
180
181 // Private Member Functions
182
183 //- No copy construct
184 ShapeSensitivitiesBase(const ShapeSensitivitiesBase&) = delete;
185
186 //- No copy assignment
187 void operator=(const ShapeSensitivitiesBase&) = delete;
188
189
190public:
191
192 //- Runtime type information
193 TypeName("ShapeSensitivitiesBase");
194
195
196 // Constructors
197
198 //- Construct from components
199 ShapeSensitivitiesBase
200 (
201 const fvMesh& mesh,
202 const dictionary& dict,
204 );
205
206
207 //- Destructor
208 virtual ~ShapeSensitivitiesBase() = default;
209
210
211 // Member Functions
212
213 //- Read dict if changed
214 virtual bool readDict(const dictionary& dict);
215
216 //- Get patch IDs on which sensitivities are computed
217 inline const labelHashSet& sensitivityPatchIDs() const
218 {
220 }
221
222 //- Overwrite sensitivityPatchIDs
223 inline void setSensitivityPatchIDs(const labelHashSet& sensPatchIDs)
224 {
225 sensitivityPatchIDs_ = sensPatchIDs;
226 }
227
228 //- Return set of patches on which to compute direct sensitivities
231 //- Accumulate sensitivity integrands
232 // Common function for the FI and E-SI approaches
233 virtual void accumulateIntegrand(const scalar dt);
234
235 //- Zero sensitivity fields and their constituents
236 void clearSensitivities();
237
238 //- Write sensitivity fields.
239 // If valid, copies boundaryFields to volFields and writes them.
240 virtual void write(const word& baseName = word::null);
241
242 //- Get wall face sensitivity vectors field
244
245 //- Get wall face sensitivity projected to normal field
247
248 //- Get wall face normal sens as vectors field
250
251 //- Get wall point sensitivity vectors field
252 // Uses volPointInterpolation
254
255 //- Get wall point sensitivity projected to normal field
256 // Uses volPointInterpolation
258
259 //- Get wall point sens as vectors field
260 // Uses volPointInterpolation
263 //- Get wall face sensitivity vectors field
264 virtual const boundaryVectorField& getWallFaceSensVecBoundary() const;
265
266 //- Get wall face sensitivity projected to normal field
267 virtual const boundaryScalarField&
269
270 //- Get wall face normal sens as vectors field
271 virtual const boundaryVectorField&
273};
274
275
276// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277
278} // End namespace Foam
279
280// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281
282#ifdef NoRepository
284#endif
285
286// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287
288#endif
289
290// ************************************************************************* //
Useful typenames for fields defined only at the boundaries.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of boundary fields.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
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.
TypeName("ShapeSensitivitiesBase")
Runtime type information.
autoPtr< pointBoundaryVectorField > wallPointSensVecPtr_
Wall point sens w.r.t. (x,y.z).
tmp< pointVectorField > getWallPointSensVec()
Get wall point sensitivity vectors field.
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.
virtual ~ShapeSensitivitiesBase()=default
Destructor.
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.
const labelHashSet & sensitivityPatchIDs() const
Get patch IDs on which sensitivities are computed.
autoPtr< boundaryScalarField > wallFaceSensNormalPtr_
Wall face sens projected to normal.
void setSensitivityPatchIDs(const labelHashSet &sensPatchIDs)
Overwrite sensitivityPatchIDs.
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.
Base class for adjoint solvers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition objective.H:58
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
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
Namespace for OpenFOAM.
volScalarField::Boundary boundaryScalarField
volFields
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
volVectorField::Boundary boundaryVectorField
runTime write()
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68