adjointSensitivityIncompressible.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-2022 PCOpt/NTUA
9 Copyright (C) 2013-2022 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::incompressible::adjointSensitivity
30
31Description
32 Abstract base class for adjoint-based sensitivities in incompressible flows
33
34 Reference:
35 \verbatim
36 For the FI and ESI formulations
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 For the SI formulation
44 Papoutsis-Kiachagias, E. M., & Giannakoglou, K. C. (2014).
45 Continuous Adjoint Methods for Turbulent Flows, Applied to Shape
46 and Topology Optimization: Industrial Applications.
47 Archives of Computational Methods in Engineering, 23(2), 255–299.
48 http://doi.org/10.1007/s11831-014-9141-9
49 \endverbatim
50
51SourceFiles
52 adjointSensitivity.C
53
54\*---------------------------------------------------------------------------*/
55
56#ifndef adjointSensitivityIncompressible_H
57#define adjointSensitivityIncompressible_H
58
59#include "sensitivity.H"
60#include "incompressibleVars.H"
62#include "wallFvPatch.H"
63
64// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65
66namespace Foam
67{
68
69// Forward declaration
70class incompressibleAdjointSolver;
71
72namespace incompressible
73{
74
75/*---------------------------------------------------------------------------*\
76 Class adjointSensitivity Declaration
77\*---------------------------------------------------------------------------*/
80:
81 public sensitivity
82{
83protected:
84
85 // Protected data
92
93
94private:
95
96 // Private Member Functions
97
98 //- No copy construct
100
101 //- No copy assignment
102 void operator=(const adjointSensitivity&) = delete;
103
104
105public:
106
107 //- Runtime type information
108 TypeName("adjointSensitivity");
109
110
111 // Declare run-time constructor selection table
114 (
115 autoPtr,
118 (
119 const fvMesh& mesh,
120 const dictionary& dict,
122 ),
123 (
124 mesh,
125 dict,
127 )
128 );
129
130
131 // Constructors
132
133 //- Construct from components
135 (
136 const fvMesh& mesh,
137 const dictionary& dict,
139 );
140
141 // Selectors
142
143 //- Return a reference to the selected turbulence model
145 (
146 const fvMesh& mesh,
147 const dictionary& dict,
149 );
150
151
152 //- Destructor
153 virtual ~adjointSensitivity() = default;
154
155
156 // Member Functions
157
158 //- Get primal variables
159 inline const incompressibleVars& primalVars() const
160 {
161 return primalVars_;
162 }
163
164 //- Get adjoint variables
165 inline const incompressibleAdjointVars& adjointVars() const
166 {
167 return adjointVars_;
168 }
169
170 //- Get adjoint solver
171 inline const incompressibleAdjointSolver& adjointSolver() const
172 {
173 return adjointSolver_;
174 }
175
176 //- Accumulate sensitivity integrands
177 // Corresponds to the flow and adjoint part of the sensitivities
178 virtual void accumulateIntegrand(const scalar dt) = 0;
179
180 //- Assemble sensitivities
181 // Adds the geometric part of the sensitivities
182 virtual void assembleSensitivities() = 0;
183
184 //- Calculates and returns sensitivity fields.
185 // Used with optimisation libraries
186 virtual const scalarField& calculateSensitivities();
187
188 //- Returns the sensitivity fields
189 // Assumes it has already been updated/computed
190 const scalarField& getSensitivities() const;
191
192 //- Zero sensitivity fields and their constituents
193 virtual void clearSensitivities();
194
195 //- Write sensitivity fields.
196 // If valid, copies boundaryFields to volFields and writes them.
197 // Virtual to be reimplemented by control points-based methods
198 // (Bezier, RBF) which do not need to write fields
199 virtual void write(const word& baseName = word::null);
200
201 //- Compute the volTensorField multiplying grad(dxdb) for
202 //- the volume-based approach to compute shape sensitivity derivatives
204
205 //- Compute source term for adjoint mesh movement equation
207};
208
209
210// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211
212} // End namespace incompressible
213} // End namespace Foam
214
215// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216
217#endif
218
219// ************************************************************************* //
Base class for adjoint solvers.
Definition: adjointSolver.H:60
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base class for incompressibleAdjoint solvers.
Class including all adjoint fields for incompressible flows.
Base class for solution control classes.
Abstract base class for adjoint-based sensitivities in incompressible flows.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
const incompressibleAdjointSolver & adjointSolver() const
Get adjoint solver.
virtual const scalarField & calculateSensitivities()
Calculates and returns sensitivity fields.
TypeName("adjointSensitivity")
Runtime type information.
virtual void assembleSensitivities()=0
Assemble sensitivities.
virtual ~adjointSensitivity()=default
Destructor.
tmp< volVectorField > adjointMeshMovementSource()
Compute source term for adjoint mesh movement equation.
declareRunTimeSelectionTable(autoPtr, adjointSensitivity, dictionary,(const fvMesh &mesh, const dictionary &dict, incompressibleAdjointSolver &adjointSolver),(mesh, dict, adjointSolver))
virtual void accumulateIntegrand(const scalar dt)=0
Accumulate sensitivity integrands.
const scalarField & getSensitivities() const
Returns the sensitivity fields.
static autoPtr< adjointSensitivity > New(const fvMesh &mesh, const dictionary &dict, incompressibleAdjointSolver &adjointSolver)
Return a reference to the selected turbulence model.
const incompressibleVars & primalVars() const
Get primal variables.
const incompressibleAdjointVars & adjointVars() const
Get adjoint variables.
class for managing incompressible objective functions.
Abstract base class for adjoint sensitivities.
Definition: sensitivity.H:64
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:57
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
dynamicFvMesh & mesh
Namespace for OpenFOAM.
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:73