sensitivitySurfacePointsIncompressible.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 -------------------------------------------------------------------------------
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::sensitivitySurfacePoints
30 
31 Description
32  Calculation of adjoint based sensitivities at wall points
33 
34 SourceFiles
35  sensitivitySurfacePoints.C
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #ifndef sensitivitySurfacePointsIncompressible_H
40 #define sensitivitySurfacePointsIncompressible_H
41 
43 #include "shapeSensitivitiesBase.H"
46 #include "deltaBoundary.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 namespace incompressible
54 {
55 
56 /*---------------------------------------------------------------------------*\
57  Class sensitivitySurfacePoints Declaration
58 \*---------------------------------------------------------------------------*/
59 
61 :
62  public adjointSensitivity,
64 {
65 protected:
66 
67  // Protected data
68 
69 
70  //- Include surface area in sens computation
72 
73  //- Include the adjoint pressure term in sens computation
75 
76  //- Include the term containing the grad of the stress at the boundary
78 
79  //- Include the transpose part of the adjoint stresses
81 
82  //- Use snGrad in the transpose part of the adjoint stresses
84 
85  //- Include the term from the deviatoric part of the stresses
86  bool includeDivTerm_;
87 
88  //- Include distance variation in sens computation
89  bool includeDistance_;
90 
91  //- Include mesh movement variation in sens computation
93 
94  //- Include terms directly emerging from the objective function
95  bool includeObjective_;
96 
98 
100 
101  //- The face-based part of the sensitivities
102  // i.e. terms that multiply dxFace/dxPoint.
103  // Sensitivities DO include locale surface area, to get
104  // the correct weighting from the contributions of various faces.
105  // Normalized at the end.
107 
108  //- Multipliers of d(Sf)/db and d(nf)/db
111 
113  // Protected Member Functions
114 
115  //- Read controls and update solver pointers if necessary
116  void read();
117 
118  //- Add terms related to post-processing PDEs
119  //- (i.e. adjoint Eikonal, adjoint mesh movement)
120  //- and add local face area
121  void finaliseFaceMultiplier();
123  //- Converts face sensitivities to point sensitivities and adds the
124  //- ones directly computed in points (i.e. dSf/db and dnf/db).
126 
127  //- Construct globally correct point normals and point areas
129  (
130  vectorField& pointNormals,
131  scalarField& pointMagSf
132  );
133 
134  //- Set suffix name for sensitivity fields
135  void setSuffixName();
136 
137 
138 private:
139 
140  // Private Member Functions
141 
142  //- No copy construct
144 
145  //- No copy assignment
146  void operator=(const sensitivitySurfacePoints&) = delete;
147 
148 
149 public:
150 
151  //- Runtime type information
152  TypeName("surfacePoints");
153 
154 
155  // Constructors
156 
157  //- Construct from components
159  (
160  const fvMesh& mesh,
161  const dictionary& dict,
163  );
164 
165 
166  //- Destructor
167  virtual ~sensitivitySurfacePoints() = default;
168 
169 
170  // Member Functions
171 
172  //- Read dict if changed
173  virtual bool readDict(const dictionary& dict);
174 
175  //- Accumulate sensitivity integrands
176  virtual void accumulateIntegrand(const scalar dt);
177 
178  //- Assemble sensitivities
179  virtual void assembleSensitivities();
180 
181  //- Zero sensitivity fields and their constituents
182  virtual void clearSensitivities();
183 
184  virtual void write(const word& baseName = word::null);
185 };
186 
187 
188 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 
190 } // End namespace incompressible
191 } // End namespace Foam
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 #endif
196 
197 // ************************************************************************* //
bool includeGradStressTerm_
Include the term containing the grad of the stress at the boundary.
virtual void accumulateIntegrand(const scalar dt)
Accumulate sensitivity integrands.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
autoPtr< boundaryVectorField > dSfdbMult_
Multipliers of d(Sf)/db and d(nf)/db.
bool includeTransposeStresses_
Include the transpose part of the adjoint stresses.
Base class for adjoint solvers.
Definition: adjointSolver.H:51
Base class for incompressibleAdjoint solvers.
void finalisePointSensitivities()
Converts face sensitivities to point sensitivities and adds the ones directly computed in points (i...
const dictionary & dict() const
Return the construction dictionary.
Definition: sensitivity.C:50
dynamicFvMesh & mesh
bool includeDivTerm_
Include the term from the deviatoric part of the stresses.
bool useSnGradInTranposeStresses_
Use snGrad in the transpose part of the adjoint stresses.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Abstract base class for adjoint-based sensitivities in incompressible flows.
static const word null
An empty word.
Definition: word.H:84
virtual void clearSensitivities()
Zero sensitivity fields and their constituents.
autoPtr< boundaryVectorField > wallFaceSens_
The face-based part of the sensitivities.
bool includeSurfaceArea_
Include surface area in sens computation.
void finaliseFaceMultiplier()
Add terms related to post-processing PDEs (i.e. adjoint Eikonal, adjoint mesh movement) and add local...
void read()
Read controls and update solver pointers if necessary.
bool includeDistance_
Include distance variation in sens computation.
bool includePressureTerm_
Include the adjoint pressure term in sens computation.
virtual void write(const word &baseName=word::null)
Write sensitivity fields.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:79
Calculation of adjoint based sensitivities at wall points.
bool includeObjective_
Include terms directly emerging from the objective function.
bool includeMeshMovement_
Include mesh movement variation in sens computation.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual bool readDict(const dictionary &dict)
Read dict if changed.
void setSuffixName()
Set suffix name for sensitivity fields.
TypeName("surfacePoints")
Runtime type information.
virtual ~sensitivitySurfacePoints()=default
Destructor.
Base class supporting shape sensitivity derivatives for incompressible flows.
void constructGlobalPointNormalsAndAreas(vectorField &pointNormals, scalarField &pointMagSf)
Construct globally correct point normals and point areas.
Namespace for OpenFOAM.