Loading...
Searching...
No Matches
adjointSolver.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
28
29Class
30 Foam::adjointSolver
31
32Description
33 Base class for adjoint solvers
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef adjointSolver_H
38#define adjointSolver_H
39
40#include "fvMesh.H"
41#include "Time.H"
42#include "IOdictionary.H"
43#include "solver.H"
44#include "objectiveManager.H"
45#include "primalSolver.H"
46#include "adjointSensitivity.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
54class designVariables;
55
56/*---------------------------------------------------------------------------*\
57 Class adjointSolver Declaration
58\*---------------------------------------------------------------------------*/
59
60class adjointSolver
61:
62 public solver
63{
64private:
65
66 // Private Member Functions
67
68 //- No copy construct
69 adjointSolver(const adjointSolver&) = delete;
70
71 //- No copy assignment
72 void operator=(const adjointSolver&) = delete;
73
74
75protected:
76
77 // Protected data
78
79 //- Name of primal solver
81
82 //- Object to manage objective functions
85 //- Sensitivities field
87
88 //- Are sensitivities computed
90
91 //- Is the adjoint solver used to tackle a constraint
92 bool isConstraint_;
93
94 //- Is the adjoint solver used to tackle a double-sided constraint
96
97 //- Sensitivity Derivatives engine
100
101 // Protected Member Functions
102
103 //- Actions to be performed before calculating sensitivities
104 // Does noting in base
105 virtual void preCalculateSensitivities()
106 {}
107
108 //- Allocate the sensitivity derivatives
109 // Since parts of the sensitivities depend on virtual functions
110 // implemented within derived classes, the actual allocation should
111 // happen there
113
114 //- Return the dictionary corresponding to the design variables
116
117
118public:
120 // Static Data Members
121
122 //- Run-time type information
123 TypeName("adjointSolver");
124
125
126 // Declare run-time constructor selection table
127
129 (
130 autoPtr,
131 adjointSolver,
132 adjointSolver,
133 (
134 fvMesh& mesh,
135 const word& managerType,
136 const dictionary& dict,
137 const word& primalSolverName,
138 const word& solverName
139 ),
141 );
142
143
144 // Constructors
145
146 //- Construct from mesh, dictionary, and primal solver name
147 adjointSolver
148 (
150 const word& managerType,
151 const dictionary& dict,
152 const word& primalSolverName,
153 const word& solverName
154 );
155
156
157 // Selectors
158
159 //- Return a reference to the selected turbulence model
161 (
162 fvMesh& mesh,
163 const word& managerType,
164 const dictionary& dict,
165 const word& primalSolverName,
166 const word& solverName
167 );
168
169
170 //- Destructor
171 virtual ~adjointSolver() = default;
172
173
174 // Member Functions
175
176 // Access
177
178 virtual bool readDict(const dictionary& dict);
179
180 //- Return the primal solver name
181 inline const word& primalSolverName() const;
182
183 //- Return a const-reference to the primal solver
184 //- corresponding to this adjoint solver
185 inline const primalSolver& getPrimalSolver() const;
186
187 //- Return a non const-reference to the primal solver
188 //- corresponding to this adjoint solver
190
191 //- Return a const reference to the objective manager
192 inline const objectiveManager& getObjectiveManager() const;
193
194 //- Return a reference to the objective manager
196
197 //- Is the solving referring to a constraint
198 inline bool isConstraint();
199
200 //- Is the solving referring to a double-sided constraint
201 inline bool isDoubleSidedConstraint();
202
203 //- Does the adjoint to an equation computing distances need to
204 //- taken into consideration
205 virtual bool includeDistance() const;
206
207 //- Return the dimensions of the adjoint distance field
208 virtual dimensionSet daDimensions() const;
209
210 //- Return the dimensions of the adjoint grid displacement variable
211 virtual dimensionSet maDimensions() const;
212
213 //- Return the source the adjoint eikonal equation
215
216 //- Return the distance field, to be used in the solution of the
217 //- adjoint eikonal PDE
218 virtual tmp<volScalarField> yWall() const;
219
220
221 // Evolution
222
223 //- Compute sensitivities of the underlaying objectives
225 (
226 autoPtr<designVariables>& designVars
227 );
228
229 //- Grab a reference to the computed sensitivities
231 (
232 autoPtr<designVariables>& designVars
233 );
234
235 //- Clears the sensitivity field known by the adjoint solver
236 virtual void clearSensitivities();
237
238 //- Update primal based quantities, e.g. the primal fields
239 //- in adjoint turbulence models
240 // Does nothing in the base
241 virtual void updatePrimalBasedQuantities();
242
243 //- Write the sensitivity derivatives
244 virtual bool writeData(Ostream& os) const;
245
246 // Functions related to the computation of sensitivity derivatives.
247 // All functions get the field to accumulate their contribution on
248 // as an argument and should be implemented by the derived classes
249
250 // Shape optimisation
251
252 //- Compute the multiplier for grad(dxdb)
253 // Used in shape sensitivity derivatives, computed with
254 // the FI and E-SI approaches
256 (
257 volTensorField& gradDxDbMult,
258 const scalar dt
259 )
260 {}
261
262 //- Compute the multiplier for div(dxdb)
263 // Used in shape sensitivity derivatives, computed with
264 // the FI and E-SI approaches
265 virtual void accumulateDivDxDbMultiplier
266 (
267 autoPtr<scalarField>& divDxDbMult,
268 const scalar dt
269 )
270 {}
271
272 //- Accumulate the multipliers of geometric quantities
273 //- defined at the boundary, usually through an objective
274 //- or constraint function
276 (
277 autoPtr<boundaryVectorField>& dSfdbMult,
278 autoPtr<boundaryVectorField>& dnfdbMult,
279 autoPtr<boundaryVectorField>& dxdbDirectMult,
280 autoPtr<pointBoundaryVectorField>& pointDxDirectDbMult,
281 const labelHashSet& sensitivityPatchIDs,
282 const scalar dt
283 )
284 {}
285
286 //- Contributions from boundary functions that inlcude
287 //- geometric aspects in them and change when the geometry
288 //- is displaced, e.g. rotationWallVelocity
290 (
291 autoPtr<boundaryVectorField>& bcDxDbMult,
292 const labelHashSet& sensitivityPatchIDs,
293 const scalar dt
294 )
295 {}
296
297 //- Contributions from fvOptions that inlcude
298 //- geometric aspects in them and change when the geometry
299 //- is displaced, e.g. MRF
301 (
302 vectorField& optionsDxDbMult,
303 const scalar dt
304 )
305 {}
306
307
308 // Topology optimisation
309
310 //- Compute the multiplier of beta
311 virtual void topOSensMultiplier
312 (
313 scalarField& betaMult,
314 const word& designVariablesName,
315 const scalar dt
316 )
317 {}
318};
319
321// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
322
323} // End namespace Foam
324
325// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
326
327#include "adjointSolverI.H"
328
329// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
330
331#endif
332
333// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual void clearSensitivities()
Clears the sensitivity field known by the adjoint solver.
virtual void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb).
objectiveManager objectiveManager_
Object to manage objective functions.
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE.
tmp< scalarField > sensitivities_
Sensitivities field.
void allocateSensitivities()
Allocate the sensitivity derivatives.
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb).
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
virtual const scalarField & getObjectiveSensitivities(autoPtr< designVariables > &designVars)
Grab a reference to the computed sensitivities.
virtual dimensionSet daDimensions() const
Return the dimensions of the adjoint distance field.
virtual bool writeData(Ostream &os) const
Write the sensitivity derivatives.
static autoPtr< adjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected turbulence model.
bool isDoubleSidedConstraint()
Is the solving referring to a double-sided constraint.
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
autoPtr< adjointSensitivity > adjointSensitivity_
Sensitivity Derivatives engine.
virtual void accumulateGeometryVariationsMultipliers(autoPtr< boundaryVectorField > &dSfdbMult, autoPtr< boundaryVectorField > &dnfdbMult, autoPtr< boundaryVectorField > &dxdbDirectMult, autoPtr< pointBoundaryVectorField > &pointDxDirectDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Accumulate the multipliers of geometric quantities defined at the boundary, usually through an object...
virtual ~adjointSolver()=default
Destructor.
virtual bool readDict(const dictionary &dict)
const primalSolver & getPrimalSolver() const
Return a const-reference to the primal solver corresponding to this adjoint solver.
dictionary designVarsDict() const
Return the dictionary corresponding to the design variables.
virtual void accumulateOptionsDxDbMultiplier(vectorField &optionsDxDbMult, const scalar dt)
Contributions from fvOptions that inlcude geometric aspects in them and change when the geometry is d...
bool computeSensitivities_
Are sensitivities computed.
virtual bool includeDistance() const
Does the adjoint to an equation computing distances need to taken into consideration.
virtual void updatePrimalBasedQuantities()
Update primal based quantities, e.g. the primal fields in adjoint turbulence models.
const word & primalSolverName() const
Return the primal solver name.
TypeName("adjointSolver")
Run-time type information.
virtual void preCalculateSensitivities()
Actions to be performed before calculating sensitivities.
bool isConstraint()
Is the solving referring to a constraint.
virtual void accumulateBCSensitivityIntegrand(autoPtr< boundaryVectorField > &bcDxDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Contributions from boundary functions that inlcude geometric aspects in them and change when the geom...
const word primalSolverName_
Name of primal solver.
declareRunTimeNewSelectionTable(autoPtr, adjointSolver, adjointSolver,(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName),(mesh, managerType, dict, primalSolverName, solverName))
bool isDoubleSidedConstraint_
Is the adjoint solver used to tackle a double-sided constraint.
bool isConstraint_
Is the adjoint solver used to tackle a constraint.
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
const objectiveManager & getObjectiveManager() const
Return a const reference to the objective manager.
virtual void computeObjectiveSensitivities(autoPtr< designVariables > &designVars)
Compute sensitivities of the underlaying objectives.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base class for defining design variables.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Class for managing objective functions.
Base class for primal solvers.
const word & managerType() const
Return the manager type.
Definition solverI.H:72
const dictionary & dict() const
Return the solver dictionary.
Definition solverI.H:54
const fvMesh & mesh() const
Return the solver mesh.
Definition solverI.H:24
const word & solverName() const
Return the solver name.
Definition solverI.H:30
solver(const solver &)=delete
No copy construct.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
OBJstream os(runTime.globalPath()/outputName)
Namespace for OpenFOAM.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Macros to ease declaration of run-time selection tables.
#define declareRunTimeNewSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection for derived classes.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68