sensitivitySurfaceIncompressible.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-2021 PCOpt/NTUA
9  Copyright (C) 2013-2021 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::incompressible::sensitivitySurface
30 
31 Description
32  Calculation of adjoint based sensitivities at wall faces
33 
34  Reference:
35  \verbatim
36  The computation of the sensitivity derivatives follows the (E-)SI
37  formulation of
38 
39  Kavvadias, I. S., Papoutsis-Kiachagias, E. M.,
40  Giannakoglou, K. C. (2015).
41  On the proper treatment of grid sensitivities in continuous adjoint
42  methods for shape optimization.
43  Journal of computational physics, 301, 1-18.
44  https://doi.org/10.1016/j.jcp.2015.08.012
45 
46  whereas their smoothing based on the computation of the 'Sobolev
47  gradient' is derived from
48 
49  Vassberg J. C., Jameson A. (2006).
50  Aerodynamic Shape Optimization Part I: Theoretical Background.
51  VKI Lecture Series, Introduction to Optimization and
52  Multidisciplinary Design, Brussels, Belgium, 8 March, 2006.
53  \endverbatim
54 
55 SourceFiles
56  sensitivitySurface.C
57 
58 \*---------------------------------------------------------------------------*/
59 
60 #ifndef sensitivitySurfaceIncompressible_H
61 #define sensitivitySurfaceIncompressible_H
62 
64 #include "shapeSensitivitiesBase.H"
67 #include "deltaBoundary.H"
68 #include "faMesh.H"
69 
70 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72 namespace Foam
73 {
74 
75 namespace incompressible
76 {
77 
78 /*---------------------------------------------------------------------------*\
79  Class sensitivitySurface Declaration
80 \*---------------------------------------------------------------------------*/
81 
83 :
84  public adjointSensitivity,
86 {
87 protected:
88 
89  // Protected data
90 
91 
92  //- Include surface area in sens computation
94 
95  //- Include the adjoint pressure term in sens computation
97 
98  //- Include the term containing the grad of the stress at the boundary
101  //- Include the transpose part of the adjoint stresses
103 
104  //- Use snGrad in the transpose part of the adjoint stresses
106 
107  //- Include the term from the deviatoric part of the stresses
108  bool includeDivTerm_;
109 
110  //- Include distance variation in sens computation
111  bool includeDistance_;
112 
113  //- Include mesh movement variation in sens computation
116  //- Include terms directly emerging from the objective function
117  bool includeObjective_;
118 
119  //- Write geometric info for use by external programs
121 
122  //- Smooth sensitivity derivatives based on the computation of the
123  //- 'Sobolev gradient'
127 
129 
130  // Export face normal and face centre for use by external users
134 
136  // Protected Member Functions
137 
138  //- Add sensitivities from dSd/db and dnf/db computed at points and
139  //- mapped to faces
140  void addGeometricSens();
142  //- Set suffix name for sensitivity fields
144 
145  //- Smooth sensitivity derivatives based on the computation of the
146  //- 'Sobolev gradient'
147  void smoothSensitivities();
149  //- Compute the physical smoothing radius based on the average boundary
150  //- face 'length'
151  scalar computeRadius(const faMesh& aMesh);
152 
153 
154 private:
155 
156  // Private Member Functions
157 
158  //- No copy construct
159  sensitivitySurface(const sensitivitySurface&) = delete;
160 
161  //- No copy assignment
162  void operator=(const sensitivitySurface&) = delete;
163 
164 
165 public:
166 
167  //- Runtime type information
168  TypeName("surface");
169 
170 
171  // Constructors
172 
173  //- Construct from components
175  (
176  const fvMesh& mesh,
177  const dictionary& dict,
179  );
180 
181 
182  //- Destructor
183  virtual ~sensitivitySurface() = default;
184 
185 
186  // Member Functions
187 
188  //- Read controls and update solver pointers if necessary
189  void read();
190 
191  //- Read dict if changed
192  virtual bool readDict(const dictionary& dict);
193 
194  //- Compute the number of faces on sensitivityPatchIDs_
195  void computeDerivativesSize();
196 
197  //- Accumulate sensitivity integrands
198  virtual void accumulateIntegrand(const scalar dt);
199 
200  //- Assemble sensitivities
201  virtual void assembleSensitivities();
202 
203  //- Zero sensitivity fields and their constituents
204  virtual void clearSensitivities();
205 
206  //- Get adjoint eikonal solver
208 
209  //- Write sensitivity maps
210  virtual void write(const word& baseName = word::null);
211 
212  // Inline getters and setters
213 
214  //- Get access to the includeObjective bool
215  inline bool getIncludeObjective() const;
216 
217  //- Get access to the includeSurfaceArea bool
218  inline bool getIncludeSurfaceArea() const;
219 
220  //- Set includeObjective bool
221  inline void setIncludeObjective(const bool includeObjective);
222 
223  //- Set includeSurfaceArea bool
224  inline void setIncludeSurfaceArea(const bool includeSurfaceArea);
225 
226 };
227 
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 } // End namespace incompressible
232 } // End namespace Foam
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 #endif
241 
242 // ************************************************************************* //
Finite area mesh (used for 2-D non-Euclidian finite area method) defined using a patch of faces on a ...
Definition: faMesh.H:88
void setIncludeObjective(const bool includeObjective)
Set includeObjective bool.
bool useSnGradInTranposeStresses_
Use snGrad in the transpose part of the adjoint stresses.
bool includeObjective_
Include terms directly emerging from the objective function.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
autoPtr< adjointMeshMovementSolver > meshMovementSolver_
virtual void write(const word &baseName=word::null)
Write sensitivity maps.
void addGeometricSens()
Add sensitivities from dSd/db and dnf/db computed at points and mapped to faces.
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
void smoothSensitivities()
Smooth sensitivity derivatives based on the computation of the &#39;Sobolev gradient&#39;.
Base class for adjoint solvers.
Definition: adjointSolver.H:51
bool includeTransposeStresses_
Include the transpose part of the adjoint stresses.
bool getIncludeSurfaceArea() const
Get access to the includeSurfaceArea bool.
void read()
Read controls and update solver pointers if necessary.
bool includeDivTerm_
Include the term from the deviatoric part of the stresses.
virtual ~sensitivitySurface()=default
Destructor.
bool includeGradStressTerm_
Include the term containing the grad of the stress at the boundary.
Base class for incompressibleAdjoint solvers.
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:50
void setIncludeSurfaceArea(const bool includeSurfaceArea)
Set includeSurfaceArea bool.
TypeName("surface")
Runtime type information.
dynamicFvMesh & mesh
autoPtr< adjointEikonalSolver > & getAdjointEikonalSolver()
Get adjoint eikonal solver.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Abstract base class for adjoint-based sensitivities in incompressible flows.
virtual bool readDict(const dictionary &dict)
Read dict if changed.
bool includeSurfaceArea_
Include surface area in sens computation.
void computeDerivativesSize()
Compute the number of faces on sensitivityPatchIDs_.
bool includeMeshMovement_
Include mesh movement variation in sens computation.
static const word null
An empty word.
Definition: word.H:84
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
virtual void assembleSensitivities()
Assemble sensitivities.
bool writeGeometricInfo_
Write geometric info for use by external programs.
bool includeDistance_
Include distance variation in sens computation.
bool includePressureTerm_
Include the adjoint pressure term in sens computation.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Calculation of adjoint based sensitivities at wall faces.
void setSuffixName()
Set suffix name for sensitivity fields.
scalar computeRadius(const faMesh &aMesh)
Compute the physical smoothing radius based on the average boundary face &#39;length&#39;.
bool smoothSensitivities_
Smooth sensitivity derivatives based on the computation of the &#39;Sobolev gradient&#39;.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
bool getIncludeObjective() const
Get access to the includeObjective bool.
faMesh aMesh(mesh)
Base class supporting shape sensitivity derivatives for incompressible flows.
Namespace for OpenFOAM.