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 -------------------------------------------------------------------------------
12 License
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 Class
29  Foam::ShapeSensitivitiesBase
30 
31 Description
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 
52 SourceFiles
53  ShapeSensitivitiesBase.C
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef ShapeSensitivitiesBase_H
58 #define ShapeSensitivitiesBase_H
59 
60 #include "adjointSensitivity.H"
61 #include "objective.H"
62 #include "volPointInterpolation.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 
73 namespace Foam
74 {
75 
76 /*---------------------------------------------------------------------------*\
77  Class ShapeSensitivitiesBase Declaration
78 \*---------------------------------------------------------------------------*/
79 
81  :
82  public adjointSensitivity
83 {
84 
85 protected:
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 
179 private:
180 
181  // Private Member Functions
182 
183  //- No copy construct
185 
186  //- No copy assignment
187  void operator=(const ShapeSensitivitiesBase&) = delete;
188 
189 
190 public:
191 
192  //- Runtime type information
193  TypeName("ShapeSensitivitiesBase");
194 
195 
196  // Constructors
197 
198  //- Construct from components
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  {
219  return sensitivityPatchIDs_;
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
229  virtual const labelHashSet& geometryVariationIntegrationPatches() const;
230 
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
243  tmp<volVectorField> getWallFaceSensVec();
244 
245  //- Get wall face sensitivity projected to normal field
246  tmp<volScalarField> getWallFaceSensNormal();
247 
248  //- Get wall face normal sens as vectors field
249  tmp<volVectorField> getWallFaceSensNormalVec();
250 
251  //- Get wall point sensitivity vectors field
252  // Uses volPointInterpolation
253  tmp<pointVectorField> getWallPointSensVec();
254 
255  //- Get wall point sensitivity projected to normal field
256  // Uses volPointInterpolation
257  tmp<pointScalarField> getWallPointSensNormal();
258 
259  //- Get wall point sens as vectors field
260  // Uses volPointInterpolation
261  tmp<pointVectorField> getWallPointSensNormalVec();
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 // ************************************************************************* //
tmp< pointScalarField > getWallPointSensNormal()
Get wall point sensitivity projected to normal field.
void writeFaceBasedSens() const
Write face-based sensitivities, if present.
tmp< volVectorField > getWallFaceSensVec()
Get wall face sensitivity vectors field.
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
autoPtr< pointBoundaryVectorField > wallPointSensVecPtr_
Wall point sens w.r.t. (x,y.z)
void allocateEikonalSolver()
Allocate the adjoint eikonal solver.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
labelHashSet sensitivityPatchIDs_
Patches on which to compute shape sensitivities.
Abstract base class for adjoint-based sensitivities.
tmp< volVectorField > getWallFaceSensNormalVec()
Get wall face normal sens as vectors field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
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.
autoPtr< boundaryVectorField > wallFaceSensNormalVecPtr_
Normal sens as vectors.
Base class for adjoint solvers.
Definition: adjointSolver.H:53
autoPtr< boundaryVectorField > wallFaceSensVecPtr_
Wall face sens w.r.t. (x,y.z)
const fvMesh & mesh() const
Return reference to mesh.
Definition: sensitivity.H:121
virtual ~ShapeSensitivitiesBase()=default
Destructor.
autoPtr< boundaryScalarField > wallFaceSensNormalPtr_
Wall face sens projected to normal.
autoPtr< pointBoundaryScalarField > wallPointSensNormalPtr_
Wall point sens projected to normal.
void clearSurfaceFields()
Clear surface/point fields.
bool writeAllSurfaceFiles_
Whether to write all surface sensitivity fields.
virtual const boundaryVectorField & getWallFaceSensNormalVecBoundary() const
Get wall face normal sens as vectors field.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
virtual const boundaryVectorField & getWallFaceSensVecBoundary() const
Get wall face sensitivity vectors field.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void constructAndWriteSensitivityField(const autoPtr< typename GeometricField< Type, fvPatchField, volMesh >::Boundary > &sensFieldPtr, const word &name) const
Constructs volField based on boundaryField and writes it.
Generic templated field type.
Definition: Field.H:62
tmp< pointVectorField > getWallPointSensNormalVec()
Get wall point sens as vectors field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
void constructAndWriteSensitivtyPointField(const autoPtr< List< Field< Type >>> &sensFieldPtr, const word &name) const
Constructs pointField based on boundaryField and writes it.
virtual const boundaryScalarField & getWallFaceSensNormalBoundary() const
Get wall face sensitivity projected to normal field.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
static const word null
An empty word.
Definition: word.H:84
void writePointBasedSens() const
Write point-based sensitivities, if present.
virtual const labelHashSet & geometryVariationIntegrationPatches() const
Return set of patches on which to compute direct sensitivities.
TypeName("ShapeSensitivitiesBase")
Runtime type information.
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: objective.H:55
tmp< volScalarField > getWallFaceSensNormal()
Get wall face sensitivity projected to normal field.
tmp< pointVectorField > getWallPointSensVec()
Get wall point sensitivity vectors field.
autoPtr< pointBoundaryVectorField > wallPointSensNormalVecPtr_
Normal sens as vectors.
Base class supporting Shape sensitivity derivatives.
Useful typenames for fields defined only at the boundaries.
Generic GeometricBoundaryField class.
Definition: areaFieldsFwd.H:46
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.H:129
void clearSensitivities()
Zero sensitivity fields and their constituents.
void setSensitivityPatchIDs(const labelHashSet &sensPatchIDs)
Overwrite sensitivityPatchIDs.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void clearMultipliers()
Clear multipliers.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const labelHashSet & sensitivityPatchIDs() const
Get patch IDs on which sensitivities are computed.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
bool hasMultiplier(bool(objective::*hasFunction)() const)
Check if any of the available objective has a certain multiplier, provided through a function object...
void allocateMultipliers()
Allocate multiplier fields.
Namespace for OpenFOAM.