Loading...
Searching...
No Matches
RASModelVariables.C
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-2023 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\*---------------------------------------------------------------------------*/
29
30#include "RASModelVariables.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
37namespace incompressible
38{
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
45
46// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47
49{
50 if (solverControl_.storeInitValues())
51 {
52 Info<< "Storing initial values of turbulence variables" << endl;
53
54 if (hasTMVar1())
55 {
56 TMVar1InitPtr_.reset
57 (
58 new volScalarField(TMVar1Inst().name()+"Init", TMVar1Inst())
59 );
60 }
61
62 if (hasTMVar2())
63 {
64 TMVar2InitPtr_.reset
65 (
66 new volScalarField(TMVar2Inst().name()+"Init", TMVar2Inst())
67 );
68 }
69
70 if (hasNut())
71 {
72 nutInitPtr_.reset
73 (
75 );
76 }
77 }
78}
79
80
82{
84 {
85 Info<< "Allocating mean values of turbulence variables" << endl;
86
87 if (hasTMVar1())
88 {
89 TMVar1MeanPtr_.reset
90 (
92 (
93 IOobject
94 (
95 TMVar1Inst().name()+"Mean",
97 mesh_,
100 ),
101 TMVar1Inst()
102 )
103 );
104 }
105
106 if (hasTMVar2())
107 {
108 TMVar2MeanPtr_.reset
109 (
111 (
112 IOobject
113 (
114 TMVar2Inst().name()+"Mean",
115 mesh_.time().timeName(),
116 mesh_,
119 ),
120 TMVar2Inst()
121 )
122 );
123 }
124
125 if (hasNut())
126 {
127 nutMeanPtr_.reset
128 (
130 (
131 IOobject
132 (
133 nutRefInst().name()+"Mean",
134 mesh_.time().timeName(),
135 mesh_,
138 ),
139 nutRefInst()
140 )
141 );
142 }
143 }
144}
145
146
148RASModelVariables::cloneRefPtr(const refPtr<volScalarField>& obj) const
149{
150 if (obj)
151 {
152 const volScalarField& sf = obj();
153
154 const word timeName = mesh_.time().timeName();
155
156 return refPtr<volScalarField>::New(sf.name() + timeName, sf);
157 }
158
159 return nullptr;
160}
161
162
164(
165 volScalarField& f1,
167)
168{
169 f1 == f2;
170 const word name1 = f1.name();
171 const word name2 = f2.name();
172
173 // Extra rename to avoid databese collision
174 f2.rename("temp");
175 f1.rename(name2);
176 f2.rename(name1);
177}
178
179
180// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
181
183(
184 const fvMesh& mesh,
185 const solverControl& SolverControl
186)
187:
188 mesh_(mesh),
189 solverControl_(SolverControl),
190
191 TMVar1BaseName_(),
192 TMVar2BaseName_(),
193 nutBaseName_("nut"),
194
195 TMVar1Ptr_(nullptr),
196 TMVar2Ptr_(nullptr),
197 nutPtr_(nullptr),
198 distPtr_(nullptr),
199
200 TMVar1InitPtr_(nullptr),
201 TMVar2InitPtr_(nullptr),
202 nutInitPtr_(nullptr),
204 TMVar1MeanPtr_(nullptr),
205 TMVar2MeanPtr_(nullptr),
206 nutMeanPtr_(nullptr)
207{}
208
209
211(
212 const RASModelVariables& rmv
213)
214:
215 mesh_(rmv.mesh_),
216 solverControl_(rmv.solverControl_),
217
218 TMVar1BaseName_(rmv.TMVar1BaseName_),
219 TMVar2BaseName_(rmv.TMVar2BaseName_),
220 nutBaseName_(rmv.nutBaseName_),
221
222 TMVar1Ptr_(cloneRefPtr(rmv.TMVar1Ptr_)),
223 TMVar2Ptr_(cloneRefPtr(rmv.TMVar2Ptr_)),
224 nutPtr_(cloneRefPtr(rmv.nutPtr_)),
225 distPtr_(cloneRefPtr(rmv.distPtr_)),
226
227 TMVar1InitPtr_(nullptr),
228 TMVar2InitPtr_(nullptr),
229 nutInitPtr_(nullptr),
231 TMVar1MeanPtr_(nullptr),
232 TMVar2MeanPtr_(nullptr),
233 nutMeanPtr_(nullptr)
234{}
235
236
240}
241
242
243// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
244
246(
247 const fvMesh& mesh,
248 const solverControl& SolverControl
249)
250{
251 const IOdictionary modelDict
252 (
254 (
256 mesh.time().constant(),
257 mesh,
261 )
262 );
263
264 word modelType("laminar"); // default to laminar
265
266 const dictionary* dictptr = modelDict.findDict("RAS");
267
268 if (dictptr)
269 {
270 // "RASModel" for v2006 and earlier
271 dictptr->readCompat("model", {{"RASModel", -2006}}, modelType);
272 }
273 else
274 {
275 dictptr = &dictionary::null;
276 }
277
278 Info<< "Creating references for RASModel variables : " << modelType << endl;
279
280 auto* ctorPtr = dictionaryConstructorTable(modelType);
281
282 if (!ctorPtr)
283 {
285 (
286 *dictptr,
287 "RASModelVariables",
288 modelType,
289 *dictionaryConstructorTablePtr_
290 ) << exit(FatalIOError);
291 }
293 return autoPtr<RASModelVariables>(ctorPtr(mesh, SolverControl));
294}
295
296
297// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
298
300(
302) const
303{
305 << "jutJacobianVar1 not implemented for the current turbulence model."
306 << "Returning zero field" << endl;
307
309 (
310 "nutJacobianVar1",
312 mesh_,
314 );
315}
316
317
319(
321) const
322{
324 << "nutJacobianVar2 not implemented for the current turbulence model."
325 << "Returning zero field" << endl;
326
328 (
329 "nutJacobianVar2",
331 mesh_,
333 );
334}
335
336
338{
340 {
341 if (hasTMVar1())
342 {
344 }
345 if (hasTMVar2())
346 {
348 }
349 if (hasNut())
351 nutRefInst() == nutInitPtr_();
352 }
353 }
354}
355
356
358{
360 {
361 Info<< "Resetting mean turbulent fields to zero" << endl;
362
363 // Reset fields to zero
364 if (TMVar1MeanPtr_)
365 {
366 TMVar1MeanPtr_.ref() == Zero;
367 }
368 if (TMVar2MeanPtr_)
369 {
370 TMVar2MeanPtr_.ref() == Zero;
371 }
372 if (nutMeanPtr_)
374 nutMeanPtr_.ref() == Zero;
375 }
376 }
377}
378
379
381{
383 {
384 const label iAverageIter = solverControl_.averageIter();
385 const scalar avIter(iAverageIter);
386 const scalar oneOverItP1 = 1./(avIter + 1);
387 const scalar mult = avIter*oneOverItP1;
388
389 if (hasTMVar1())
390 {
391 TMVar1MeanPtr_.ref() ==
392 (TMVar1MeanPtr_()*mult + TMVar1Inst()*oneOverItP1);
393 }
394 if (hasTMVar2())
395 {
396 TMVar2MeanPtr_.ref() ==
397 (TMVar2MeanPtr_()*mult + TMVar2Inst()*oneOverItP1);
398 }
399 if (hasNut())
400 {
401 nutMeanPtr_.ref() ==
402 (nutMeanPtr_()*mult + nutRefInst()*oneOverItP1);
403 }
404 }
405}
406
407
408tmp<volSymmTensorField> RASModelVariables::devReff
409(
410 const singlePhaseTransportModel& laminarTransport,
411 const volVectorField& U
412) const
413{
415 (
416 "devRhoReff",
419 );
420}
421
422
424(
426)
427{
428 if (hasTMVar1())
429 {
432 {
433 TMVar1MeanPtr_.ref().correctBoundaryConditions();
434 }
435 }
436
437 if (hasTMVar2())
438 {
441 {
442 TMVar2MeanPtr_.ref().correctBoundaryConditions();
443 }
444 }
445
446 if (hasNut())
447 {
449 if (solverControl_.average())
451 nutMeanPtr_.ref().correctBoundaryConditions();
452 }
453 }
454}
455
456
457void RASModelVariables::transfer(RASModelVariables& rmv)
458{
459 if (rmv.hasTMVar1() && hasTMVar1())
460 {
461 copyAndRename(TMVar1Inst(), rmv.TMVar1Inst());
462 }
463
464 if (rmv.hasTMVar2() && hasTMVar2())
465 {
466 copyAndRename(TMVar2Inst(), rmv.TMVar2Inst());
467 }
468
469 if (rmv.hasNut() && hasNut())
470 {
471 copyAndRename(nutRefInst(), rmv.nutRefInst());
472 }
473
474 if (rmv.hasDist() && hasDist())
475 {
476 copyAndRename(d(), rmv.d());
477 }
478}
479
480
481// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
482
483} // End namespace incompressible
484} // End namespace Foam
485
486// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
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())
void correctBoundaryConditions()
Correct boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
bool readCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val using any compatibility names if needed. FatalIOError if there are exc...
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition dictionary.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
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 & TMVar1Inst() const
Return references to instantaneous turbulence fields.
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
correct bounday conditions of turbulent fields
autoPtr< RASModelVariables > clone() const
Clone.
void restoreInitValues()
Restore turbulent fields to their initial values.
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
RASModelVariables(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
void copyAndRename(volScalarField &f1, volScalarField &f2)
const volScalarField & nutRefInst() const
virtual void computeMeanFields()
Compute mean fields on the fly.
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.
virtual bool hasTMVar1() const
Bools to identify which turbulent fields are present.
const volScalarField & TMVar2Inst() const
A class for managing references or pointers (no reference counting).
Definition refPtr.H:54
static refPtr< T > New(Args &&... args)
Construct refPtr with forwarding arguments.
Definition refPtr.H:187
virtual void rename(const word &newName)
Rename.
A simple single-phase transport model based on viscosityModel.
Base class for solver control classes.
bool doAverageIter() const
Whether or not to add fields of the current iteration to the average fields.
bool storeInitValues() const
Re-initialize.
label & averageIter()
Return average iteration index reference.
bool average() const
Whether averaging is enabled or not.
A class for managing temporary objects.
Definition tmp.H:75
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition error.H:637
auto & name
word timeName
Definition getTimeIndex.H:3
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
singlePhaseTransportModel laminarTransport(U, phi)