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  obr_.store
53  (
54  new FieldType
55  (
56  IOobject
57  (
58  localName,
59  obr_.time().timeName(),
60  obr_,
61  IOobject::READ_IF_PRESENT,
62  IOobject::NO_WRITE
63  ),
64  tvalue
65  )
66  );
67  }
68 }
69 
70 
71 template<class Model>
74 (
75  const Model& model
76 ) const
77 {
78  const dimensionedScalar omega0(dimless/dimTime, SMALL);
79 
81  (
82  "nuTilda.tmp",
83  model.k()/(model.omega() + omega0)
84  );
85 }
86 
87 
88 template<class Model>
91 (
92  const Model& model
93 ) const
94 {
95  // (P:Eq. 10.37)
96  const scalar Cmu = 0.09;
97  const dimensionedScalar eps0(sqr(dimVelocity)/dimTime, SMALL);
98 
100  (
101  "L.tmp",
102  pow(Cmu, 0.75)*pow(model.k(), 1.5)/(model.epsilon() + eps0)
103  );
104 }
105 
106 
107 template<class Model>
110 (
111  const Model& model
112 ) const
113 {
114  // (P:p. 183)
115  // root-mean-square of velocity fluctuations - isotropic turbulence
116  const volScalarField uPrime(sqrt((2.0/3.0)*model.k()));
117  const dimensionedScalar U0(dimVelocity, SMALL);
118 
120  (
121  "I.tmp",
122  uPrime/max(max(uPrime, mag(model.U())), U0)
123  );
124 }
125 
126 
127 // ************************************************************************* //
volFields.H
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::functionObjects::turbulenceFields::processField
void processField(const word &fieldName, const tmp< GeometricField< Type, fvPatchField, volMesh >> &tvalue)
Process the turbulence field.
Definition: turbulenceFieldsTemplates.C:35
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::functionObjects::turbulenceFields::L
tmp< volScalarField > L(const Model &model) const
Return integral length scale, L, calculated from k and epsilon.
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::dimensioned< scalar >
Foam::functionObjects::turbulenceFields::nuTilda
tmp< volScalarField > nuTilda(const Model &model) const
Return nuTilda calculated from k and omega.
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::functionObjects::turbulenceFields::I
tmp< volScalarField > I(const Model &model) const
Return turbulence intensity, I, calculated from k and U.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189