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-2021 PCOpt/NTUA
9 Copyright (C) 2013-2021 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
65namespace incompressibleAdjoint
66{
67
68/*---------------------------------------------------------------------------*\
69 Class adjointTurbulenceModel Declaration
70\*---------------------------------------------------------------------------*/
73:
74 public regIOobject
75{
76private:
77
78 // Private Member Functions
79
80 //- No copy construct
82
83 //- No copy assignment
84 void operator=(const adjointTurbulenceModel&) = delete;
85
86
87protected:
88
89 // Protected data
93 const Time& runTime_;
94 const fvMesh& mesh_;
95
96
97public:
98
99 //- Runtime type information
100 TypeName("adjointTurbulenceModel");
101
102
103 // Declare run-time New selection table
106 (
107 autoPtr,
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
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 {
159 }
160
161 //- Return the turbulence viscosity
162 virtual const volScalarField& nut() const
163 {
164 return primalVars_.RASModelVariables()().nutRef();
165 }
166
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 =
175 turbVars = primalVars_.RASModelVariables();
176
177 return
179 (
180 "nuEff",
181 lamTrans.nu()() + turbVars().nutRef()
182 );
183 //return primalVars_.turbulence()().nuEff();
184 }
185
186 //- Return the effective viscosity on a given patch
187 virtual tmp<scalarField> nuEff(const label patchI) const
188 {
189 // Go through RASModelVariables::nutRef in order to obtain
190 // the mean field, if present
191 const singlePhaseTransportModel& lamTrans =
194 turbVars = primalVars_.RASModelVariables();
195
196 return
197 (
198 lamTrans.nu()().boundaryField()[patchI]
199 + turbVars().nutRef().boundaryField()[patchI]
200 );
201 }
202
203 //- Return the effective stress tensor including the laminar stress
204 virtual tmp<volSymmTensorField> devReff() const = 0;
205
206 //- Return the effective stress tensor based on a given velocity field
208 (
209 const volVectorField& U
210 ) const = 0;
211
212 //- Return the diffusion term for the momentum equation
214
215 //- Source term added to the adjoint mean flow due to the
216 // differentiation of the turbulence model
218
219 //- Solve the adjoint turbulence equations
220 virtual void correct() = 0;
221
222 //- Read adjointLESProperties or adjointRASProperties dictionary
223 virtual bool read() = 0;
224
225 //- Default dummy write function
226 virtual bool writeData(Ostream&) const
227 {
228 return true;
229 }
230
231 //- Nullify all adjoint turbulence model fields and their old times
232 virtual void nullify() = 0;
233};
234
235
236// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237
238} // End namespace incompressibleAdjoint
239} // End namespace Foam
240
241// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
242
243#endif
244
245// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Manages the adjoint mean flow fields and their mean values.
Abstract base class for incompressible adjoint turbulence models (RAS, LES and laminar).
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 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 incompressible objective functions.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:76
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:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
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:73