shapeDesignVariables.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 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::shapeDesignVariables
29 
30 Description
31  Abstract base class for defining design variables for shape optimisation.
32 
33 SourceFiles
34  shapeDesignVariables.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef shapeDesignVariables_H
39 #define shapeDesignVariables_H
40 
41 #include "designVariables.H"
42 #include "displacementMethod.H"
43 #include "runTimeSelectionTables.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class shapeDesignVariables Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 :
56  public designVariables
57 {
58 protected:
59 
60  // Protected data
61 
62  //- Patches to be moved by the design variables
64 
65  //- Mesh movement mechanism
67 
68  //- Store old points. Useful for line search
70 
71  //- Write the mesh points irrespective of whether this is a write time
72  bool writeEachMesh_;
73 
74  // Auxiliary fields keeping track of the various sensitiity components
75 
76  //- Flow related term.
77  // Term including grad(dxdb) in the volume integrals of the FI
78  // formulation
80 
81  //- Flow related term.
82  // Main term in the E-SI formulation. Is the surface intergral
83  // emerging after performing the Gauss-divergence theorem on the
84  // FI-based sensitivities
86 
87  //- Term depending on delta(n dS)/delta b
89 
90  //- Term depending on delta(n)/delta b
92 
93  //- Term depending on delta(x)/delta b for objectives that directly
94  //- depend on x
96 
97  //- Term depending on delta(V)/delta b
99 
100  //- Term depending on distance differentiation
102 
103  //- Term depending on fvOptions
105 
106  //- Term depending on the differenation of boundary conditions
108 
109  //- Name of the sensitivity derivatives folder
111 
112 
113  // Protected Member Functions
114 
115  //- Size of the sensitivity derivatives.
116  // Might be different than this->size() in some cases
117  virtual label sensSize() const;
118 
119  //- Active variables for which to compute sensitivities
120  // Might be different than this->activeDesignVariables_ in some cases
121  virtual const labelList& activeSensitivities() const;
122 
123  //- Compute dxdb at the mesh cell centers by solving a Laplace PDE
125  (
126  const label patchI,
127  const label varID
128  ) const;
129 
130  //- Allocate the fields assosiated with the computation of sensitivities
131  // Not allocated in the constructor since the size of the design
132  // variables is usually not known there
133  void allocateSensFields();
134 
135  //- Zero the fields assosiated with the computation of sensitivities
136  void zeroSensFields();
137 
138 
139 private:
140 
141  // Private Member Functions
142 
143  //- Disallow default bitwise copy construct
145 
146  //- Disallow default bitwise assignment
147  void operator=(const shapeDesignVariables&) = delete;
148 
149 
150 public:
151 
152  //- Runtime type information
153  TypeName("shape");
154 
155 
156  // Declare run-time constructor selection table
157 
159  (
160  autoPtr,
162  dictionary,
163  (
164  fvMesh& mesh,
165  const dictionary& dict
166  ),
167  (mesh, dict)
168  );
169 
170 
171  // Constructors
172 
173  //- Construct from components
175  (
176  fvMesh& mesh,
177  const dictionary& dict
178  );
179 
180 
181  // Selectors
182 
183  //- Construct and return the selected shapeDesignVariables
185  (
186  fvMesh& mesh,
187  const dictionary& dict
188  );
189 
190 
191  //- Destructor
192  virtual ~shapeDesignVariables() = default;
193 
194 
195  // Member Functions
196 
197  //- Read dictionary if changed
198  virtual bool readDict(const dictionary& dict);
199 
200  //- Update design variables based on a given correction.
201  // Translates the scalarField of corrections to a meaningful
202  // update of the design variables
203  virtual void update(scalarField& correction) = 0;
204 
205  //- Store design variables, as the starting point for line search
206  virtual void storeDesignVariables();
207 
208  //- Reset to starting point of line search
209  virtual void resetDesignVariables();
210 
211  //- Compute eta if not set in the first step
212  virtual scalar computeEta(scalarField& correction) = 0;
213 
214  //- Whether to use global sum when computing matrix-vector products
215  //- in update methods
216  // Depends on whether the design variables are common for all
217  // processors (e.g. volumetric B-Splines control points) or distributed
218  // across them (e.g. topology optimisation)
219  virtual bool globalSum() const = 0;
220 
221  // Functions related to mesh movement
222 
223  //- Move mesh based on displacementMethod
224  virtual void moveMesh();
225 
226  //- Patches affected by the parameterisation
227  inline const labelHashSet& getPatches() const
228  {
229  return parametertisedPatches_;
230  }
231 
232  //- Return displacementMethod
234  {
235  return displMethodPtr_;
236  }
237 
238 
239  // Fields related to sensitivity computations
240 
241  //- Add part of sensitivity derivatives related to geometry
242  //- variations
243  virtual tmp<scalarField> assembleSensitivities
244  (
245  adjointSensitivity& adjointSens
246  );
247 
248  //- Write final sensitivity derivatives to files
249  virtual void writeSensitivities
250  (
251  const scalarField& sens,
252  const adjointSensitivity& adjointSens
253  );
254 
255  //- Get dxdb for all mesh points
256  virtual tmp<vectorField> dxdbVol
257  (
258  const label varID
259  ) const;
260 
261  //- Get dxdb for a given design variable and patch
262  virtual tmp<vectorField> dxdbFace
263  (
264  const label patchI,
265  const label varID
266  ) const;
267 
268  //- Get dndb for a given design variable and patch
269  virtual tmp<vectorField> dndb
270  (
271  const label patchI,
272  const label varID
273  ) const;
274 
275  //- Get dSdb for a given design variable and patch
276  virtual tmp<vectorField> dSdb
277  (
278  const label patchI,
279  const label varID
280  ) const;
281 
282  //- Get dCdb for a given design variable.
283  // Used for FI-based sensitivities
284  virtual tmp<volVectorField> dCdb(const label varID) const;
285 };
286 
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 } // End namespace Foam
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 #endif
296 // ************************************************************************* //
virtual scalar computeEta(scalarField &correction)=0
Compute eta if not set in the first step.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix by subtracting the matrix multiplied by the current fi...
scalarField optionsSens_
Term depending on fvOptions.
virtual const labelList & activeSensitivities() const
Active variables for which to compute sensitivities.
bool writeEachMesh_
Write the mesh points irrespective of whether this is a write time.
dictionary dict
virtual tmp< vectorField > dxdbFace(const label patchI, const label varID) const
Get dxdb for a given design variable and patch.
A class for handling file names.
Definition: fileName.H:72
virtual void writeSensitivities(const scalarField &sens, const adjointSensitivity &adjointSens)
Write final sensitivity derivatives to files.
autoPtr< displacementMethod > displMethodPtr_
Mesh movement mechanism.
virtual void storeDesignVariables()
Store design variables, as the starting point for line search.
scalarField dxdbDirectSens_
Term depending on delta(x)/delta b for objectives that directly depend on x.
virtual tmp< volVectorField > dCdb(const label varID) const
Get dCdb for a given design variable.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
void zeroSensFields()
Zero the fields assosiated with the computation of sensitivities.
virtual void resetDesignVariables()
Reset to starting point of line search.
virtual void moveMesh()
Move mesh based on displacementMethod.
virtual tmp< scalarField > assembleSensitivities(adjointSensitivity &adjointSens)
Add part of sensitivity derivatives related to geometry variations.
virtual label sensSize() const
Size of the sensitivity derivatives.
virtual bool globalSum() const =0
Whether to use global sum when computing matrix-vector products in update methods.
virtual tmp< volVectorField > solveMeshMovementEqn(const label patchI, const label varID) const
Compute dxdb at the mesh cell centers by solving a Laplace PDE.
static autoPtr< shapeDesignVariables > New(fvMesh &mesh, const dictionary &dict)
Construct and return the selected shapeDesignVariables.
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< vectorField > dndb(const label patchI, const label varID) const
Get dndb for a given design variable and patch.
declareRunTimeSelectionTable(autoPtr, shapeDesignVariables, dictionary,(fvMesh &mesh, const dictionary &dict),(mesh, dict))
scalarField dndbSens_
Term depending on delta(n)/delta b.
scalarField dVdbSens_
Term depending on delta(V)/delta b.
scalarField bcSens_
Term depending on the differenation of boundary conditions.
scalarField dSdbSens_
Term depending on delta(n dS)/delta b.
virtual void update(scalarField &correction)=0
Update design variables based on a given correction.
scalarField dxdbSurfSens_
Flow related term.
fileName derivativesFolder_
Name of the sensitivity derivatives folder.
virtual bool readDict(const dictionary &dict)
Read dictionary if changed.
const labelHashSet & getPatches() const
Patches affected by the parameterisation.
Abstract base class for defining design variables.
TypeName("shape")
Runtime type information.
Abstract base class for defining design variables for shape optimisation.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual tmp< vectorField > dSdb(const label patchI, const label varID) const
Get dSdb for a given design variable and patch.
void allocateSensFields()
Allocate the fields assosiated with the computation of sensitivities.
autoPtr< pointField > pointsInit_
Store old points. Useful for line search.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Macros to ease declaration of run-time selection tables.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
scalarField dxdbVolSens_
Flow related term.
scalarField distanceSens_
Term depending on distance differentiation.
autoPtr< displacementMethod > & returnDisplacementMethod()
Return displacementMethod.
virtual ~shapeDesignVariables()=default
Destructor.
labelHashSet parametertisedPatches_
Patches to be moved by the design variables.
virtual tmp< vectorField > dxdbVol(const label varID) const
Get dxdb for all mesh points.
Namespace for OpenFOAM.