adjointSensitivity.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::adjointSensitivity
30 
31 Description
32  Abstract base class for adjoint-based sensitivities
33 
34 SourceFiles
35  adjointSensitivity.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef adjointSensitivityIncompressible_H
40 #define adjointSensitivityIncompressible_H
41 
42 #include "boundaryFieldsFwd.H"
43 #include "adjointEikonalSolver.H"
45 #include "sensitivity.H"
46 #include "volFieldsFwd.H"
47 #include "wallFvPatch.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 // Forward declaration
55 class adjointSolver;
56 
57 /*---------------------------------------------------------------------------*\
58  Class adjointSensitivity Declaration
59 \*---------------------------------------------------------------------------*/
60 
62 :
63  public sensitivity
64 {
65 protected:
66 
67 
68  // Protected data
69 
70  //- Reference to the underlaying adjoint solver
72 
73  //- The sensitivity derivative values
75 
76  //- Append this word to files related to the sensitivities
77  word suffix_;
78 
79  //- Include distance variation in sensitivity computations
80  bool includeDistance_;
81 
82  //- Adjoint eikonal equation solver
84 
85  //- Adjoint grid displacement solver
87 
88  // Fields to accumulated through the adjoint solver
89 
90  // Shape optimisation
91 
92  //- Multiplier of grad(dx/b)
94 
95  //- Multiplier of div(dx/db)
97 
98  //- Multiplier of face dx/db
99  // The term that multiplies the adjoint-related part of the
100  // sensitivities in the (E)SI approach
103  //- Multiplier of dSf/db
105 
106  //- Multiplier of dnf/db
108 
109  //- Multiplier of dCf/db, found in the objective function
111 
112  //- Multiplier of dx/db computed at points,
113  //- found in the objective function
116  //- Multiplier of dx/db, coming from boundary conditions that
117  //- depend on the geometry, like rotatingWallVelocity
119 
120  //- dx/db multiplier coming from fvOptions
122 
123 
124 private:
126  // Private Member Functions
127 
128  //- No copy construct
129  adjointSensitivity(const adjointSensitivity&) = delete;
131  //- No copy assignment
132  void operator=(const adjointSensitivity&) = delete;
133 
134 
135 public:
137  //- Runtime type information
138  TypeName("adjointSensitivity");
139 
140 
141  // Declare run-time constructor selection table
144  (
145  autoPtr,
148  (
149  const fvMesh& mesh,
150  const dictionary& dict,
152  ),
153  (
154  mesh,
155  dict,
157  )
158  );
159 
160 
161  // Constructors
162 
163  //- Construct from components
165  (
166  const fvMesh& mesh,
167  const dictionary& dict,
169  );
170 
171  // Selectors
172 
173  //- Return a reference to the selected turbulence model
175  (
176  const fvMesh& mesh,
177  const dictionary& dict,
179  );
180 
181 
182  //- Destructor
183  virtual ~adjointSensitivity() = default;
184 
185 
186  // Member Functions
187 
188  //- Read dictionary if changed
189  virtual bool readDict(const dictionary& dict);
190 
191  //- Const access to adjoint solver
192  inline const adjointSolver& getAdjointSolver() const
193  {
194  return adjointSolver_;
195  }
196 
197  //- Non-const access to adjoint solver
199  {
200  return adjointSolver_;
201  }
202 
203  //- Should the adjoint eikonal PDE should be solved
204  inline bool includeDistance() const
205  {
206  return includeDistance_;
207  }
208 
209  //- Return the adjoint eikonal solver
210  inline autoPtr<adjointEikonalSolver>& getAdjointEikonalSolver()
211  {
212  return eikonalSolver_;
213  }
214 
215  //- Return the adjoint eikonal solver
216  inline autoPtr<adjointMeshMovementSolver>&
218  {
220  }
221 
222  //- Set suffix
223  inline void setSuffix(const word& suffix)
224  {
225  suffix_ = suffix;
226  }
227 
228  //- Get suffix
229  inline const word& getSuffix() const
230  {
231  return suffix_;
232  }
233 
234  //- Should the parameterization compute the internalField of dxdb
235  virtual bool computeDxDbInternalField() const;
236 
237  //- Accumulate sensitivity integrands
238  // Corresponds to the flow and adjoint part of the sensitivities
239  virtual void accumulateIntegrand(const scalar dt) = 0;
240 
241  //- Assemble sensitivities
242  // Adds the geometric part of the sensitivities
243  virtual void assembleSensitivities
244  (
245  autoPtr<designVariables>& designVars
246  );
247 
248  //- Calculates and returns sensitivity fields.
249  // Used with optimisation libraries
251  (
252  autoPtr<designVariables>& designVars
253  );
254 
255  //- Returns the sensitivity fields
256  // Assumes it has already been updated/computed
257  const scalarField& getSensitivities() const;
259  //- Zero sensitivity fields and their constituents
260  virtual void clearSensitivities();
261 
262  //- Write sensitivity fields.
263  // If valid, copies boundaryFields to volFields and writes them.
264  // Virtual to be reimplemented by control points-based methods
265  // (Bezier, RBF) which do not need to write fields
266  virtual void write(const word& baseName = word::null);
268  // Access functions to multipliers
269 
270  // Shape optimisation
271 
272  inline const autoPtr<volTensorField>& gradDxDbMult() const;
274  inline const autoPtr<scalarField>& divDxDbMult() const;
275  inline const autoPtr<boundaryVectorField>& dxdbMult() const;
276  inline const autoPtr<boundaryVectorField>& dSfdbMult() const;
277  inline const autoPtr<boundaryVectorField>& dnfdbMult() const;
278  inline const autoPtr<boundaryVectorField>&
279  dxdbDirectMult() const;
281  pointDxDbDirectMult() const;
282  inline const autoPtr<boundaryVectorField>& bcDxDbMult() const;
283  inline const autoPtr<vectorField>& optionsDxDbMult() const;
284 };
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 #include "adjointSensitivityI.H"
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 #endif
298 
299 // ************************************************************************* //
autoPtr< boundaryVectorField > dSfdbMult_
Multiplier of dSf/db.
adjointSolver & adjointSolver_
Reference to the underlaying adjoint solver.
Forwards and collection of common volume field types.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Abstract base class for adjoint-based sensitivities.
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
Base class for adjoint solvers.
Definition: adjointSolver.H:53
Abstract base class for adjoint sensitivities.
Definition: sensitivity.H:62
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
virtual bool computeDxDbInternalField() const
Should the parameterization compute the internalField of dxdb.
autoPtr< scalarField > divDxDbMult_
Multiplier of div(dx/db)
const fvMesh & mesh() const
Return reference to mesh.
Definition: sensitivity.H:121
const autoPtr< scalarField > & divDxDbMult() const
autoPtr< pointBoundaryVectorField > pointDxDbDirectMult_
Multiplier of dx/db computed at points, found in the objective function.
const autoPtr< volTensorField > & gradDxDbMult() const
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
Return a reference to the selected turbulence model.
autoPtr< volTensorField > gradDxDbMult_
Multiplier of grad(dx/b)
A class for handling words, derived from Foam::string.
Definition: word.H:63
declareRunTimeSelectionTable(autoPtr, adjointSensitivity, dictionary,(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver),(mesh, dict, adjointSolver))
autoPtr< adjointMeshMovementSolver > adjointMeshMovementSolver_
Adjoint grid displacement solver.
const autoPtr< boundaryVectorField > & bcDxDbMult() const
static const word null
An empty word.
Definition: word.H:84
const autoPtr< boundaryVectorField > & dxdbMult() const
void setSuffix(const word &suffix)
Set suffix.
virtual ~adjointSensitivity()=default
Destructor.
autoPtr< boundaryVectorField > dxdbDirectMult_
Multiplier of dCf/db, found in the objective function.
virtual void accumulateIntegrand(const scalar dt)=0
Accumulate sensitivity integrands.
const autoPtr< boundaryVectorField > & dnfdbMult() const
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver()
Return the adjoint eikonal solver.
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
const scalarField & getSensitivities() const
Returns the sensitivity fields.
autoPtr< boundaryVectorField > bcDxDbMult_
Multiplier of dx/db, coming from boundary conditions that depend on the geometry, like rotatingWallVe...
virtual const scalarField & calculateSensitivities(autoPtr< designVariables > &designVars)
Calculates and returns sensitivity fields.
const word & getSuffix() const
Get suffix.
Useful typenames for fields defined only at the boundaries.
autoPtr< adjointMeshMovementSolver > & getAdjointMeshMovementSolver()
Return the adjoint eikonal solver.
word suffix_
Append this word to files related to the sensitivities.
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.H:129
autoPtr< vectorField > optionsDxDbMult_
dx/db multiplier coming from fvOptions
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
bool includeDistance_
Include distance variation in sensitivity computations.
const autoPtr< vectorField > & optionsDxDbMult() const
scalarField derivatives_
The sensitivity derivative values.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const autoPtr< boundaryVectorField > & dxdbDirectMult() const
TypeName("adjointSensitivity")
Runtime type information.
autoPtr< adjointEikonalSolver > eikonalSolver_
Adjoint eikonal equation solver.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool includeDistance() const
Should the adjoint eikonal PDE should be solved.
const autoPtr< boundaryVectorField > & dSfdbMult() const
const autoPtr< pointBoundaryVectorField > & pointDxDbDirectMult() const
Namespace for OpenFOAM.
autoPtr< boundaryVectorField > dxdbMult_
Multiplier of face dx/db.
autoPtr< boundaryVectorField > dnfdbMult_
Multiplier of dnf/db.