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-2020 OpenCFD Ltd.
11 -------------------------------------------------------------------------------
12 License
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 
35 namespace Foam
36 {
37 namespace incompressible
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(RASModelVariables, 0);
43 defineRunTimeSelectionTable(RASModelVariables, dictionary);
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 {
83  if (solverControl_.average())
84  {
85  Info<< "Allocating mean values of turbulence variables" << endl;
86 
87  if (hasTMVar1())
88  {
89  TMVar1MeanPtr_.reset
90  (
91  new volScalarField
92  (
93  IOobject
94  (
95  TMVar1Inst().name()+"Mean",
96  mesh_.time().timeName(),
97  mesh_,
100  ),
101  TMVar1Inst()
102  )
103  );
104  }
105 
106  if (hasTMVar2())
107  {
108  TMVar2MeanPtr_.reset
109  (
110  new volScalarField
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  (
129  new volScalarField
130  (
131  IOobject
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,
166  volScalarField& f2
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 {
239  return autoPtr<RASModelVariables>::New(*this);
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
244 
246 (
247  const fvMesh& mesh,
248  const solverControl& SolverControl
249 )
250 {
251  const IOdictionary modelDict
252  (
253  IOobject
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 cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
281 
282  if (!cstrIter.found())
283  {
285  (
286  *dictptr,
287  "RASModelVariables",
288  modelType,
289  *dictionaryConstructorTablePtr_
290  ) << exit(FatalIOError);
291  }
292 
293  return autoPtr<RASModelVariables>(cstrIter()(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  IOobject
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  (
335  IOobject
336  (
337  "nutJacobianVar2",
338  mesh_.time().timeName(),
339  mesh_,
342  ),
343  mesh_,
345  );
346 }
347 
348 
350 {
352  {
353  if (hasTMVar1())
354  {
355  TMVar1Inst() == TMVar1InitPtr_();
356  }
357  if (hasTMVar2())
358  {
359  TMVar2Inst() == TMVar2InitPtr_();
360  }
361  if (hasNut())
362  {
363  nutRefInst() == nutInitPtr_();
364  }
365  }
366 }
367 
368 
370 {
371  if (solverControl_.average())
372  {
373  Info<< "Resetting mean turbulent fields to zero" << endl;
374 
375  // Reset fields to zero
376  if (TMVar1Ptr_)
377  {
378  TMVar1MeanPtr_.ref() ==
379  dimensionedScalar(TMVar1Inst().dimensions(), Zero);
380  }
381  if (TMVar2Ptr_)
382  {
383  TMVar2MeanPtr_.ref() ==
384  dimensionedScalar(TMVar2Inst().dimensions(), Zero);
385  }
386  if (nutPtr_)
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  (
431  IOobject
432  (
433  "devRhoReff",
434  mesh_.time().timeName(),
435  mesh_,
438  ),
439  -(laminarTransport.nu() + nutRef())*dev(twoSymm(fvc::grad(U)))
440  );
441 }
442 
443 
445 (
447 )
448 {
449  if (hasTMVar1())
450  {
451  TMVar1Inst().correctBoundaryConditions();
452  if (solverControl_.average())
453  {
454  TMVar1MeanPtr_.ref().correctBoundaryConditions();
455  }
456  }
457 
458  if (hasTMVar2())
459  {
460  TMVar2Inst().correctBoundaryConditions();
461  if (solverControl_.average())
462  {
463  TMVar2MeanPtr_.ref().correctBoundaryConditions();
464  }
465  }
466 
467  if (hasNut())
468  {
469  nutRefInst().correctBoundaryConditions();
470  if (solverControl_.average())
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 // ************************************************************************* //
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:508
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::refPtr::New
static refPtr< T > New(Args &&... args)
Construct refPtr of T with forwarding arguments.
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::solverControl::doAverageIter
bool doAverageIter() const
Definition: solverControlI.H:81
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::incompressible::RASModelVariables::d
const volScalarField & d() const
Definition: RASModelVariablesI.H:144
Foam::incompressible::RASModelVariables::TMVar1Inst
const volScalarField & TMVar1Inst() const
return references to instantaneous turbulence fields
Definition: RASModelVariablesI.H:156
Foam::incompressible::RASModelVariables::solverControl_
const solverControl & solverControl_
Definition: RASModelVariables.H:66
Foam::incompressible::defineTypeNameAndDebug
defineTypeNameAndDebug(adjointEikonalSolver, 0)
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
Foam::incompressible::RASModelVariables
Abstract base class for objective functions. No point in making this runTime selectable since its chi...
Definition: RASModelVariables.H:59
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::incompressible::RASModelVariables::TMVar1InitPtr_
refPtr< volScalarField > TMVar1InitPtr_
Definition: RASModelVariables.H:80
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::incompressible::RASModelVariables::TMVar2Ptr_
refPtr< volScalarField > TMVar2Ptr_
Definition: RASModelVariables.H:74
Foam::incompressible::RASModelVariables::nutInitPtr_
refPtr< volScalarField > nutInitPtr_
Definition: RASModelVariables.H:82
Foam::incompressible::RASModelVariables::restoreInitValues
void restoreInitValues()
Restore turbulent fields to their initial values.
Definition: RASModelVariables.C:349
Foam::FatalIOError
IOerror FatalIOError
Foam::solverControl
Base class for solver control classes.
Definition: solverControl.H:51
Foam::incompressible::RASModelVariables::allocateInitValues
void allocateInitValues()
Definition: RASModelVariables.C:48
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::solverControl::storeInitValues
bool storeInitValues() const
Re-initialize.
Definition: solverControlI.H:51
Foam::incompressible::RASModelVariables::nutMeanPtr_
refPtr< volScalarField > nutMeanPtr_
Definition: RASModelVariables.H:87
Foam::incompressible::RASModelVariables::TMVar2BaseName_
word TMVar2BaseName_
Definition: RASModelVariables.H:70
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:57
Foam::incompressible::RASModelVariables::distPtr_
refPtr< volScalarField > distPtr_
Definition: RASModelVariables.H:76
Foam::incompressible::RASModelVariables::correctBoundaryConditions
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
correct bounday conditions of turbulent fields
Definition: RASModelVariables.C:445
Foam::incompressible::defineRunTimeSelectionTable
defineRunTimeSelectionTable(adjointSensitivity, dictionary)
Foam::incompressible::RASModelVariables::RASModelVariables
RASModelVariables(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
Definition: RASModelVariables.C:183
Foam::incompressible::RASModelVariables::computeMeanFields
void computeMeanFields()
Compute mean fields on the fly.
Definition: RASModelVariables.C:395
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:385
Foam::incompressible::RASModelVariables::transfer
virtual void transfer(RASModelVariables &rmv)
Transfer turbulence fields from an another object.
Definition: RASModelVariables.C:478
Foam::incompressible::RASModelVariables::nutRefInst
const volScalarField & nutRefInst() const
Definition: RASModelVariablesI.H:180
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:406
Foam::incompressible::RASModelVariables::hasDist
bool hasDist() const
Definition: RASModelVariablesI.H:74
Foam::incompressible::RASModelVariables::allocateMeanFields
void allocateMeanFields()
Definition: RASModelVariables.C:81
Foam::incompressible::RASModelVariables::resetMeanFields
void resetMeanFields()
Reset mean fields to zero.
Definition: RASModelVariables.C:369
RASModelVariables.H
Foam::incompressible::RASModelVariables::hasTMVar2
bool hasTMVar2() const
Definition: RASModelVariablesI.H:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::incompressible::RASModelVariables::nutBaseName_
word nutBaseName_
Definition: RASModelVariables.H:71
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::incompressible::RASModelVariables::mesh_
const fvMesh & mesh_
Definition: RASModelVariables.H:65
Foam::incompressible::RASModelVariables::nutPtr_
refPtr< volScalarField > nutPtr_
Definition: RASModelVariables.H:75
Foam::incompressible::RASModelVariables::nutJacobianVar2
virtual tmp< volScalarField > nutJacobianVar2(const singlePhaseTransportModel &laminarTransport) const
Definition: RASModelVariables.C:325
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::incompressible::RASModelVariables::New
static autoPtr< RASModelVariables > New(const fvMesh &mesh, const solverControl &SolverControl)
Return a reference to the selected turbulence model.
Definition: RASModelVariables.C:246
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::incompressible::RASModelVariables::TMVar1Ptr_
refPtr< volScalarField > TMVar1Ptr_
Definition: RASModelVariables.H:73
Foam::incompressible::RASModelVariables::TMVar1BaseName_
word TMVar1BaseName_
Definition: RASModelVariables.H:69
Foam::solverControl::averageIter
label & averageIter()
Return average iteration index reference.
Definition: solverControlI.H:63
Foam::incompressible::RASModelVariables::TMVar2Inst
const volScalarField & TMVar2Inst() const
Definition: RASModelVariablesI.H:168
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::incompressible::RASModelVariables::hasTMVar1
bool hasTMVar1() const
Bools to idenify which turbulent fields are present.
Definition: RASModelVariablesI.H:56
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::incompressible::RASModelVariables::devReff
tmp< volSymmTensorField > devReff(const singlePhaseTransportModel &laminarTransport, const volVectorField &U) const
Return stress tensor based on the mean flow variables.
Definition: RASModelVariables.C:424
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::incompressible::RASModelVariables::TMVar2InitPtr_
refPtr< volScalarField > TMVar2InitPtr_
Definition: RASModelVariables.H:81
Foam::incompressible::RASModelVariables::copyAndRename
void copyAndRename(volScalarField &f1, volScalarField &f2)
Definition: RASModelVariables.C:164
Foam::dictionary::readCompat
bool readCompat(const word &keyword, std::initializer_list< std::pair< const char *, int >> compat, T &val, enum keyType::option=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:384
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:121
Foam::incompressible::RASModelVariables::TMVar1MeanPtr_
refPtr< volScalarField > TMVar1MeanPtr_
Definition: RASModelVariables.H:85
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:275
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::incompressible::RASModelVariables::TMVar2MeanPtr_
refPtr< volScalarField > TMVar2MeanPtr_
Definition: RASModelVariables.H:86
Foam::incompressible::RASModelVariables::nutJacobianVar1
virtual tmp< volScalarField > nutJacobianVar1(const singlePhaseTransportModel &laminarTransport) const
Return nut Jacobian wrt the TM vars.
Definition: RASModelVariables.C:300
Foam::incompressible::RASModelVariables::hasNut
bool hasNut() const
Definition: RASModelVariablesI.H:68
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::solverControl::average
bool average() const
Whether averaging is enabled or not.
Definition: solverControlI.H:107
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
Foam::incompressible::RASModelVariables::clone
autoPtr< RASModelVariables > clone() const
Clone.
Definition: RASModelVariables.C:237
Foam::incompressible::RASModelVariables::cloneRefPtr
refPtr< volScalarField > cloneRefPtr(const refPtr< volScalarField > &obj) const
Definition: RASModelVariables.C:148
Foam::refPtr
A class for managing references or pointers (no reference counting)
Definition: PtrList.H:60
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106