RASModelVariables.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-2019, 2022 PCOpt/NTUA
9 Copyright (C) 2013-2019, 2022 FOSS GP
10 Copyright (C) 2019-2020 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
28Class
29 Foam::incompressible::RASModelVariables
30
31Description
32 Abstract base class for objective functions. No point in making this
33 runTime selectable since its children will have different constructors.
34
35SourceFiles
36 RASModelVariables.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef RASModelVariables_H
41#define RASModelVariables_H
42
43#include "solverControl.H"
47#include "refPtr.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51namespace Foam
52{
53namespace incompressible
54{
55
56/*---------------------------------------------------------------------------*\
57 Class RASModelVariables Declaration
58\*---------------------------------------------------------------------------*/
61{
62protected:
63
64 // Protected Data
66 const fvMesh& mesh_;
68
69 // Base names of the turbulent fields
78
79 // Initial values (demand-driven)
80 // For finite differences and optimisation runs
84
85 // Mean values (demand-driven)
89
90
91 // Protected Member functions
92
93 virtual void allocateInitValues();
94 virtual void allocateMeanFields();
95
97 cloneRefPtr(const refPtr<volScalarField>& obj) const;
98
100
101 //- No copy assignment
102 void operator=(const RASModelVariables&) = delete;
103
104
105public:
106
107
108 //- Runtime type information
109 TypeName("RASModelVariables");
110
111 // Declare run-time constructor selection table
114 (
115 autoPtr,
118 (
119 const fvMesh& mesh,
120 const solverControl& SolverControl
121 ),
122 (mesh, SolverControl)
123 );
124
125
126 // Constructors
127
128 //- Construct from components
130 (
131 const fvMesh& mesh,
132 const solverControl& SolverControl
133 );
134
135 //- Copy constructor
136 // Will allocate new fields (instead of referencing the ones in the
137 // turbulence model), so cannot be used directly to access the fields
138 // of the turbulence model.
139 // Mainly used for checkpointing in unsteady adjoint
141 (
142 const RASModelVariables& rmv
143 );
144
145 //- Clone
146 // Will allocate new fields (instead of referencing the ones in the
147 // turbulence model), so cannot be used directly to access the fields
148 // of the turbulence model.
149 // Mainly used for checkpointing in unsteady adjoint
151
152
153 // Selectors
154
155 //- Return a reference to the selected turbulence model
157 (
158 const fvMesh& mesh,
159 const solverControl& SolverControl
160 );
161
162
163 // Destructor
164
165 // Destructor does nothing on base since depending on the case new
166 // fields might be allocated
167 // MUST be overloaded in inherited classes
168 virtual ~RASModelVariables() = default;
169
170
171 // Member Functions
172
173 //- Turbulence field names
174 inline const word& TMVar1BaseName() const;
175 inline const word& TMVar2BaseName() const;
176 inline const word& nutBaseName() const;
177
178 //- Bools to identify which turbulent fields are present
179 // Apart from the distance pointer, all other pointers are
180 // allocated even if the the corresponding field does not exist.
181 // Hence, the pointer itself cannot be used to determine the
182 // existance of the field
183 inline virtual bool hasTMVar1() const;
184 inline virtual bool hasTMVar2() const;
185 inline virtual bool hasNut() const;
186 inline bool hasDist() const;
187
188 //- Return references to turbulence fields
189 // will return the mean field if it exists, otherwise the
190 // instantaneous one
191 inline const volScalarField& TMVar1() const;
192 inline volScalarField& TMVar1();
193 inline const volScalarField& TMVar2() const;
194 inline volScalarField& TMVar2();
195 inline const volScalarField& nutRef() const;
196 inline volScalarField& nutRef();
197 inline const volScalarField& d() const;
198 inline volScalarField& d();
199
200 //- return references to instantaneous turbulence fields
201 inline const volScalarField& TMVar1Inst() const;
202 inline volScalarField& TMVar1Inst();
203 inline const volScalarField& TMVar2Inst() const;
204 inline volScalarField& TMVar2Inst();
205 inline const volScalarField& nutRefInst() const;
206 inline volScalarField& nutRefInst();
207
208
209 //- Return nut Jacobian wrt the TM vars
211 (
213 ) const;
214
216 (
218 ) const;
219
220 //- Return the turbulence production term
222 {
224 return nullptr;
225 }
226
227 //- Restore turbulent fields to their initial values
228 void restoreInitValues();
229
230 //- Reset mean fields to zero
231 void resetMeanFields();
232
233 //- Compute mean fields on the fly
234 virtual void computeMeanFields();
235
236 //- Return stress tensor based on the mean flow variables
238 (
240 const volVectorField& U
241 ) const;
242
243 //- correct bounday conditions of turbulent fields
244 virtual void correctBoundaryConditions
245 (
246 const incompressible::turbulenceModel& turbulence
247 );
248
249 //- Transfer turbulence fields from an another object
250 // Copies values since the ownership of the original fields is held
251 // by the turbulence model
252 virtual void transfer(RASModelVariables& rmv);
253};
254
255
256// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257
258} // End namespace incompressible
259} // End namespace Foam
260
261// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262
263#include "RASModelVariablesI.H"
264
265// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266
267#endif
268
269// ************************************************************************* //
Templated abstract base class for single-phase incompressible turbulence models.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
virtual void transfer(RASModelVariables &rmv)
Transfer turbulence fields from an another object.
static autoPtr< RASModelVariables > New(const fvMesh &mesh, const solverControl &SolverControl)
Return a reference to the selected turbulence model.
const volScalarField & TMVar2() const
const volScalarField & TMVar1Inst() const
return references to instantaneous turbulence fields
declareRunTimeSelectionTable(autoPtr, RASModelVariables, dictionary,(const fvMesh &mesh, const solverControl &SolverControl),(mesh, SolverControl))
autoPtr< RASModelVariables > clone() const
Clone.
void restoreInitValues()
Restore turbulent fields to their initial values.
const volScalarField & d() const
refPtr< volScalarField > cloneRefPtr(const refPtr< volScalarField > &obj) const
virtual tmp< volScalarField > nutJacobianVar1(const singlePhaseTransportModel &laminarTransport) const
Return nut Jacobian wrt the TM vars.
virtual tmp< volScalarField > nutJacobianVar2(const singlePhaseTransportModel &laminarTransport) const
void copyAndRename(volScalarField &f1, volScalarField &f2)
const volScalarField & nutRefInst() const
virtual void computeMeanFields()
Compute mean fields on the fly.
const word & TMVar1BaseName() const
Turbulence field names.
tmp< volSymmTensorField > devReff(const singlePhaseTransportModel &laminarTransport, const volVectorField &U) const
Return stress tensor based on the mean flow variables.
void resetMeanFields()
Reset mean fields to zero.
const volScalarField & TMVar1() const
Return references to turbulence fields.
TypeName("RASModelVariables")
Runtime type information.
virtual bool hasTMVar1() const
Bools to identify which turbulent fields are present.
virtual tmp< volScalarField::Internal > G()
Return the turbulence production term.
const volScalarField & nutRef() const
const volScalarField & TMVar2Inst() const
void operator=(const RASModelVariables &)=delete
No copy assignment.
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
A simple single-phase transport model based on viscosityModel.
Base class for solver control classes.
Definition: solverControl.H:52
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
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Namespace for OpenFOAM.
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)
cellMask correctBoundaryConditions()
singlePhaseTransportModel laminarTransport(U, phi)
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73