Loading...
Searching...
No Matches
adjointTurbulenceModel.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::incompressibleAdjoint
31
32Description
33 Namespace for incompressible adjoint turbulence models.
34
35Class
36 Foam::incompressibleAdjoint::adjointTurbulenceModel
37
38Description
39 Abstract base class for incompressible adjoint turbulence models
40 (RAS, LES and laminar).
41
42SourceFiles
43 adjointTurbulenceModel.C
44 newTurbulenceModel.C
45
46\*---------------------------------------------------------------------------*/
47
48#ifndef adjointTurbulenceModel_H
49#define adjointTurbulenceModel_H
50
51#include "incompressibleVars.H"
53#include "objectiveManager.H"
54#include "Time.H"
56
57// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58
59namespace Foam
60{
61
62// Forward declarations
63class fvMesh;
64
66{
67
68/*---------------------------------------------------------------------------*\
69 Class adjointTurbulenceModel Declaration
70\*---------------------------------------------------------------------------*/
71
72class adjointTurbulenceModel
73:
74 public regIOobject
75{
76private:
77
78 // Private Member Functions
79
80 //- No copy construct
81 adjointTurbulenceModel(const adjointTurbulenceModel&) = delete;
82
83 //- No copy assignment
84 void operator=(const adjointTurbulenceModel&) = delete;
85
86
87protected:
88
89 // Protected data
90
93 const Time& runTime_;
94 const fvMesh& mesh_;
95
96
97public:
98
99 //- Runtime type information
100 TypeName("adjointTurbulenceModel");
101
103 // Declare run-time New selection table
104
106 (
108 adjointTurbulenceModel,
109 adjointTurbulenceModel,
110 (
111 incompressibleVars& primalVars,
113 objectiveManager& objManager,
114 const word& adjointTurbulenceModelName
115 ),
116 (
117 primalVars,
118 adjointVars,
119 objManager,
120 adjointTurbulenceModelName
121 )
122 );
123
124
125 // Constructors
126
127 //- Construct from components
128 adjointTurbulenceModel
129 (
130 incompressibleVars& primalVars,
132 objectiveManager& objManager,
133 const word& adjointTurbulenceModelName = typeName
134 );
135
136
137 // Selectors
138
139 //- Return a reference to the selected turbulence model
141 (
142 incompressibleVars& primalVars,
144 objectiveManager& objManager,
145 const word& adjointTurbulenceModelName = typeName
146 );
147
148
149 //- Destructor
150 virtual ~adjointTurbulenceModel() = default;
151
152
153 // Member Functions
154
155 //- Return the laminar viscosity
156 inline tmp<volScalarField> nu() const
157 {
158 return primalVars_.laminarTransport().nu();
159 }
160
161 //- Return the turbulence viscosity
162 virtual const tmp<volScalarField> nut() const
163 {
164 return primalVars_.RASModelVariables()().nut();
165 }
167 //- Return the effective viscosity
168 virtual tmp<volScalarField> nuEff() const
169 {
170 // Go through RASModelVariables::nutRef in order to obtain
171 // the mean field, if present
172 const singlePhaseTransportModel& lamTrans =
173 primalVars_.laminarTransport();
175 turbVars = primalVars_.RASModelVariables();
176
178 (
179 "nuEff",
181 lamTrans.nu() + turbVars().nut()
182 );
183 }
184
185 //- Return the effective viscosity on a given patch
186 virtual tmp<scalarField> nuEff(const label patchI) const
187 {
188 // Go through RASModelVariables::nutRef in order to obtain
189 // the mean field, if present
190 const singlePhaseTransportModel& lamTrans =
193 turbVars = primalVars_.RASModelVariables();
194
195 return (lamTrans.nu(patchI) + turbVars().nut(patchI));
196 }
197
198 //- Return the effective stress tensor including the laminar stress
199 virtual tmp<volSymmTensorField> devReff() const = 0;
200
201 //- Return the effective stress tensor based on a given velocity field
203 (
204 const volVectorField& U
205 ) const = 0;
206
207 //- Return the diffusion term for the momentum equation
209
210 //- Source term added to the adjoint mean flow due to the
211 // differentiation of the turbulence model
213
214 //- Solve the adjoint turbulence equations
215 virtual void correct() = 0;
216
217 //- Read adjointLESProperties or adjointRASProperties dictionary
218 virtual bool read() = 0;
219
220 //- Default dummy write function
221 virtual bool writeData(Ostream&) const
223 return true;
224 }
225
226 //- Nullify all adjoint turbulence model fields and their old times
227 virtual void nullify() = 0;
228};
229
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233} // End namespace incompressibleAdjoint
234} // End namespace Foam
235
236// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238#endif
239
240// ************************************************************************* //
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_REGISTER
Do not request registration (bool: false).
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Manages the adjoint mean flow fields and their mean values.
virtual tmp< volVectorField > adjointMeanFlowSource()=0
Source term added to the adjoint mean flow due to the.
virtual tmp< volScalarField > nuEff() const
Return the effective viscosity.
virtual bool read()=0
Read adjointLESProperties or adjointRASProperties dictionary.
virtual tmp< volSymmTensorField > devReff() const =0
Return the effective stress tensor including the laminar stress.
virtual tmp< volSymmTensorField > devReff(const volVectorField &U) const =0
Return the effective stress tensor based on a given velocity field.
virtual tmp< fvVectorMatrix > divDevReff(volVectorField &U) const =0
Return the diffusion term for the momentum equation.
virtual ~adjointTurbulenceModel()=default
Destructor.
tmp< volScalarField > nu() const
Return the laminar viscosity.
static autoPtr< adjointTurbulenceModel > New(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName=typeName)
Return a reference to the selected turbulence model.
TypeName("adjointTurbulenceModel")
Runtime type information.
declareRunTimeNewSelectionTable(autoPtr, adjointTurbulenceModel, adjointTurbulenceModel,(incompressibleVars &primalVars, incompressibleAdjointMeanFlowVars &adjointVars, objectiveManager &objManager, const word &adjointTurbulenceModelName),(primalVars, adjointVars, objManager, adjointTurbulenceModelName))
virtual void correct()=0
Solve the adjoint turbulence equations.
virtual tmp< scalarField > nuEff(const label patchI) const
Return the effective viscosity on a given patch.
virtual bool writeData(Ostream &) const
Default dummy write function.
virtual void nullify()=0
Nullify all adjoint turbulence model fields and their old times.
virtual const tmp< volScalarField > nut() const
Return the turbulence viscosity.
Base class for solution control classes.
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
Class for managing objective functions.
regIOobject(const IOobject &io, const bool isTimeObject=false)
Construct from IOobject. The optional flag adds special handling if the object is the top-level regIO...
Definition regIOobject.C:43
A simple single-phase transport model based on viscosityModel.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
Namespace for incompressible adjoint turbulence models.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Macros to ease declaration of run-time selection tables.
#define declareRunTimeNewSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection for derived classes.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68