incompressibleVars.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
28Class
29 Foam::incompressibleVars
30
31Description
32 Base class for solution control classes
33
34\*---------------------------------------------------------------------------*/
35
36#ifndef incompressibleVars_H
37#define incompressibleVars_H
38
39#include "variablesSet.H"
40#include "fvMesh.H"
43#include "RASModelVariables.H"
44#include "solverControl.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51/*---------------------------------------------------------------------------*\
52 Class incompressibleVars Declaration
53\*---------------------------------------------------------------------------*/
56:
57 public variablesSet
58{
59protected:
60
61 // Protected data
62
63 //- Reference to the solverControl of the solver allocating the fields
65
66 //- Fields involved in the solution of the incompressible NS equations
73 //autoPtr<volScalarField> TPtr_;
74
75 //- Keep a copy of the initial field values if necessary
79
80 //- Manage mean fields. Turbulent mean fields are managed
81 //- in RASModelVariables
85
86 //- Update boundary conditions upon construction.
87 // Useful for cases in which information has been lost on boundary
88 // (e.g. averaging, redistribution)
90
91
92 // Protected Member Functions
93
94 //- Read fields and set turbulence
95 void setFields();
96
97 //- Set initial fields if necessary
98 void setInitFields();
99
100 //- Set mean fields if necessary
101 void setMeanFields();
102
103 //- Rename turbulence fields if necessary
105
106 //- Update boundary conditions of mean-flow
108
109 //- Update boundary conditions of turbulent fields
111
112 //- No copy assignment
113 void operator=(const incompressibleVars&);
114
115
116public:
117
118 // Static Data Members
119
120 //- Run-time type information
121 TypeName("incompressibleVars");
122
123
124 // Constructors
125
126 //- Construct from mesh
128 (
129 fvMesh& mesh,
130 solverControl& SolverControl
131 );
132
133 //- Copy constructor
134 // turbulence_ and laminarTransport_ are uninitialized since there is
135 // no clear way (and possibly no usage) of cloning them.
136 // The resulting incompressibleVars is hence incomplete.
137 // Additionally, copied fields have the original name, appended with
138 // the timeName, to avoid database collisions.
139 // To be used as storing space during the solutiuon of the unsteady
140 // adjoint equations
142
143 //- Clone the incompressibleVars
144 virtual autoPtr<variablesSet> clone() const;
145
146
147
148 //- Destructor
149 virtual ~incompressibleVars() = default;
150
151
152 // Member Functions
153
154 // Access
155
156 // Access to fields that will be later used for the solution of
157 // the adjoint equations. Might be averaged or not, depending on
158 // the corresponding switch
159 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
160
161 //- Return const reference to pressure
162 const volScalarField& p() const;
163
164 //- Return reference to pressure
165 volScalarField& p();
166
167 //- Return const reference to velocity
168 const volVectorField& U() const;
169
170 //- Return reference to velocity
171 volVectorField& U();
172
173 //- Return const reference to volume flux
174 const surfaceScalarField& phi() const;
175
176 //- Return reference to volume flux
178
179
180 // Access to instantaneous fields. Solvers should use these fields
181 // to execute a solver iteration
182 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
183
184 //- Return const reference to pressure
185 inline const volScalarField& pInst() const;
186
187 //- Return reference to pressure
188 inline volScalarField& pInst();
189
190 //- Return const reference to velocity
191 inline const volVectorField& UInst() const;
192
193 //- Return reference to velocity
194 inline volVectorField& UInst();
195
196 //- Return const reference to volume flux
197 inline const surfaceScalarField& phiInst() const;
198
199 //- Return reference to volume flux
200 inline surfaceScalarField& phiInst();
201
202
203 //- Return const reference to transport model
204 inline const singlePhaseTransportModel& laminarTransport() const;
205
206 //- Return reference to transport model
208
209 //- Return const reference to the turbulence model
211 turbulence() const;
212
213 //- Return reference to the turbulence model
215
216 //- Return const reference to the turbulence model variables
218 RASModelVariables() const;
219
220 //- Return reference to the turbulence model variables
223
224 //- Restore field values to the initial ones
225 void restoreInitValues();
226
227 //- Reset mean fields to zero
228 void resetMeanFields();
229
230 //- Compute mean fields on the fly
231 void computeMeanFields();
232
233 //- correct boundaryconditions for all volFields
235
236 //- Return storeInitValues bool
237 bool storeInitValues() const;
238
239 //- Return computeMeanFields bool
240 bool computeMeanFields() const;
241
242 /*
243 //- Return const reference to the temperature
244 const volScalarField& T() const;
245
246 //- Return reference to the temperature
247 volScalarField& T();
248 */
249
250 //- Transfer the fields of another variablesSet to this
251 virtual void transfer(variablesSet& vars);
252
253 //- Write dummy turbulent fields to allow for continuation in
254 //- multi-point, turbulent runs
255 bool write() const;
256};
257
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261} // End namespace Foam
262
263// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264
265#include "incompressibleVarsI.H"
266
267// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268
269#endif
270
271// ************************************************************************* //
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
Base class for solution control classes.
const volVectorField & UInst() const
Return const reference to velocity.
TypeName("incompressibleVars")
Run-time type information.
solverControl & solverControl_
Reference to the solverControl of the solver allocating the fields.
bool storeInitValues() const
Return storeInitValues bool.
const volVectorField & U() const
Return const reference to velocity.
autoPtr< volVectorField > UMeanPtr_
autoPtr< volVectorField > UInitPtr_
const volScalarField & p() const
Return const reference to pressure.
autoPtr< volVectorField > UPtr_
const surfaceScalarField & phi() const
Return const reference to volume flux.
autoPtr< surfaceScalarField > phiPtr_
void restoreInitValues()
Restore field values to the initial ones.
const surfaceScalarField & phiInst() const
Return const reference to volume flux.
autoPtr< singlePhaseTransportModel > laminarTransportPtr_
autoPtr< volScalarField > pPtr_
Fields involved in the solution of the incompressible NS equations.
void correctTurbulentBoundaryConditions()
Update boundary conditions of turbulent fields.
const autoPtr< incompressible::turbulenceModel > & turbulence() const
Return const reference to the turbulence model.
const singlePhaseTransportModel & laminarTransport() const
Return const reference to transport model.
void setFields()
Read fields and set turbulence.
const volScalarField & pInst() const
Return const reference to pressure.
bool correctBoundaryConditions_
Update boundary conditions upon construction.
void computeMeanFields()
Compute mean fields on the fly.
void resetMeanFields()
Reset mean fields to zero.
virtual ~incompressibleVars()=default
Destructor.
autoPtr< surfaceScalarField > phiMeanPtr_
void renameTurbulenceFields()
Rename turbulence fields if necessary.
autoPtr< incompressible::RASModelVariables > RASModelVariables_
void operator=(const incompressibleVars &)
No copy assignment.
void correctNonTurbulentBoundaryConditions()
Update boundary conditions of mean-flow.
virtual void transfer(variablesSet &vars)
Transfer the fields of another variablesSet to this.
void correctBoundaryConditions()
correct boundaryconditions for all volFields
autoPtr< surfaceScalarField > phiInitPtr_
autoPtr< incompressible::turbulenceModel > turbulence_
virtual autoPtr< variablesSet > clone() const
Clone the incompressibleVars.
const autoPtr< incompressible::RASModelVariables > & RASModelVariables() const
Return const reference to the turbulence model variables.
autoPtr< volScalarField > pMeanPtr_
void setInitFields()
Set initial fields if necessary.
autoPtr< volScalarField > pInitPtr_
Keep a copy of the initial field values if necessary.
void setMeanFields()
Set mean fields if necessary.
A simple single-phase transport model based on viscosityModel.
Base class for solver control classes.
Definition: solverControl.H:52
Base class for creating a set of variables.
Definition: variablesSet.H:50
dynamicFvMesh & mesh
Namespace for OpenFOAM.
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73