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 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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
50  {
51  Info<< "Storing initial values of turbulence variables" << endl;
52  if (hasTMVar1_)
53  {
54  TMVar1InitPtr_.reset
55  (
56  new volScalarField
57  (
58  TMVar1Inst().name()+"Init",TMVar1Inst()
59  )
60  );
61  }
62 
63  if (hasTMVar2_)
64  {
65  TMVar2InitPtr_.reset
66  (
67  new volScalarField
68  (
69  TMVar2Inst().name()+"Init",TMVar2Inst()
70  )
71  );
72  }
73 
74  if (hasNut_)
75  {
76  nutInitPtr_.reset
77  (
78  new volScalarField
79  (
80  nutRefInst().name()+"Init",nutRefInst()
81  )
82  );
83  }
84  }
85 }
86 
87 
89 {
90  if (solverControl_.average())
91  {
92  Info<< "Allocating mean values of turbulence variables" << endl;
93  if (hasTMVar1_)
94  {
95  TMVar1MeanPtr_.reset
96  (
97  new volScalarField
98  (
99  IOobject
100  (
101  TMVar1Inst().name()+"Mean",
102  mesh_.time().timeName(),
103  mesh_,
106  ),
107  TMVar1Inst()
108  )
109  );
110  }
111  if (hasTMVar2_)
112  {
113  TMVar2MeanPtr_.reset
114  (
115  new volScalarField
116  (
117  IOobject
118  (
119  TMVar2Inst().name()+"Mean",
120  mesh_.time().timeName(),
121  mesh_,
124  ),
125  TMVar2Inst()
126  )
127  );
128  }
129 
130  if (hasNut_)
131  {
132  nutMeanPtr_.reset
133  (
134  new volScalarField
135  (
136  IOobject
137  (
138  nutRefInst().name()+"Mean",
139  mesh_.time().timeName(),
140  mesh_,
143  ),
144  nutRefInst()
145  )
146  );
147  }
148  }
149 }
150 
151 
154 {
155  autoTmp returnField(nullptr);
156  if (source.valid() && source().valid())
157  {
158  const volScalarField& sf = source()();
159  DebugInfo
160  << "Cloning " << sf.name() << endl;
161  const word timeName = mesh_.time().timeName();
162  returnField.reset
163  (
165  (
166  new volScalarField(sf.name() + timeName, sf)
167  )
168  );
169  }
170  return returnField;
171 }
172 
173 
175 (
176  volScalarField& f1,
177  volScalarField& f2
178 )
179 {
180  f1 == f2;
181  const word name1 = f1.name();
182  const word name2 = f2.name();
183 
184  // Extra rename to avoid databese collision
185  f2.rename("temp");
186  f1.rename(name2);
187  f2.rename(name1);
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
194 (
195  const fvMesh& mesh,
196  const solverControl& SolverControl
197 )
198 :
199  mesh_(mesh),
200  solverControl_(SolverControl),
201  hasTMVar1_(false),
202  hasTMVar2_(false),
203  hasNut_(false),
204  hasDist_(false),
205  TMVar1Ptr_(nullptr),
206  TMVar2Ptr_(nullptr),
207  nutPtr_(nullptr),
208  dPtr_(nullptr),
209  TMVar1BaseName_(word::null),
210  TMVar2BaseName_(word::null),
211  nutBaseName_("nut"),
212  TMVar1InitPtr_(nullptr),
213  TMVar2InitPtr_(nullptr),
214  nutInitPtr_(nullptr),
215  TMVar1MeanPtr_(nullptr),
216  TMVar2MeanPtr_(nullptr),
217  nutMeanPtr_(nullptr)
218 {}
219 
220 
222 (
223  const RASModelVariables& rmv
224 )
225 :
226  mesh_(rmv.mesh_),
227  solverControl_(rmv.solverControl_),
228  hasTMVar1_(rmv.hasTMVar1_),
229  hasTMVar2_(rmv.hasTMVar2_),
230  hasNut_(rmv.hasNut_),
231  hasDist_(rmv.hasDist_),
232  TMVar1Ptr_(cloneAutoTmp(rmv.TMVar1Ptr_)),
233  TMVar2Ptr_(cloneAutoTmp(rmv.TMVar2Ptr_)),
234  nutPtr_(cloneAutoTmp(rmv.nutPtr_)),
235  dPtr_(cloneAutoTmp(rmv.dPtr_)),
236  TMVar1BaseName_(rmv.TMVar1BaseName_),
237  TMVar2BaseName_(rmv.TMVar2BaseName_),
238  nutBaseName_(rmv.nutBaseName_),
239  TMVar1InitPtr_(nullptr),
240  TMVar2InitPtr_(nullptr),
241  nutInitPtr_(nullptr),
242  TMVar1MeanPtr_(nullptr),
243  TMVar2MeanPtr_(nullptr),
244  nutMeanPtr_(nullptr)
245 {
246 }
247 
248 
250 {
251  return autoPtr<RASModelVariables>::New(*this);
252 }
253 
254 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
255 
257 (
258  const fvMesh& mesh,
259  const solverControl& SolverControl
260 )
261 {
262  const IOdictionary modelDict
263  (
264  IOobject
265  (
267  mesh.time().constant(),
268  mesh,
271  false // Do not register
272  )
273  );
274 
275  const dictionary dict(modelDict.subOrEmptyDict("RAS"));
276 
277  const word modelType(dict.getOrDefault<word>("RASModel", "laminar"));
278 
279  Info<< "Creating references for RASModel variables : " << modelType << endl;
280 
281  auto cstrIter = dictionaryConstructorTablePtr_->cfind(modelType);
282 
283  if (!cstrIter.found())
284  {
286  (
287  dict,
288  "RASModelVariables",
289  modelType,
290  *dictionaryConstructorTablePtr_
291  ) << exit(FatalIOError);
292  }
293 
294  return autoPtr<RASModelVariables>(cstrIter()(mesh, SolverControl));
295 }
296 
297 
298 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
299 
301 {
302  return hasTMVar1_;
303 }
304 
305 
307 {
308  return hasTMVar2_;
309 }
310 
311 
313 {
314  return hasNut_;
315 }
316 
317 
319 {
320  return hasDist_;
321 }
322 
323 
325 {
326  return TMVar1BaseName_;
327 }
328 
329 
331 {
332  return TMVar2BaseName_;
333 }
334 
335 
337 {
338  return nutBaseName_;
339 }
340 
341 
343 {
345  {
346  return TMVar1MeanPtr_();
347  }
348  else
349  {
350  return TMVar1Ptr_()();
351  }
352 }
353 
354 
356 {
358  {
359  return TMVar1MeanPtr_();
360  }
361  else
362  {
363  return TMVar1Ptr_().constCast();
364  }
365 }
366 
367 
369 {
371  {
372  return TMVar2MeanPtr_();
373  }
374  else
375  {
376  return TMVar2Ptr_()();
377  }
378 }
379 
381 {
383  {
384  return TMVar2MeanPtr_();
385  }
386  else
387  {
388  return TMVar2Ptr_().constCast();
389  }
390 }
391 
393 {
395  {
396  return nutMeanPtr_();
397  }
398  else
399  {
400  return nutPtr_()();
401  }
402 }
403 
404 
406 {
408  {
409  return nutMeanPtr_();
410  }
411  else
412  {
413  return nutPtr_().constCast();
414  }
415 }
416 
417 
419 {
420  return dPtr_()();
421 }
422 
423 
425 {
426  return dPtr_().constCast();
427 }
428 
429 
431 {
432  return TMVar1Ptr_()();
433 }
434 
435 
437 {
438  return TMVar1Ptr_().constCast();
439 }
440 
441 
443 {
444  return TMVar2Ptr_()();
445 }
446 
447 
449 {
450  return TMVar2Ptr_().constCast();
451 }
452 
453 
455 {
456  return nutPtr_()();
457 }
458 
459 
461 {
462  return nutPtr_().constCast();
463 }
464 
465 
467 (
469 ) const
470 {
472  << "jutJacobianVar1 not implemented for the current turbulence model."
473  << "Returning zero field" << endl;
474 
475  tmp<volScalarField> nutJacobian
476  (
477  new volScalarField
478  (
479  IOobject
480  (
481  "nutJacobianVar1",
482  mesh_.time().timeName(),
483  mesh_,
486  ),
487  mesh_,
489  )
490  );
491  return nutJacobian;
492 }
493 
494 
496 (
498 ) const
499 {
501  << "nutJacobianVar2 not implemented for the current turbulence model."
502  << "Returning zero field" << endl;
503 
504  tmp<volScalarField> nutJacobian
505  (
506  new volScalarField
507  (
508  IOobject
509  (
510  "nutJacobianVar2",
511  mesh_.time().timeName(),
512  mesh_,
515  ),
516  mesh_,
518  )
519  );
520  return nutJacobian;
521 }
522 
524 {
526  {
527  if (hasTMVar1_)
528  {
529  TMVar1Inst() == TMVar1InitPtr_();
530  }
531  if (hasTMVar2_)
532  {
533  TMVar2Inst() == TMVar2InitPtr_();
534  }
535  if (hasNut_)
536  {
537  nutRefInst() == nutInitPtr_();
538  }
539  }
540 }
541 
542 
544 {
545  if (solverControl_.average())
546  {
547  Info<< "Reseting mean turbulent fields to zero" << endl;
548 
549  // Reset fields to zero
550  if (hasTMVar1_)
551  {
552  TMVar1MeanPtr_() ==
553  dimensionedScalar(TMVar1Inst().dimensions(), Zero);
554  }
555  if (hasTMVar2_)
556  {
557  TMVar2MeanPtr_() ==
558  dimensionedScalar(TMVar2Inst().dimensions(), Zero);
559  }
560  if (hasNut_)
561  {
562  nutMeanPtr_() == dimensionedScalar(nutRefInst().dimensions(), Zero);
563  }
564  }
565 }
566 
567 
569 {
571  {
572  const label iAverageIter = solverControl_.averageIter();
573  scalar avIter(iAverageIter);
574  scalar oneOverItP1 = 1./(avIter + 1);
575  scalar mult = avIter*oneOverItP1;
576  if (hasTMVar1_)
577  {
578  TMVar1MeanPtr_() ==
579  TMVar1MeanPtr_()*mult + TMVar1Inst()*oneOverItP1;
580  }
581  if (hasTMVar2_)
582  {
583  TMVar2MeanPtr_() ==
584  TMVar2MeanPtr_()*mult + TMVar2Inst()*oneOverItP1;
585  }
586  if (hasNut_)
587  {
588  nutMeanPtr_() == nutMeanPtr_()*mult + nutRefInst()*oneOverItP1;
589  }
590  }
591 }
592 
593 
595 (
597  const volVectorField& U
598 ) const
599 {
601  (
603  (
604  IOobject
605  (
606  "devRhoReff",
607  mesh_.time().timeName(),
608  mesh_,
611  ),
612  -(laminarTransport.nu() + nutRef())*dev(twoSymm(fvc::grad(U)))
613  )
614  );
615 }
616 
617 
619 (
621 )
622 {
623  if (hasTMVar1())
624  {
625  TMVar1Ptr_().constCast().correctBoundaryConditions();
626  if (solverControl_.average())
627  {
628  TMVar1MeanPtr_().correctBoundaryConditions();
629  }
630  }
631 
632  if (hasTMVar2())
633  {
634  TMVar2Ptr_().constCast().correctBoundaryConditions();
635  if (solverControl_.average())
636  {
637  TMVar2MeanPtr_().correctBoundaryConditions();
638  }
639  }
640 
641  if (hasNut())
642  {
643  nutPtr_().constCast().correctBoundaryConditions();
644  if (solverControl_.average())
645  {
646  nutMeanPtr_().correctBoundaryConditions();
647  }
648  }
649 }
650 
651 
653 {
654  if (rmv.hasTMVar1() && hasTMVar1_)
655  {
657  }
658 
659  if (rmv.hasTMVar2() && hasTMVar2_)
660  {
662  }
663 
664  if (rmv.hasNut() && hasNut_)
665  {
666  copyAndRename(nutRef(), rmv.nutRef());
667  }
668 
669  if (rmv.hasDist() && hasDist_)
670  {
671  copyAndRename(d(), rmv.d());
672  }
673 }
674 
675 
676 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
677 
678 } // End namespace incompressible
679 } // End namespace Foam
680 
681 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::autoPtr::New
static autoPtr< T > New(Args &&... args)
Construct autoPtr of T with forwarding arguments.
Foam::incompressible::RASModelVariables::nutInitPtr_
autoPtr< volScalarField > nutInitPtr_
Definition: RASModelVariables.H:95
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::autoPtr::reset
void reset(T *p=nullptr) noexcept
Delete managed object and set to new given pointer.
Definition: autoPtrI.H:158
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::incompressible::RASModelVariables::nutBaseName
const word & nutBaseName() const
Definition: RASModelVariables.C:336
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: RASModelVariables.C:418
Foam::incompressible::RASModelVariables::TMVar2Ptr_
autoTmp TMVar2Ptr_
Definition: RASModelVariables.H:83
Foam::incompressible::RASModelVariables::TMVar1Inst
const volScalarField & TMVar1Inst() const
return references to instantaneous turbulence fields
Definition: RASModelVariables.C:430
Foam::incompressible::RASModelVariables::solverControl_
const solverControl & solverControl_
Definition: RASModelVariables.H:73
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:58
Foam::incompressible::RASModelVariables::TMVar1InitPtr_
autoPtr< volScalarField > TMVar1InitPtr_
Definition: RASModelVariables.H:93
Foam::tmp< volScalarField >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::incompressible::RASModelVariables::TMVar1MeanPtr_
autoPtr< volScalarField > TMVar1MeanPtr_
Definition: RASModelVariables.H:98
Foam::incompressible::RASModelVariables::dPtr_
autoTmp dPtr_
Definition: RASModelVariables.H:85
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::incompressible::RASModelVariables::restoreInitValues
void restoreInitValues()
Restore turbulent fields to their initial values.
Definition: RASModelVariables.C:523
Foam::autoPtr::valid
bool valid() const noexcept
True if the managed pointer is non-null.
Definition: autoPtrI.H:107
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:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::solverControl::storeInitValues
bool storeInitValues() const
Re-initialize.
Definition: solverControlI.H:51
Foam::incompressible::RASModelVariables::TMVar2BaseName_
word TMVar2BaseName_
Definition: RASModelVariables.H:88
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:57
Foam::incompressible::RASModelVariables::correctBoundaryConditions
virtual void correctBoundaryConditions(const incompressible::turbulenceModel &turbulence)
correct bounday conditions of turbulent fields
Definition: RASModelVariables.C:619
Foam::incompressible::defineRunTimeSelectionTable
defineRunTimeSelectionTable(adjointSensitivity, dictionary)
Foam::incompressible::RASModelVariables::RASModelVariables
RASModelVariables(const fvMesh &mesh, const solverControl &SolverControl)
Construct from components.
Definition: RASModelVariables.C:194
Foam::incompressible::RASModelVariables::computeMeanFields
void computeMeanFields()
Compute mean fields on the fly.
Definition: RASModelVariables.C:568
Foam::incompressible::RASModelVariables::nutMeanPtr_
autoPtr< volScalarField > nutMeanPtr_
Definition: RASModelVariables.H:100
Foam::incompressible::RASModelVariables::transfer
virtual void transfer(RASModelVariables &rmv)
Transfer turbulence fields from an another object.
Definition: RASModelVariables.C:652
Foam::incompressible::RASModelVariables::nutRefInst
const volScalarField & nutRefInst() const
Definition: RASModelVariables.C:454
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:380
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::incompressible::RASModelVariables::hasDist
bool hasDist() const
Definition: RASModelVariables.C:318
Foam::incompressible::RASModelVariables::allocateMeanFields
void allocateMeanFields()
Definition: RASModelVariables.C:88
Foam::incompressible::RASModelVariables::resetMeanFields
void resetMeanFields()
Reset mean fields to zero.
Definition: RASModelVariables.C:543
RASModelVariables.H
Foam::incompressible::RASModelVariables::hasTMVar2
bool hasTMVar2() const
Definition: RASModelVariables.C:306
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::incompressible::RASModelVariables::nutBaseName_
word nutBaseName_
Definition: RASModelVariables.H:89
Foam::incompressible::RASModelVariables::hasTMVar2_
bool hasTMVar2_
Definition: RASModelVariables.H:79
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:72
Foam::incompressible::RASModelVariables::nutJacobianVar2
virtual tmp< volScalarField > nutJacobianVar2(const singlePhaseTransportModel &laminarTransport) const
Definition: RASModelVariables.C:496
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::solverControl::useAveragedFields
bool useAveragedFields() const
Definition: solverControlI.H:94
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:257
Foam::incompressible::RASModelVariables::TMVar2MeanPtr_
autoPtr< volScalarField > TMVar2MeanPtr_
Definition: RASModelVariables.H:99
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::TMVar1BaseName_
word TMVar1BaseName_
Definition: RASModelVariables.H:87
Foam::solverControl::averageIter
label & averageIter()
Return average iteration index reference.
Definition: solverControlI.H:63
Foam::incompressible::RASModelVariables::TMVar2Inst
const volScalarField & TMVar2Inst() const
Definition: RASModelVariables.C:442
timeName
word timeName
Definition: getTimeIndex.H:3
dict
dictionary dict
Definition: searchingEngine.H:14
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::nutRef
const volScalarField & nutRef() const
Definition: RASModelVariables.C:392
Foam::incompressible::RASModelVariables::hasTMVar1
bool hasTMVar1() const
Bools to idenify which turbulent fields are present.
Definition: RASModelVariables.C:300
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::incompressible::RASModelVariables::hasDist_
bool hasDist_
Definition: RASModelVariables.H:81
Foam::dictionary::subOrEmptyDict
dictionary subOrEmptyDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX, const bool mandatory=false) const
Definition: dictionary.C:603
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:595
Foam::autoPtr< tmp< volScalarField > >
U
U
Definition: pEqn.H:72
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::incompressible::RASModelVariables::cloneAutoTmp
autoTmp cloneAutoTmp(const autoTmp &source) const
Definition: RASModelVariables.C:153
Foam::incompressible::RASModelVariables::hasTMVar1_
bool hasTMVar1_
Definition: RASModelVariables.H:78
Foam::incompressible::RASModelVariables::TMVar1BaseName
const word & TMVar1BaseName() const
Turbulence field names.
Definition: RASModelVariables.C:324
Foam::incompressible::RASModelVariables::copyAndRename
void copyAndRename(volScalarField &f1, volScalarField &f2)
Definition: RASModelVariables.C:175
Foam::incompressible::RASModelVariables::TMVar1Ptr_
autoTmp TMVar1Ptr_
Definition: RASModelVariables.H:82
Foam::incompressible::RASModelVariables::TMVar2BaseName
const word & TMVar2BaseName() const
Definition: RASModelVariables.C:330
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:121
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::incompressible::RASModelVariables::TMVar1
const volScalarField & TMVar1() const
Return references to turbulence fields.
Definition: RASModelVariables.C:342
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
Foam::incompressible::RASModelVariables::hasNut_
bool hasNut_
Definition: RASModelVariables.H:80
Foam::incompressible::RASModelVariables::nutJacobianVar1
virtual tmp< volScalarField > nutJacobianVar1(const singlePhaseTransportModel &laminarTransport) const
Return nut Jacobian wrt the TM vars.
Definition: RASModelVariables.C:467
Foam::incompressible::RASModelVariables::hasNut
bool hasNut() const
Definition: RASModelVariables.C:312
Foam::incompressible::RASModelVariables::TMVar2InitPtr_
autoPtr< volScalarField > TMVar2InitPtr_
Definition: RASModelVariables.H:94
Foam::incompressible::RASModelVariables::TMVar2
const volScalarField & TMVar2() const
Definition: RASModelVariables.C:368
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::incompressible::RASModelVariables::nutPtr_
autoTmp nutPtr_
Definition: RASModelVariables.H:84
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
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:249
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106