turbulenceFieldsTemplates.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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2018-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const word& fieldName,
38 )
39 {
41 
42  const word localName(IOobject::scopedName(prefix_, fieldName));
43 
44  FieldType* fldPtr = obr_.getObjectPtr<FieldType>(localName);
45 
46  if (fldPtr)
47  {
48  (*fldPtr) == tvalue();
49  }
50  else
51  {
52  fldPtr = new FieldType
53  (
54  IOobject
55  (
56  localName,
57  obr_.time().timeName(),
58  obr_,
62  ),
63  tvalue
64  );
65 
66  obr_.store(fldPtr);
67  }
68 }
69 
70 
71 template<class Model>
74 (
75  const Model& model
76 ) const
77 {
78  const dimensionedScalar omega0(dimless/dimTime, SMALL);
79 
80  return volScalarField::New
81  (
82  "nuTilda.tmp",
84  model.k()/(model.omega() + omega0)
85  );
86 }
87 
88 
89 template<class Model>
92 (
93  const Model& model
94 ) const
95 {
96  // (P:Eq. 10.37)
97  const scalar Cmu = 0.09;
98  const dimensionedScalar eps0(sqr(dimVelocity)/dimTime, SMALL);
99 
100  return volScalarField::New
101  (
102  "L.tmp",
104  pow(Cmu, 0.75)*pow(model.k(), 1.5)/(model.epsilon() + eps0)
105  );
106 }
107 
108 
109 template<class Model>
112 (
113  const Model& model
114 ) const
115 {
116  // (P:p. 183)
117  // root-mean-square of velocity fluctuations - isotropic turbulence
118  const volScalarField uPrime(sqrt((2.0/3.0)*model.k()));
119  const dimensionedScalar U0(dimVelocity, SMALL);
120 
121  return volScalarField::New
122  (
123  "I.tmp",
125  uPrime/max(max(uPrime, mag(model.U())), U0)
126  );
127 }
128 
129 
130 // ************************************************************************* //
tmp< volScalarField > I(const Model &model) const
Return turbulence intensity, I, calculated from k and U.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
Generic GeometricField class.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
word prefix_
Name of output-field prefix.
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
const Time & time() const noexcept
Return time registry.
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< volScalarField > nuTilda(const Model &model) const
Return nuTilda calculated from k and omega.
void processField(const word &fieldName, const tmp< GeometricField< Type, fvPatchField, volMesh >> &tvalue)
Process the turbulence field.
const objectRegistry & obr_
Reference to the region objectRegistry.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Reading is optional [identical to READ_IF_PRESENT].
tmp< volScalarField > L(const Model &model) const
Return integral length scale, L, calculated from k and epsilon.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Request registration (bool: true)
Do not request registration (bool: false)
const dimensionSet dimVelocity