adjointRASModel.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 
29 Namespace
30  Foam::incompressible::adjointRASModels
31 
32 Description
33  Namespace for incompressible adjointRAS turbulence models.
34 
35 Class
36  Foam::incompressibleAdjoint::adjointRASModel
37 
38 Description
39  Abstract base class for incompressible turbulence models.
40 
41 SourceFiles
42  adjointRASModel.C
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #ifndef adjointRASModel_H
47 #define adjointRASModel_H
48 
49 #include "adjointTurbulenceModel.H"
50 #include "volFields.H"
51 #include "surfaceFields.H"
52 #include "nearWallDist.H"
53 #include "fvm.H"
54 #include "fvc.H"
55 #include "fvMatrices.H"
57 #include "IOdictionary.H"
58 #include "Switch.H"
59 #include "bound.H"
60 #include "autoPtr.H"
61 #include "runTimeSelectionTables.H"
62 #include "objectiveManager.H"
63 #include "boundaryFieldsFwd.H"
64 #include "createZeroField.H"
65 #include "solverControl.H"
66 
67 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69 namespace Foam
70 {
71 namespace incompressibleAdjoint
72 {
73 
74 /*---------------------------------------------------------------------------*\
75  Class adjointRASModel Declaration
76 \*---------------------------------------------------------------------------*/
77 
78 class adjointRASModel
79 :
81  public IOdictionary
82 {
83 private:
84 
85  // Private Member Functions
86 
87  //- No copy construct
88  adjointRASModel(const adjointRASModel&) = delete;
89 
90  //- No copy assignment
91  void operator=(const adjointRASModel&) = delete;
92 
93 
94 protected:
95 
96  // Protected data
97 
98  //- Reference to the objectiveManager
100 
101  //- Turbulence on/off flag
103 
104  //- Flag to print the model coeffs at run-time
106 
107  //- Model coefficients dictionary
109 
110  //- Near wall distance boundary field
112 
113  //- Adjoint turbulence model variable 1
115 
116  //- Adjoint turbulence model variable 2
118 
119  //- Adjoint turbulence model variable 1, mean value
121 
122  //- Adjoint turbulence model variable 2, mean value
124 
125  //- Source to the adjoint momentum BC emerging
126  //- from differentiating the turbulence model
128 
129  //- Wall sensitivity term for shape optimisation
131 
132  //- Wall sensitivity term for flow control optimisation
134 
135  //- Does the turbulence model include distances and should the
136  //- adjoint to the distance field be computed
137  bool includeDistance_;
138 
139  //- Has the primal solution changed?
141 
142 
143  // Protected Member Functions
144 
145  //- Print model coefficients
146  virtual void printCoeffs();
147 
148  //- Set mean fields
149  void setMeanFields();
150 
151 
152 public:
153 
154  //- Runtime type information
155  TypeName("adjointRASModel");
156 
157 
158  // Declare run-time constructor selection table
159 
161  (
162  autoPtr,
164  dictionary,
165  (
166  incompressibleVars& primalVars,
168  objectiveManager& objManager,
169  const word& adjointTurbulenceModelName
170  ),
171  (
172  primalVars,
173  adjointVars,
174  objManager,
175  adjointTurbulenceModelName
176  )
177  );
178 
179 
180  // Constructors
181 
182  //- Construct from components
184  (
185  const word& type,
186  incompressibleVars& primalVars,
188  objectiveManager& objManager,
189  const word& adjointTurbulenceModelName =
190  adjointTurbulenceModel::typeName
191  );
192 
193 
194  // Selectors
195 
196  //- Return a reference to the selected adjointRAS model
198  (
199  incompressibleVars& primalVars,
201  objectiveManager& objManager,
202  const word& adjointTurbulenceModelName =
203  adjointTurbulenceModel::typeName
204  );
205 
206 
207  //- Destructor
208  virtual ~adjointRASModel() = default;
209 
210 
211  // Member Functions
212 
213  //- Return the near wall distances
214  const nearWallDist& y() const
215  {
216  return y_;
217  }
218 
219  //- Const access to the coefficients dictionary
220  const dictionary& coeffDict() const
221  {
222  return coeffDict_;
223  }
224 
225  //- Const access to the primal solver name
226  const word& primalSolverName() const
227  {
228  return primalVars_.solverName();
229  }
230 
231  //- Const access to the adjoint solver name
232  const word& adjointSolverName() const
233  {
234  return adjointVars_.solverName();
235  }
236 
237  //- Return non-constant reference to adjoint turbulence model variable 1
238  // Will allocate and return a zero field in case it does not exist
240 
241  //- Return non-constant reference to adjoint turbulence model variable 2
242  // Will allocate and return a zero field in case it does not exist
244 
245  //- Return non-constant reference to adjoint turbulence model variable 1
246  // Will return the mean value if present,
247  // otherwise the instantaneous value
249 
250  //- Return non-constant reference to adjoint turbulence model variable 2
251  // Will return the mean value if present,
252  // otherwise the instantaneous value
254 
255  //- Return non-constant autoPtr to adjoint turbulence model variable 1
257 
258  //- Return non-constant autoPtr to adjoint turbulence model variable 2
260 
261  //- Return the effective stress tensor including the laminar stress
262  virtual tmp<volSymmTensorField> devReff() const = 0;
263 
264  //- Return the effective stress tensor based on a given velocity field
266  (
267  const volVectorField& U
268  ) const = 0;
269 
270  //- Return the diffusion term for the momentum equation
271  virtual tmp<fvVectorMatrix> divDevReff(volVectorField& U) const = 0;
272 
273  //- Source terms to the adjoint momentum equation due to
274  //- the differentiation of the turbulence model
276 
277  //- Jacobian of nut wrt the first turbulence model variable
278  // Needed for objective functions that depend on nut. Defaults to zero
279  virtual tmp<volScalarField> nutJacobianTMVar1() const;
280 
281  //- Jacobian of nut wrt the second turbulence model variable
282  // Needed for objective functions that depend on nut. Defaults to zero
283  virtual tmp<volScalarField> nutJacobianTMVar2() const;
284 
285  //- Diffusion coefficient of the first primal and adjoint turbulence
286  //- model equation. Needed for some adjoint BCs. Defaults to zero
287  virtual tmp<scalarField> diffusionCoeffVar1(label patchI) const;
288 
289  //- Diffusion coefficient of the second primal and adjoint turbulence
290  //- model equation. Needed for some adjoint BCs. Defaults to zero
291  virtual tmp<scalarField> diffusionCoeffVar2(label patchI) const;
292 
293  //- Source for the outlet adjoint momentum BC coming from
294  //- differentiating the turbulence model
295  virtual const boundaryVectorField& adjointMomentumBCSource() const = 0;
296 
297  //- Sensitivity terms for shape optimisation, emerging from
298  // the turbulence model differentiation.
299  // Misses dxdb, to be added by the classes assembling the sensitivities
300  virtual const boundaryVectorField& wallShapeSensitivities() = 0;
301 
302  //- Sensitivity terms for flow control, emerging from the
303  // turbulence model differentiation
304  virtual const boundaryVectorField& wallFloCoSensitivities() = 0;
305 
306  //- Sensitivity terms resulting from the differentiation of the
307  //- distance field. Misses dxdb, to be added by the classes
308  //- assembling the sensitivities
310 
311  //- Term contributing to the computation of FI-based sensitivities
312  // Misses grad(dxdb), to be added by the assembling the sensitivities
314 
315  //- Solve the adjoint turbulence equations
316  virtual void correct();
317 
318  //- Read adjointRASProperties dictionary
319  virtual bool read();
320 
321  //- Set flag of changed primal solution to true
323 
324  //- Reset mean fields to zero
325  void resetMeanFields();
326 
327  //- Average adjoint fields on the fly
328  void computeMeanFields();
329 
330  //- Should the adjoint to the eikonal equation be computed
331  bool includeDistance() const;
332 
333  //- Nullify all adjoint turbulence model fields and their old times
334  virtual void nullify() = 0;
335 };
336 
337 
338 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
339 
340 } // End namespace incompressible
341 } // End namespace Foam
342 
343 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
344 
345 #endif
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::incompressibleAdjointMeanFlowVars
Manages the adjoint mean flow fields and their mean values.
Definition: incompressibleAdjointMeanFlowVars.H:51
Foam::objectiveManager
class for managing incompressible objective functions.
Definition: objectiveManager.H:54
Foam::variablesSet::solverName
const word & solverName() const
Return solver name.
Definition: variablesSet.C:84
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:77
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable1
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
Definition: adjointRASModel.C:315
Foam::incompressibleAdjoint::adjointRASModel::nutJacobianTMVar2
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
Definition: adjointRASModel.C:377
Foam::incompressibleAdjoint::adjointRASModel::wallShapeSensitivitiesPtr_
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
Definition: adjointRASModel.H:129
Foam::incompressibleAdjoint::adjointRASModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: adjointRASModel.H:213
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
adjointTurbulenceModel.H
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::incompressibleAdjoint::adjointRASModel::nutJacobianTMVar1
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
Definition: adjointRASModel.C:354
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable2MeanPtr_
autoPtr< volScalarField > adjointTMVariable2MeanPtr_
Adjoint turbulence model variable 2, mean value.
Definition: adjointRASModel.H:122
Foam::incompressibleAdjoint::adjointRASModel::coeffDict_
dictionary coeffDict_
Model coefficients dictionary.
Definition: adjointRASModel.H:107
Foam::incompressibleAdjoint::adjointRASModel::wallFloCoSensitivities
virtual const boundaryVectorField & wallFloCoSensitivities()=0
Sensitivity terms for flow control, emerging from the.
Foam::incompressibleAdjoint::adjointRASModel::correct
virtual void correct()
Solve the adjoint turbulence equations.
Definition: adjointRASModel.C:215
Foam::incompressibleAdjoint::adjointRASModel::wallFloCoSensitivitiesPtr_
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
Definition: adjointRASModel.H:132
surfaceFields.H
Foam::surfaceFields.
Foam::incompressibleAdjoint::adjointRASModel::changedPrimalSolution_
bool changedPrimalSolution_
Has the primal solution changed?
Definition: adjointRASModel.H:139
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable1MeanPtr_
autoPtr< volScalarField > adjointTMVariable1MeanPtr_
Adjoint turbulence model variable 1, mean value.
Definition: adjointRASModel.H:119
objectiveManager.H
Foam::incompressibleAdjoint::adjointRASModel::diffusionCoeffVar2
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Definition: adjointRASModel.C:406
Foam::incompressibleAdjoint::adjointRASModel::adjointSolverName
const word & adjointSolverName() const
Const access to the adjoint solver name.
Definition: adjointRASModel.H:231
Foam::incompressibleAdjoint::adjointRASModel::y_
nearWallDist y_
Near wall distance boundary field.
Definition: adjointRASModel.H:110
Foam::incompressibleAdjoint::adjointRASModel::read
virtual bool read()
Read adjointRASProperties dictionary.
Definition: adjointRASModel.C:226
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable2Ptr_
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
Definition: adjointRASModel.H:116
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable2Inst
volScalarField & getAdjointTMVariable2Inst()
Return non-constant reference to adjoint turbulence model variable 2.
Definition: adjointRASModel.C:288
nearWallDist.H
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable2
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
Definition: adjointRASModel.C:329
Foam::incompressibleAdjoint::adjointRASModel::nullify
virtual void nullify()=0
Nullify all adjoint turbulence model fields and their old times.
createZeroField.H
Foam::incompressibleAdjoint::adjointRASModel::adjMomentumBCSourcePtr_
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Definition: adjointRASModel.H:126
bound.H
Bound the given scalar field if it has gone unbounded.
Switch.H
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable2InstPtr
autoPtr< volScalarField > & getAdjointTMVariable2InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 2.
Definition: adjointRASModel.C:348
Foam::incompressibleAdjoint::adjointRASModel::declareRunTimeSelectionTable
declareRunTimeSelectionTable(autoPtr, adjointRASModel, dictionary,(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName),(primalVars, adjointVars, objManager, adjointTurbulenceModelName))
Foam::incompressibleAdjoint::adjointTurbulenceModel
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
Definition: adjointTurbulenceModel.H:71
transportModel.H
Foam::incompressibleAdjoint::adjointRASModel::setChangedPrimalSolution
void setChangedPrimalSolution()
Set flag of changed primal solution to true.
Definition: adjointRASModel.C:412
Foam::incompressibleAdjoint::adjointRASModel::primalSolverName
const word & primalSolverName() const
Const access to the primal solver name.
Definition: adjointRASModel.H:225
Foam::incompressibleAdjoint::adjointRASModel::FISensitivityTerm
virtual tmp< volTensorField > FISensitivityTerm()=0
Term contributing to the computation of FI-based sensitivities.
Foam::incompressibleAdjoint::adjointRASModel::distanceSensitivities
virtual tmp< volScalarField > distanceSensitivities()=0
fvm.H
Foam::incompressibleAdjoint::adjointRASModel::adjointTurbulence_
Switch adjointTurbulence_
Turbulence on/off flag.
Definition: adjointRASModel.H:101
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::incompressibleAdjoint::adjointRASModel::includeDistance
bool includeDistance() const
Should the adjoint to the eikonal equation be computed.
Definition: adjointRASModel.C:462
Foam::incompressibleAdjoint::adjointRASModel
Abstract base class for incompressible turbulence models.
Definition: adjointRASModel.H:77
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressibleAdjoint::adjointRASModel::objectiveManager_
objectiveManager & objectiveManager_
Reference to the objectiveManager.
Definition: adjointRASModel.H:98
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable1Inst
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
Definition: adjointRASModel.C:261
IOdictionary.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::incompressibleAdjoint::adjointRASModel::TypeName
TypeName("adjointRASModel")
Runtime type information.
Foam::GeometricField::Boundary
The boundary fields.
Definition: GeometricField.H:115
boundaryFieldsFwd.H
Useful typenames for fields defined only at the boundaries.
runTimeSelectionTables.H
Macros to ease declaration of run-time selection tables.
Foam::incompressibleAdjoint::adjointRASModel::adjointMeanFlowSource
virtual tmp< volVectorField > adjointMeanFlowSource()=0
Foam::incompressibleAdjoint::adjointRASModel::adjointTMVariable1Ptr_
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
Definition: adjointRASModel.H:113
Foam::incompressibleAdjoint::adjointRASModel::printCoeffs
virtual void printCoeffs()
Print model coefficients.
Definition: adjointRASModel.C:55
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::incompressibleAdjoint::adjointRASModel::getAdjointTMVariable1InstPtr
autoPtr< volScalarField > & getAdjointTMVariable1InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 1.
Definition: adjointRASModel.C:342
Foam::incompressibleAdjoint::adjointRASModel::diffusionCoeffVar1
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Definition: adjointRASModel.C:400
Foam::incompressibleAdjoint::adjointTurbulenceModel::adjointVars_
incompressibleAdjointMeanFlowVars & adjointVars_
Definition: adjointTurbulenceModel.H:91
Foam::incompressibleAdjoint::adjointRASModel::coeffDict
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: adjointRASModel.H:219
Foam::incompressibleAdjoint::adjointTurbulenceModel::primalVars_
incompressibleVars & primalVars_
Definition: adjointTurbulenceModel.H:90
Foam::incompressibleAdjoint::adjointRASModel::adjointMomentumBCSource
virtual const boundaryVectorField & adjointMomentumBCSource() const =0
Foam::incompressibleAdjoint::adjointRASModel::~adjointRASModel
virtual ~adjointRASModel()=default
Destructor.
Foam::nearWallDist
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:53
fvc.H
Foam::incompressibleAdjoint::adjointRASModel::resetMeanFields
void resetMeanFields()
Reset mean fields to zero.
Definition: adjointRASModel.C:418
Foam::incompressibleAdjoint::adjointRASModel::includeDistance_
bool includeDistance_
Definition: adjointRASModel.H:136
Foam::incompressibleAdjoint::adjointRASModel::setMeanFields
void setMeanFields()
Set mean fields.
Definition: adjointRASModel.C:64
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::incompressibleAdjoint::adjointRASModel::printCoeffs_
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Definition: adjointRASModel.H:104
Foam::incompressibleVars
Base class for solution control classes.
Definition: incompressibleVars.H:54
Foam::incompressibleAdjoint::adjointRASModel::wallShapeSensitivities
virtual const boundaryVectorField & wallShapeSensitivities()=0
Sensitivity terms for shape optimisation, emerging from.
solverControl.H
Foam::incompressibleAdjoint::adjointRASModel::computeMeanFields
void computeMeanFields()
Average adjoint fields on the fly.
Definition: adjointRASModel.C:437
Foam::incompressibleAdjoint::adjointRASModel::devReff
virtual tmp< volSymmTensorField > devReff() const =0
Return the effective stress tensor including the laminar stress.
Foam::incompressibleAdjoint::adjointRASModel::divDevReff
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const =0
Return the diffusion term for the momentum equation.
autoPtr.H
Foam::incompressibleAdjoint::adjointRASModel::New
static autoPtr< adjointRASModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=adjointTurbulenceModel::typeName)
Return a reference to the selected adjointRAS model.
Definition: adjointRASModel.C:163