Loading...
Searching...
No Matches
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-2023 PCOpt/NTUA
9 Copyright (C) 2013-2023 FOSS GP
10 Copyright (C) 2019 OpenCFD Ltd.
11-------------------------------------------------------------------------------
12License
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
29Namespace
30 Foam::incompressible::adjointRASModels
31
32Description
33 Namespace for incompressible adjointRAS turbulence models.
34
35Class
36 Foam::incompressibleAdjoint::adjointRASModel
37
38Description
39 Abstract base class for incompressible turbulence models.
40
41SourceFiles
42 adjointRASModel.C
43
44\*---------------------------------------------------------------------------*/
45
46#ifndef adjointRASModel_H
47#define adjointRASModel_H
48
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"
62#include "objectiveManager.H"
63#include "boundaryFieldsFwd.H"
64#include "createZeroField.H"
65#include "solverControl.H"
66
67// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69namespace Foam
70{
72{
73
74/*---------------------------------------------------------------------------*\
75 Class adjointRASModel Declaration
76\*---------------------------------------------------------------------------*/
77
78class adjointRASModel
79:
80 public adjointTurbulenceModel,
81 public IOdictionary
82{
83private:
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
94protected:
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
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 //- Base names of the adjoint fields
122 //- Adjoint turbulence model variable 1, mean value
124
125 //- Adjoint turbulence model variable 2, mean value
127
128 //- Source to the adjoint momentum BC emerging
129 //- from differentiating the turbulence model
132 //- Wall sensitivity term for shape optimisation
134
135 //- Wall sensitivity term for flow control optimisation
137
138 //- Does the turbulence model include distances and should the
139 //- adjoint to the distance field be computed
140 bool includeDistance_;
142 //- Has the primal solution changed?
144
145
146 // Protected Member Functions
147
148 //- Print model coefficients
149 virtual void printCoeffs();
150
151 //- Set mean fields
153
154
155public:
156
157 //- Runtime type information
158 TypeName("adjointRASModel");
159
160
161 // Declare run-time constructor selection table
164 (
165 autoPtr,
166 adjointRASModel,
169 incompressibleVars& primalVars,
171 objectiveManager& objManager,
172 const word& adjointTurbulenceModelName
173 ),
174 (
175 primalVars,
176 adjointVars,
177 objManager,
178 adjointTurbulenceModelName
179 )
180 );
181
182
183 // Constructors
184
185 //- Construct from components
186 adjointRASModel
187 (
188 const word& type,
189 incompressibleVars& primalVars,
191 objectiveManager& objManager,
192 const word& adjointTurbulenceModelName =
193 adjointTurbulenceModel::typeName
194 );
195
196
197 // Selectors
198
199 //- Return a reference to the selected adjointRAS model
201 (
202 incompressibleVars& primalVars,
204 objectiveManager& objManager,
205 const word& adjointTurbulenceModelName =
206 adjointTurbulenceModel::typeName
207 );
208
209
210 //- Destructor
211 virtual ~adjointRASModel() = default;
212
213
214 // Member Functions
215
216 //- Return the near wall distances
217 const nearWallDist& y() const
218 {
219 return y_;
220 }
221
222 //- Const access to the coefficients dictionary
223 const dictionary& coeffDict() const
224 {
225 return coeffDict_;
226 }
227
228 //- Const access to the primal solver name
229 const word& primalSolverName() const
230 {
231 return primalVars_.solverName();
232 }
233
234 //- Const access to the adjoint solver name
235 const word& adjointSolverName() const
236 {
237 return adjointVars_.solverName();
238 }
239
240 //- Return non-constant reference to adjoint turbulence model variable 1
241 // Will allocate and return a zero field in case it does not exist
243
244 //- Return non-constant reference to adjoint turbulence model variable 2
245 // Will allocate and return a zero field in case it does not exist
247
248 //- Return non-constant reference to adjoint turbulence model variable 1
249 // Will return the mean value if present,
250 // otherwise the instantaneous value
252
253 //- Return non-constant reference to adjoint turbulence model variable 2
254 // Will return the mean value if present,
255 // otherwise the instantaneous value
257
258 //- Return non-constant autoPtr to adjoint turbulence model variable 1
260
261 //- Return non-constant autoPtr to adjoint turbulence model variable 2
263
264 //- Return reference to the adjoint turbulence model variables base
265 //- names
267
268 //- Return the effective stress tensor including the laminar stress
269 virtual tmp<volSymmTensorField> devReff() const = 0;
270
271 //- Return the effective stress tensor based on a given velocity field
273 (
274 const volVectorField& U
275 ) const = 0;
276
277 //- Return the diffusion term for the momentum equation
279
280 //- Source terms to the adjoint momentum equation due to
281 //- the differentiation of the turbulence model
283
284 //- Jacobian of nut wrt the first turbulence model variable
285 // Needed for objective functions that depend on nut. Defaults to zero
287
288 //- Jacobian of nut wrt the second turbulence model variable
289 // Needed for objective functions that depend on nut. Defaults to zero
291
292 //- Jacobian of nut wrt the flow velocity
293 // Assumes we want to get contributions of mult*d(nut)/dU
294 // Since the dependency of nut to U is usually through a differential
295 // operator, the term multiplying d(nut)/dU is passed as an argument
296 // to this function; the latter should then compute the
297 // contribution of the afforementioned term to the adjoint mean flow
298 // equations
300 (
301 tmp<volScalarField>& dNutdUMult
302 ) const;
303
304 //- Diffusion coefficient of the first primal and adjoint turbulence
305 //- model equation. Needed for some adjoint BCs. Defaults to zero
306 virtual tmp<scalarField> diffusionCoeffVar1(label patchI) const;
307
308 //- Diffusion coefficient of the second primal and adjoint turbulence
309 //- model equation. Needed for some adjoint BCs. Defaults to zero
310 virtual tmp<scalarField> diffusionCoeffVar2(label patchI) const;
311
312 //- Source for the outlet adjoint momentum BC coming from
313 //- differentiating the turbulence model
314 virtual const boundaryVectorField& adjointMomentumBCSource() const = 0;
315
316 //- Sensitivity terms for shape optimisation, emerging from
317 //- the turbulence model differentiation.
318 // Misses dxdb, to be added by the classes assembling the sensitivities
319 virtual const boundaryVectorField& wallShapeSensitivities() = 0;
320
321 //- Sensitivity terms for flow control, emerging from the
322 //- turbulence model differentiation
323 virtual const boundaryVectorField& wallFloCoSensitivities() = 0;
324
325 //- Sensitivity terms resulting from the differentiation of the
326 //- distance field. Misses dxdb, to be added by the classes
327 //- assembling the sensitivities
329
330 //- Term contributing to the computation of FI-based sensitivities
331 // Misses grad(dxdb), to be added by the assembling the sensitivities
333
334 //- Term contributing to the computation of topology optimisation
335 //- sensitivities
336 // Misses betaMax, dBeta/dAlpha and the mesh volume, to be added
337 // during the assembly of the sensitivities
340 const word& designVarsName
341 ) const = 0;
342
343 //- Solve the adjoint turbulence equations
344 virtual void correct();
345
346 //- Read adjointRASProperties dictionary
347 virtual bool read();
348
349 //- Set flag of changed primal solution to true
351
352 //- Restore field values to the initial ones
353 void restoreInitValues();
354
355 //- Reset mean fields to zero
356 void resetMeanFields();
357
358 //- Average adjoint fields on the fly
359 void computeMeanFields();
360
361 //- Should the adjoint to the eikonal equation be computed
362 bool includeDistance() const;
363
364 //- Nullify all adjoint turbulence model fields and their old times
365 virtual void nullify() = 0;
366};
367
368
369// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
370
371} // End namespace incompressible
372} // End namespace Foam
373
374// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
375
376#endif
Bound the given scalar field if it has gone unbounded.
Useful typenames for fields defined only at the boundaries.
IOdictionary(const IOobject &io, const dictionary *fallback=nullptr)
Construct given an IOobject and optional fallback dictionary content.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition Switch.H:81
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
dictionary()
Default construct, a top-level empty dictionary.
Definition dictionary.C:68
Manages the adjoint mean flow fields and their mean values.
virtual tmp< scalarField > diffusionCoeffVar1(label patchI) const
Diffusion coefficient of the first primal and adjoint turbulence model equation. Needed for some adjo...
virtual tmp< volVectorField > adjointMeanFlowSource()=0
Source terms to the adjoint momentum equation due to the differentiation of the turbulence model.
virtual tmp< volScalarField > nutJacobianTMVar2() const
Jacobian of nut wrt the second turbulence model variable.
volScalarField & getAdjointTMVariable1Inst()
Return non-constant reference to adjoint turbulence model variable 1.
virtual const boundaryVectorField & wallShapeSensitivities()=0
Sensitivity terms for shape optimisation, emerging from the turbulence model differentiation.
dictionary coeffDict_
Model coefficients dictionary.
virtual const boundaryVectorField & adjointMomentumBCSource() const =0
Source for the outlet adjoint momentum BC coming from differentiating the turbulence model.
virtual tmp< volScalarField > distanceSensitivities()=0
Sensitivity terms resulting from the differentiation of the distance field. Misses dxdb,...
bool includeDistance_
Does the turbulence model include distances and should the adjoint to the distance field be computed.
virtual tmp< volTensorField > FISensitivityTerm()=0
Term contributing to the computation of FI-based sensitivities.
autoPtr< boundaryVectorField > wallShapeSensitivitiesPtr_
Wall sensitivity term for shape optimisation.
autoPtr< volScalarField > adjointTMVariable1MeanPtr_
Adjoint turbulence model variable 1, mean value.
autoPtr< volScalarField > adjointTMVariable2MeanPtr_
Adjoint turbulence model variable 2, mean value.
volScalarField & getAdjointTMVariable2Inst()
Return non-constant reference to adjoint turbulence model variable 2.
virtual void correct()
Solve the adjoint turbulence equations.
virtual tmp< volSymmTensorField > devReff() const =0
Return the effective stress tensor including the laminar stress.
const word & adjointSolverName() const
Const access to the adjoint solver name.
nearWallDist y_
Near wall distance boundary field.
Switch adjointTurbulence_
Turbulence on/off flag.
virtual tmp< volSymmTensorField > devReff(const volVectorField &U) const =0
Return the effective stress tensor based on a given velocity field.
void restoreInitValues()
Restore field values to the initial ones.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const =0
Return the diffusion term for the momentum equation.
virtual tmp< scalarField > topologySensitivities(const word &designVarsName) const =0
Term contributing to the computation of topology optimisation sensitivities.
bool changedPrimalSolution_
Has the primal solution changed?
autoPtr< boundaryVectorField > adjMomentumBCSourcePtr_
Source to the adjoint momentum BC emerging from differentiating the turbulence model.
virtual tmp< volScalarField > nutJacobianTMVar1() const
Jacobian of nut wrt the first turbulence model variable.
const nearWallDist & y() const
Return the near wall distances.
virtual void printCoeffs()
Print model coefficients.
bool includeDistance() const
Should the adjoint to the eikonal equation be computed.
autoPtr< volScalarField > & getAdjointTMVariable2InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 2.
autoPtr< volScalarField > adjointTMVariable1Ptr_
Adjoint turbulence model variable 1.
autoPtr< volScalarField > adjointTMVariable2Ptr_
Adjoint turbulence model variable 2.
void computeMeanFields()
Average adjoint fields on the fly.
declareRunTimeSelectionTable(autoPtr, adjointRASModel, dictionary,(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName),(primalVars, adjointVars, objManager, adjointTurbulenceModelName))
void resetMeanFields()
Reset mean fields to zero.
void setChangedPrimalSolution()
Set flag of changed primal solution to true.
objectiveManager & objectiveManager_
Reference to the objectiveManager.
volScalarField & getAdjointTMVariable2()
Return non-constant reference to adjoint turbulence model variable 2.
volScalarField & getAdjointTMVariable1()
Return non-constant reference to adjoint turbulence model variable 1.
TypeName("adjointRASModel")
Runtime type information.
const wordList & getAdjointTMVariablesBaseNames() const
Return reference to the adjoint turbulence model variables base names.
wordList adjointTMVariablesBaseNames_
Base names of the adjoint fields.
virtual tmp< volVectorField > nutJacobianU(tmp< volScalarField > &dNutdUMult) const
Jacobian of nut wrt the flow velocity.
virtual tmp< scalarField > diffusionCoeffVar2(label patchI) const
Diffusion coefficient of the second primal and adjoint turbulence model equation. Needed for some adj...
autoPtr< boundaryVectorField > wallFloCoSensitivitiesPtr_
Wall sensitivity term for flow control optimisation.
autoPtr< volScalarField > & getAdjointTMVariable1InstPtr()
Return non-constant autoPtr to adjoint turbulence model variable 1.
virtual void nullify()=0
Nullify all adjoint turbulence model fields and their old times.
virtual ~adjointRASModel()=default
Destructor.
virtual const boundaryVectorField & wallFloCoSensitivities()=0
Sensitivity terms for flow control, emerging from the turbulence model differentiation.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
virtual bool read()
Read adjointRASProperties dictionary.
static autoPtr< adjointRASModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=adjointTurbulenceModel::typeName)
Return a reference to the selected adjointRAS model.
const word & primalSolverName() const
Const access to the primal solver name.
Base class for solution control classes.
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Class for managing objective functions.
A class for managing temporary objects.
Definition tmp.H:75
const word & solverName() const
Return solver name.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for incompressible adjoint turbulence models.
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
GeometricField< vector, fvPatchField, volMesh > volVectorField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
volVectorField::Boundary boundaryVectorField
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
Foam::surfaceFields.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68