Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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::adjointSensitivity
30
31Description
32 Abstract base class for adjoint-based sensitivities
33
34SourceFiles
35 adjointSensitivity.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef adjointSensitivityIncompressible_H
40#define adjointSensitivityIncompressible_H
41
42#include "boundaryFieldsFwd.H"
45#include "sensitivity.H"
46#include "volFieldsFwd.H"
47#include "wallFvPatch.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53
54// Forward declaration
55class adjointSolver;
57/*---------------------------------------------------------------------------*\
58 Class adjointSensitivity Declaration
59\*---------------------------------------------------------------------------*/
60
61class adjointSensitivity
62:
63 public sensitivity
64{
65protected:
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
79 //- Include distance variation in sensitivity computations
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
124private:
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
135public:
137 //- Runtime type information
138 TypeName("adjointSensitivity");
139
140
141 // Declare run-time constructor selection table
144 (
145 autoPtr,
146 adjointSensitivity,
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
164 adjointSensitivity
165 (
166 const fvMesh& mesh,
167 const dictionary& dict,
169 );
171 // Selectors
172
173 //- Return a reference to the selected turbulence model
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 }
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;
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;
284};
285
286
287// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288
289} // End namespace Foam
290
291// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292
293#include "adjointSensitivityI.H"
294
295// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296
297#endif
299// ************************************************************************* //
Useful typenames for fields defined only at the boundaries.
autoPtr< scalarField > divDxDbMult_
Multiplier of div(dx/db).
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver()
Return the adjoint eikonal solver.
bool includeDistance_
Include distance variation in sensitivity computations.
const autoPtr< volTensorField > & gradDxDbMult() const
autoPtr< boundaryVectorField > dxdbDirectMult_
Multiplier of dCf/db, found in the objective function.
virtual const scalarField & calculateSensitivities(autoPtr< designVariables > &designVars)
Calculates and returns sensitivity fields.
const autoPtr< vectorField > & optionsDxDbMult() const
declareRunTimeSelectionTable(autoPtr, adjointSensitivity, dictionary,(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver),(mesh, dict, adjointSolver))
autoPtr< pointBoundaryVectorField > pointDxDbDirectMult_
Multiplier of dx/db computed at points, found in the objective function.
void setSuffix(const word &suffix)
Set suffix.
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...
TypeName("adjointSensitivity")
Runtime type information.
autoPtr< adjointMeshMovementSolver > & getAdjointMeshMovementSolver()
Return the adjoint eikonal solver.
adjointSolver & adjointSolver_
Reference to the underlaying adjoint solver.
autoPtr< vectorField > optionsDxDbMult_
dx/db multiplier coming from fvOptions
virtual ~adjointSensitivity()=default
Destructor.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
bool includeDistance() const
Should the adjoint eikonal PDE should be solved.
const autoPtr< boundaryVectorField > & dxdbDirectMult() const
autoPtr< boundaryVectorField > dxdbMult_
Multiplier of face dx/db.
const adjointSolver & getAdjointSolver() const
Const access to adjoint solver.
virtual void accumulateIntegrand(const scalar dt)=0
Accumulate sensitivity integrands.
adjointSolver & getAdjointSolver()
Non-const access to adjoint solver.
const autoPtr< boundaryVectorField > & dnfdbMult() const
virtual bool computeDxDbInternalField() const
Should the parameterization compute the internalField of dxdb.
scalarField derivatives_
The sensitivity derivative values.
const autoPtr< boundaryVectorField > & dxdbMult() const
const scalarField & getSensitivities() const
Returns the sensitivity fields.
const autoPtr< scalarField > & divDxDbMult() const
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, adjointSolver &adjointSolver)
Return a reference to the selected turbulence model.
const autoPtr< pointBoundaryVectorField > & pointDxDbDirectMult() const
autoPtr< boundaryVectorField > dnfdbMult_
Multiplier of dnf/db.
const autoPtr< boundaryVectorField > & dSfdbMult() const
autoPtr< adjointMeshMovementSolver > adjointMeshMovementSolver_
Adjoint grid displacement solver.
autoPtr< boundaryVectorField > dSfdbMult_
Multiplier of dSf/db.
const word & getSuffix() const
Get suffix.
autoPtr< volTensorField > gradDxDbMult_
Multiplier of grad(dx/b).
virtual void assembleSensitivities(autoPtr< designVariables > &designVars)
Assemble sensitivities.
const autoPtr< boundaryVectorField > & bcDxDbMult() const
word suffix_
Append this word to files related to the sensitivities.
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
const fvMesh & mesh() const
Return reference to mesh.
const dictionary & dict() const
Return the construction dictionary.
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.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
runTime write()
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68
Forwards and collection of common volume field types.