Loading...
Searching...
No Matches
adjointEikonalSolver.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-2020 PCOpt/NTUA
9 Copyright (C) 2013-2020 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::adjointEikonalSolver
31
32Description
33 Solver of the adjoint to the eikonal PDE
34
35 Reference:
36 \verbatim
37 For the development of the adjoint eikonal PDE and its boundary
38 conditions
39
40 Papoutsis-Kiachagias, E. M., & Giannakoglou, K. C. (2014).
41 Continuous Adjoint Methods for Turbulent Flows, Applied to Shape
42 and Topology Optimization: Industrial Applications.
43 Archives of Computational Methods in Engineering, 23(2), 255-299.
44 http://doi.org/10.1007/s11831-014-9141-9
45 \endverbatim
46
47 To be as consistent as possible, it is recommended to use the
48 advectionDiffusion wallDist method in fvSchemes, instead of the more widely
49 used meshWave
50
51 Example of the adjointEikonalSolver specification in optimisationDict:
52 \verbatim
53 optimisation
54 {
55 sensitivities
56 {
57 includeDistance true;
58 adjointEikonalSolver
59 {
60 // epsilon should be the same as the one used
61 // in fvSchemes/wallDist/advectionDiffusionCoeffs
62 epsilon 0.1;
63 iters 1000;
64 tolerance 1e-6;
65 }
66 }
67 }
68 \endverbatim
69
70 Example of the entries in fvSchemes:
71 \verbatim
72 divSchemes
73 {
74 .
75 .
76 // avoid bounded schemes since yPhi is not conservative
77 div(-yPhi,da) Gauss linearUpwind grad(da);
78 .
79 .
80 }
81 laplacianSchemes
82 {
83 .
84 .
85 laplacian(yWall,da) Gauss linear corrected;
86 .
87 .
88 }
89 \endverbatim
90
91 Also, the solver specification and a relaxation factor for da are required
92 in fvSolution
93
94 \verbatim
95 da
96 {
97 solver PBiCGStab;
98 preconditioner DILU;
99 tolerance 1e-9;
100 relTol 0.1;
101 }
102
103 relaxationFactors
104 {
105 equations
106 {
107 .
108 .
109 da 0.5;
110 .
111 .
112 }
113 }
114 \endverbatim
115
116See also
117 Foam::patchDistMethod::advectionDiffusion
118 Foam::wallDist
119
120
121SourceFiles
122 adjointEikonalSolver.C
123
124\*---------------------------------------------------------------------------*/
125
126#ifndef adjointEikonalSolver_H
127#define adjointEikonalSolver_H
128
129#include "IOdictionary.H"
130#include "volFieldsFwd.H"
131#include "fvMesh.H"
133#include "createZeroField.H"
134#include "boundaryFieldsFwd.H"
135
136// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137
138namespace Foam
139{
140
141// Forward Declaration
142class adjointSolver;
144/*---------------------------------------------------------------------------*\
145 Class adjointEikonalSolver Declaration
146\*---------------------------------------------------------------------------*/
147
148class adjointEikonalSolver
149{
150private:
151
152 // Private Member Functions
153
154 //- No copy construct
155 adjointEikonalSolver(const adjointEikonalSolver&) = delete;
156
157 //- No copy assignment
158 void operator=(const adjointEikonalSolver&) = delete;
159
160
161protected:
162
163 // Protected data
165 const fvMesh& mesh_;
173 label nEikonalIters_;
175 scalar tolerance_;
177 scalar epsilon_;
184
185 //- Wall face sens w.r.t. (x,y.z)
188
189 // Protected functions
190
191 //- Return the boundary condition types for da
192 wordList patchTypes() const;
193
194 //- Compute convecting velocity
196
197 //- Read options each time a new solution is found
198 void read();
199
200
201
202public:
203
204 //- Runtime type information
205 TypeName("adjointEikonalSolver");
206
207
208 // Constructors
209
210 //- Construct from components
211 adjointEikonalSolver
212 (
213 const fvMesh& mesh,
216 const labelHashSet& sensitivityPatchIDs
217 );
218
219 // Destructor
220
221 virtual ~adjointEikonalSolver() = default;
222
223
224 // Member Functions
225
226 //- Read dict if changed
227 virtual bool readDict(const dictionary& dict);
228
229 //- Accumulate source term
230 void accumulateIntegrand(const scalar dt);
231
232 //- Calculate the adjoint distance field
233 void solve();
234
235 //- Reset source term
236 void reset();
237
238 //- Return the sensitivity term depending on da
240
241 //- Return the volume-based sensitivity term depending on da
243
244 //- Return sensitivity contribution to topology optimisation
245 tmp<scalarField> topologySensitivities(const word& designVarsName) const;
246
247 //- Return the adjoint distance field
248 const volScalarField& da();
249
250 //- Return the distance field
251 const volScalarField& d();
252
253 //- Return the gradient of the eikonal equation
255};
256
257
258// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259
260} // End namespace Foam
261
262// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
263
264#endif
265
266// ************************************************************************* //
Useful typenames for fields defined only at the boundaries.
tmp< volVectorField > gradEikonal()
Return the gradient of the eikonal equation.
boundaryVectorField & distanceSensitivities()
Return the sensitivity term depending on da.
const volScalarField & da()
Return the adjoint distance field.
autoPtr< boundaryVectorField > distanceSensPtr_
Wall face sens w.r.t. (x,y.z).
tmp< surfaceScalarField > computeYPhi()
Compute convecting velocity.
const volScalarField & d()
Return the distance field.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
virtual ~adjointEikonalSolver()=default
const labelHashSet & sensitivityPatchIDs_
tmp< volTensorField > getFISensitivityTerm() const
Return the volume-based sensitivity term depending on da.
void read()
Read options each time a new solution is found.
void accumulateIntegrand(const scalar dt)
Accumulate source term.
TypeName("adjointEikonalSolver")
Runtime type information.
void reset()
Reset source term.
tmp< scalarField > topologySensitivities(const word &designVarsName) const
Return sensitivity contribution to topology optimisation.
wordList patchTypes() const
Return the boundary condition types for da.
void solve()
Calculate the adjoint distance field.
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
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
volVectorField::Boundary boundaryVectorField
dictionary dict
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68
Forwards and collection of common volume field types.