incompressibleAdjointSolver.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  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 
29 Class
30  Foam::incompressibleAdjointSolver
31 
32 Description
33  Base class for incompressibleAdjoint solvers
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef incompressibleAdjointSolver_H
38 #define incompressibleAdjointSolver_H
39 
40 #include "adjointSolver.H"
41 #include "incompressibleVars.H"
43 #include "ATCModel.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class incompressibleAdjointSolver Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 :
56  public adjointSolver
57 {
58 private:
59 
60  // Private Member Functions
61 
62  //- No copy construct
64  (
66  ) = delete;
67 
68  //- No copy assignment
69  void operator=(const incompressibleAdjointSolver&) = delete;
70 
71 
72 protected:
73 
74  // Protected data
75 
76  //- Primal variable set
78 
79  //- Adjoint Transpose Convection options
81 
82  //- Auxiliary bool to avoid a potentially expensive
83  //- part of the sensitivity computation
84  // As a Switch with delayed evaluation since adjointBoundaryConditions
85  // have not been allocated at the time of construction
87 
88 
89  // Protected Member Functions
90 
91  //- Compute, if necessary, and return the hasBCdxdbMult_ bool
92  bool hasBCdxdbMult(const labelHashSet& sensitivityPatchIDs);
93 
94 
95 public:
96 
97 
98  // Static Data Members
99 
100  //- Run-time type information
101  TypeName("incompressible");
102 
103 
104  // Declare run-time constructor selection table
105 
107  (
108  autoPtr,
110  dictionary,
111  (
112  fvMesh& mesh,
113  const word& managerType,
114  const dictionary& dict,
115  const word& primalSolverName,
116  const word& solverName
117  ),
119  );
120 
121 
122  // Constructors
123 
124  //- Construct from mesh and dictionary
126  (
127  fvMesh& mesh,
128  const word& managerType,
129  const dictionary& dict,
130  const word& primalSolverName,
131  const word& solverName
132  );
133 
134 
135  // Selectors
136 
137  //- Return a reference to the selected incompressible adjoint solver
139  (
140  fvMesh& mesh,
141  const word& managerType,
142  const dictionary& dict,
143  const word& primalSolverName,
144  const word& solverName
145  );
146 
147 
148  //- Destructor
149  virtual ~incompressibleAdjointSolver() = default;
150 
151 
152  // Member Functions
153 
154  // Access
155 
156  //- Access to the incompressible primal variables set
157  const incompressibleVars& getPrimalVars() const;
158 
159  //- Access to the incompressible adjoint variables set
160  virtual const incompressibleAdjointVars& getAdjointVars() const;
161 
162  //- Access to the incompressible adjoint variables set
164 
165  //- Access to the ATC model
166  const autoPtr<ATCModel>& getATCModel() const;
167 
168  //- Access to the ATC model
170 
171  //- Should the adjoint to the eikonal equation be solved
172  virtual bool includeDistance() const;
173 
174  //- Return the dimensions of the adjoint distance field
175  virtual dimensionSet daDimensions() const;
176 
177  //- Return the dimensions of the adjoint grid displacement variable
178  virtual dimensionSet maDimensions() const;
179 
180  //- Return the source the adjoint eikonal equation
182 
183  //- Return the distance field, to be used in the solution of the
184  //- adjoint eikonal PDE
185  virtual tmp<volScalarField> yWall() const;
186 
187 
188  // Evolution
189 
190  //- Update primal based quantities, e.g. the primal fields
191  //- in adjoint turbulence models
192  virtual void updatePrimalBasedQuantities();
193 
194 
195  // Sensitivity related functions
196 
197  // Functions related to the computation of sensitivity derivatives
198  // All functions get the field to accumulate their contibution on
199  // as an argument. Should be implemented by the derived classes
200  // if the physics there adds terms to the sensitivity derivatives
201 
202  // Shape optimisation
203 
204  //- Compute the multiplier for grad(dxdb)
205  // Used in shape sensitivity derivatives, computed with
206  // the FI and E-SI approaches
207  virtual void accumulateGradDxDbMultiplier
208  (
209  volTensorField& gradDxDbMult,
210  const scalar dt
211  );
212 
213  //- Compute the multiplier for div(dxdb)
214  // Used in shape sensitivity derivatives, computed with
215  // the FI and E-SI approaches
216  virtual void accumulateDivDxDbMultiplier
217  (
218  autoPtr<scalarField>& divDxDbMult,
219  const scalar dt
220  );
221 
222  //- Accumulate the multipliers of geometric quantities
223  //- defined at the boundary, usually through an objective
224  //- or constraint function
226  (
227  autoPtr<boundaryVectorField>& dSfdbMult,
228  autoPtr<boundaryVectorField>& dnfdbMult,
229  autoPtr<boundaryVectorField>& dxdbDirectMult,
230  autoPtr<pointBoundaryVectorField>& pointDxDirectDbMult,
231  const labelHashSet& sensitivityPatchIDs,
232  const scalar dt
233  );
234 
235  //- Contributions from boundary functions that inlcude
236  //- geometric aspects in them and change when the geometry
237  //- is displaced, e.g. rotationWallVelocity
239  (
240  autoPtr<boundaryVectorField>& bcDxDbMult,
241  const labelHashSet& sensitivityPatchIDs,
242  const scalar dt
243  );
244 
245  //- Contributions from fvOptions that inlcude
246  //- geometric aspects in them and change when the geometry
247  //- is displaced, e.g. MRF
249  (
250  vectorField& optionsDxDbMult,
251  const scalar dt
252  );
253 
254 
255  // Topology optimisation
256 
257  //- Compute the multiplier of beta
258  virtual void topOSensMultiplier
259  (
260  scalarField& betaMult,
261  const word& designVariablesName,
262  const scalar dt
263  );
264 
265 
266  // IO
267 
268  //- In case of multi-point runs with turbulent flows,
269  //- output dummy turbulence fields with the base names, to allow
270  //- continuation
271  virtual bool write(const bool valid = true) const
272  {
273  if (mesh_.time().writeTime())
274  {
275  return primalVars_.write();
276  }
277 
278  return false;
279  }
280 
281  //- In case of multi-point runs with turbulent flows,
282  //- output dummy turbulence fields with the base names, to allow
283  //- continuation
284  virtual bool writeNow() const
285  {
286  return primalVars_.write();
287  }
288 };
289 
290 
291 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
292 
293 } // End namespace Foam
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 #endif
298 
299 // ************************************************************************* //
virtual bool write(const bool valid=true) const
In case of multi-point runs with turbulent flows, output dummy turbulence fields with the base names...
bool hasBCdxdbMult(const labelHashSet &sensitivityPatchIDs)
Compute, if necessary, and return the hasBCdxdbMult_ bool.
const word & primalSolverName() const
Return the primal solver name.
virtual ~incompressibleAdjointSolver()=default
Destructor.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual void accumulateDivDxDbMultiplier(autoPtr< scalarField > &divDxDbMult, const scalar dt)
Compute the multiplier for div(dxdb)
virtual const incompressibleAdjointVars & getAdjointVars() const
Access to the incompressible adjoint variables set.
virtual void accumulateGradDxDbMultiplier(volTensorField &gradDxDbMult, const scalar dt)
Compute the multiplier for grad(dxdb)
const fvMesh & mesh() const
Return the solver mesh.
Definition: solverI.H:24
const incompressibleVars & getPrimalVars() const
Access to the incompressible primal variables set.
autoPtr< ATCModel > ATCModel_
Adjoint Transpose Convection options.
virtual void updatePrimalBasedQuantities()
Update primal based quantities, e.g. the primal fields in adjoint turbulence models.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
Base class for adjoint solvers.
Definition: adjointSolver.H:53
bool writeTime() const noexcept
True if this is a write interval.
Definition: TimeStateI.H:73
virtual bool includeDistance() const
Should the adjoint to the eikonal equation be solved.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class including all adjoint fields for incompressible flows.
Base class for incompressibleAdjoint solvers.
Base class for solution control classes.
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Definition: dimensionSet.H:105
declareRunTimeSelectionTable(autoPtr, incompressibleAdjointSolver, dictionary,(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName),(mesh, managerType, dict, primalSolverName, solverName))
incompressibleVars & primalVars_
Primal variable set.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual void accumulateBCSensitivityIntegrand(autoPtr< boundaryVectorField > &bcDxDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Contributions from boundary functions that inlcude geometric aspects in them and change when the geom...
static autoPtr< incompressibleAdjointSolver > New(fvMesh &mesh, const word &managerType, const dictionary &dict, const word &primalSolverName, const word &solverName)
Return a reference to the selected incompressible adjoint solver.
virtual dimensionSet maDimensions() const
Return the dimensions of the adjoint grid displacement variable.
Switch hasBCdxdbMult_
Auxiliary bool to avoid a potentially expensive part of the sensitivity computation.
const word & managerType() const
Return the manager type.
Definition: solverI.H:72
const dictionary & dict() const
Return the solver dictionary.
Definition: solverI.H:54
const autoPtr< ATCModel > & getATCModel() const
Access to the ATC model.
virtual tmp< volScalarField > yWall() const
Return the distance field, to be used in the solution of the adjoint eikonal PDE. ...
virtual void topOSensMultiplier(scalarField &betaMult, const word &designVariablesName, const scalar dt)
Compute the multiplier of beta.
const word & solverName() const
Return the solver name.
Definition: solverI.H:30
virtual tmp< volScalarField > adjointEikonalSource()
Return the source the adjoint eikonal equation.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual dimensionSet daDimensions() const
Return the dimensions of the adjoint distance field.
virtual bool writeNow() const
In case of multi-point runs with turbulent flows, output dummy turbulence fields with the base names...
fvMesh & mesh_
Reference to the mesh database.
Definition: solver.H:56
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
virtual void accumulateGeometryVariationsMultipliers(autoPtr< boundaryVectorField > &dSfdbMult, autoPtr< boundaryVectorField > &dnfdbMult, autoPtr< boundaryVectorField > &dxdbDirectMult, autoPtr< pointBoundaryVectorField > &pointDxDirectDbMult, const labelHashSet &sensitivityPatchIDs, const scalar dt)
Accumulate the multipliers of geometric quantities defined at the boundary, usually through an object...
bool write() const
Write dummy turbulent fields to allow for continuation in multi-point, turbulent runs.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual void accumulateOptionsDxDbMultiplier(vectorField &optionsDxDbMult, const scalar dt)
Contributions from fvOptions that inlcude geometric aspects in them and change when the geometry is d...
Namespace for OpenFOAM.
TypeName("incompressible")
Run-time type information.