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-2019 PCOpt/NTUA
9 Copyright (C) 2013-2019 FOSS GP
10 Copyright (C) 2019-2021 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
36{
37namespace incompressible
38{
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
44
45
46// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47
49{
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 (
74 new volScalarField(nutRefInst().name()+"Init", nutRefInst())
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 (
94 (
95 TMVar1Inst().name()+"Mean",
97 mesh_,
100 ),
101 TMVar1Inst()
102 )
103 );
104 }
105
106 if (hasTMVar2())
107 {
108 TMVar2MeanPtr_.reset
109 (
111 (
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 (
132 (
133 nutRefInst().name()+"Mean",
134 mesh_.time().timeName(),
135 mesh_,
138 ),
139 nutRefInst()
140 )
141 );
142 }
143 }
144}
145
146
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),
203
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),
230
231 TMVar1MeanPtr_(nullptr),
232 TMVar2MeanPtr_(nullptr),
233 nutMeanPtr_(nullptr)
234{}
235
236
238{
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,
260 false // Do not register
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 }
292
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 (
311 (
312 "nutJacobianVar1",
313 mesh_.time().timeName(),
314 mesh_,
317 ),
318 mesh_,
320 );
321}
322
323
325(
327) const
328{
330 << "nutJacobianVar2 not implemented for the current turbulence model."
331 << "Returning zero field" << endl;
332
334 (
336 (
337 "nutJacobianVar2",
338 mesh_.time().timeName(),
339 mesh_,
342 ),
343 mesh_,
345 );
346}
347
348
350{
352 {
353 if (hasTMVar1())
354 {
356 }
357 if (hasTMVar2())
358 {
360 }
361 if (hasNut())
362 {
363 nutRefInst() == nutInitPtr_();
364 }
365 }
366}
367
368
370{
372 {
373 Info<< "Resetting mean turbulent fields to zero" << endl;
374
375 // Reset fields to zero
376 if (hasTMVar1())
377 {
378 TMVar1MeanPtr_.ref() ==
379 dimensionedScalar(TMVar1Inst().dimensions(), Zero);
380 }
381 if (hasTMVar2())
382 {
383 TMVar2MeanPtr_.ref() ==
384 dimensionedScalar(TMVar2Inst().dimensions(), Zero);
385 }
386 if (hasNut())
387 {
388 nutMeanPtr_.ref() ==
389 dimensionedScalar(nutRefInst().dimensions(), Zero);
390 }
391 }
392}
393
394
396{
398 {
399 const label iAverageIter = solverControl_.averageIter();
400 const scalar avIter(iAverageIter);
401 const scalar oneOverItP1 = 1./(avIter + 1);
402 const scalar mult = avIter*oneOverItP1;
403
404 if (hasTMVar1())
405 {
406 TMVar1MeanPtr_.ref() ==
407 (TMVar1MeanPtr_()*mult + TMVar1Inst()*oneOverItP1);
408 }
409 if (hasTMVar2())
410 {
411 TMVar2MeanPtr_.ref() ==
412 (TMVar2MeanPtr_()*mult + TMVar2Inst()*oneOverItP1);
413 }
414 if (hasNut())
415 {
416 nutMeanPtr_.ref() ==
417 (nutMeanPtr_()*mult + nutRefInst()*oneOverItP1);
418 }
419 }
420}
421
422
424(
426 const volVectorField& U
427) const
428{
430 (
432 (
433 "devRhoReff",
434 mesh_.time().timeName(),
435 mesh_,
438 ),
440 );
441}
442
443
445(
446 const incompressible::turbulenceModel& turbulence
447)
448{
449 if (hasTMVar1())
450 {
453 {
454 TMVar1MeanPtr_.ref().correctBoundaryConditions();
455 }
456 }
457
458 if (hasTMVar2())
459 {
462 {
463 TMVar2MeanPtr_.ref().correctBoundaryConditions();
464 }
465 }
466
467 if (hasNut())
468 {
471 {
472 nutMeanPtr_.ref().correctBoundaryConditions();
473 }
474 }
475}
476
477
479{
480 if (rmv.hasTMVar1() && hasTMVar1())
481 {
483 }
484
485 if (rmv.hasTMVar2() && hasTMVar2())
486 {
488 }
489
490 if (rmv.hasNut() && hasNut())
491 {
493 }
494
495 if (rmv.hasDist() && hasDist())
496 {
497 copyAndRename(d(), rmv.d());
498 }
499}
500
501
502// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
503
504} // End namespace incompressible
505} // End namespace Foam
506
507// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
void correctBoundaryConditions()
Correct boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:180
Templated abstract base class for single-phase incompressible turbulence models.
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
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
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
bool readCompat(const word &keyword, std::initializer_list< std::pair< const char *, int > > compat, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:394
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
virtual tmp< volSymmTensorField > devReff() const
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
const volScalarField & TMVar1Inst() const
return references to instantaneous turbulence fields
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.
void resetMeanFields()
Reset mean fields to zero.
virtual bool hasTMVar1() const
Bools to identify which turbulent fields are present.
const volScalarField & nutRef() const
const volScalarField & TMVar2Inst() const
static const word propertiesName
Default name of the phase properties dictionary.
Definition: phaseSystem.H:290
A class for managing references or pointers (no reference counting)
Definition: refPtr.H:58
virtual void rename(const word &newName)
Rename.
Definition: regIOobject.C:417
A simple single-phase transport model based on viscosityModel.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Base class for solver control classes.
Definition: solverControl.H:52
bool doAverageIter() const
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:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
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:54
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
singlePhaseTransportModel laminarTransport(U, phi)