atmBuoyancyTurbSource.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) 2020 ENERCON GmbH
9  Copyright (C) 2020-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 "atmBuoyancyTurbSource.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
38  defineTypeNameAndDebug(atmBuoyancyTurbSource, 0);
39  addToRunTimeSelectionTable(option, atmBuoyancyTurbSource, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 
46 void Foam::fv::atmBuoyancyTurbSource::calcB()
47 {
48  //- Temperature field [K]
50 
51  //- Kinematic turbulent thermal conductivity field [m2/s]
52  const volScalarField& alphat =
54 
55  // (ARAL:Eq. 7)
56  B_ = beta_*alphat()*(fvc::grad(T) & g_)();
57 }
58 
59 
61 Foam::fv::atmBuoyancyTurbSource::calcC3
62 (
66 ) const
67 {
68  // Gradient Richardson number (ARAL:p. 4)
69  const volScalarField::Internal Rig
70  (
71  -B_/(G + dimensionedScalar(G.dimensions(), SMALL))
72  );
73 
74  // Mixing-length scale estimation (P:Eq. 10.37 & p. 374) normalised by Lmax_
75  const volScalarField::Internal LbyLmax
76  (
77  (pow(Cmu_, 0.75)/Lmax_)*pow(k, 1.5)/epsilon
78  );
79 
80  // (ARAL:Eq. 10), with a typo of (C2_) instead of using (C2_ - 1.0)
81  volScalarField::Internal alphaB(1.0 - LbyLmax);
82 
83  alphaB =
84  neg0(Rig)*(1.0 - (1.0 + (C2_ - 1.0)/(C2_ - C1_))*LbyLmax)
85  + pos(Rig)*(1.0 - LbyLmax);
86 
87  // (SKL:Eq. 18, rhs-term:3); (ARAL:Eq. 5, rhs-term:3) has a typo
88  return (C1_ - C2_)*alphaB + 1.0;
89 }
90 
91 
93 Foam::fv::atmBuoyancyTurbSource::calcC3
94 (
96  const volScalarField::Internal& omega,
100 ) const
101 {
102  // Gradient Richardson number (ARAL:p. 4)
103  const volScalarField::Internal Rig
104  (
105  -B_/(G + dimensionedScalar(G.dimensions(), SMALL))
106  );
107 
108  // Mixing-length scale estimation (L:Eq. 3.20) normalised by Lmax_
109  const volScalarField::Internal LbyLmax
110  (
111  (1.0/(pow025(Cmu_)*Lmax_))*sqrt(k)/omega
112  );
113 
114  // (ARAL:Eq. 10)
115  volScalarField::Internal alphaB(1.0 - LbyLmax);
116 
117  alphaB =
118  neg0(Rig)*(1.0 - (1.0 + beta/(beta - gamma))*LbyLmax)
119  + pos(Rig)*(1.0 - LbyLmax);
120 
121  // (SKL:Eq. 19, rhs-term:3); (ARAL:Eq. 5, rhs-term:3) has a typo
122  return (gamma - beta)*alphaB;
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127 
129 (
130  const word& sourceName,
131  const word& modelType,
132  const dictionary& dict,
133  const fvMesh& mesh
134 )
135 :
136  fv::cellSetOption(sourceName, modelType, dict, mesh),
137  isEpsilon_(true),
138  rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
139  Lmax_
140  (
142  (
143  dimLength,
144  coeffs_.getCheckOrDefault<scalar>
145  (
146  "Lmax",
147  41.575,
148  [&](const scalar Lmax){ return Lmax > SMALL; }
149  )
150  )
151  ),
152  beta_
153  (
155  (
157  coeffs_.getCheckOrDefault<scalar>
158  (
159  "beta",
160  3.3e-3,
161  [&](const scalar x){ return x > SMALL; }
162  )
163  )
164  ),
165  Cmu_(Zero),
166  C1_(Zero),
167  C2_(Zero),
168  g_
169  (
170  "g",
172  meshObjects::gravity::New(mesh_.time()).value()
173  ),
174  B_
175  (
176  IOobject
177  (
178  "B",
179  mesh.time().timeName(),
180  mesh,
183  ),
184  mesh,
186  )
187 {
188  const auto* turbPtr =
189  mesh_.findObject<turbulenceModel>
190  (
192  );
193 
194  if (!turbPtr)
195  {
197  << "Unable to find a turbulence model."
198  << abort(FatalError);
199  }
200 
201  fieldNames_.resize(2);
202 
203  tmp<volScalarField> tepsilon = turbPtr->epsilon();
204  tmp<volScalarField> tomega = turbPtr->omega();
205 
206  if (!tepsilon.isTmp())
207  {
208  fieldNames_[0] = tepsilon().name();
209 
210  const dictionary& turbDict = turbPtr->coeffDict();
211 
212  Cmu_.read("Cmu", turbDict);
213  C1_.read("C1", turbDict);
214  C2_.read("C2", turbDict);
215  }
216  else if (!tomega.isTmp())
217  {
218  isEpsilon_ = false;
219  fieldNames_[0] = tomega().name();
220 
221  const dictionary& turbDict = turbPtr->coeffDict();
222 
223  Cmu_.read("betaStar", turbDict);
224  }
225  else
226  {
228  << "Unable to find neither epsilon nor omega field." << nl
229  << "atmBuoyancyTurbSource needs either epsilon or omega field."
230  << abort(FatalError);
231  }
232 
233  fieldNames_[1] = turbPtr->k()().name();
234 
236 
237  Log << " Applying atmBuoyancyTurbSource to: "
238  << fieldNames_[0] << " and " << fieldNames_[1]
239  << endl;
240 }
241 
242 
243 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244 
246 (
247  fvMatrix<scalar>& eqn,
248  const label fieldi
249 )
250 {
251  if (fieldi == 1)
252  {
253  atmBuoyancyTurbSourceK
254  (
257  eqn,
258  fieldi
259  );
260  return;
261  }
262 
263  calcB();
264 
265  if (isEpsilon_)
266  {
267  atmBuoyancyTurbSourceEpsilon
268  (
271  eqn,
272  fieldi
273  );
274  }
275  else
276  {
277  atmBuoyancyTurbSourceOmega
278  (
281  eqn,
282  fieldi
283  );
284  }
285 }
286 
287 
289 (
290  const volScalarField& rho,
291  fvMatrix<scalar>& eqn,
292  const label fieldi
293 )
294 {
295  if (fieldi == 1)
296  {
297  atmBuoyancyTurbSourceK(geometricOneField(), rho, eqn, fieldi);
298  return;
299  }
300 
301  calcB();
302 
303  if (isEpsilon_)
304  {
305  atmBuoyancyTurbSourceEpsilon(geometricOneField(), rho, eqn, fieldi);
306  }
307  else
308  {
309  atmBuoyancyTurbSourceOmega(geometricOneField(), rho, eqn, fieldi);
310  }
311 }
312 
313 
315 (
316  const volScalarField& alpha,
317  const volScalarField& rho,
318  fvMatrix<scalar>& eqn,
319  const label fieldi
320 )
321 {
322  if (fieldi == 1)
323  {
324  atmBuoyancyTurbSourceK(alpha, rho, eqn, fieldi);
325  return;
326  }
327 
328  calcB();
329 
330  if (isEpsilon_)
331  {
332  atmBuoyancyTurbSourceEpsilon(alpha, rho, eqn, fieldi);
333  }
334  else
335  {
336  atmBuoyancyTurbSourceOmega(alpha, rho, eqn, fieldi);
337  }
338 }
339 
340 
341 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Log
#define Log
Definition: PDRblock.C:35
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::cellSetOption
Intermediate abstract class for handling cell-set options for the derived fvOptions.
Definition: cellSetOption.H:163
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
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::fv::option::resetApplied
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:48
Foam::fv::atmBuoyancyTurbSource::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Definition: atmBuoyancyTurbSource.C:246
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::neg0
dimensionedScalar neg0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:210
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
Foam::tmp::isTmp
bool isTmp() const noexcept
Identical to is_pointer()
Definition: tmp.H:295
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rho
rho
Definition: readInitialConditions.H:88
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:139
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::GeometricField< scalar, fvPatchField, volMesh >::Internal
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Definition: GeometricField.H:107
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
atmBuoyancyTurbSource.H
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::dictionary::read
bool read(Istream &is)
Read dictionary from Istream. Discards the header.
Definition: dictionaryIO.C:141
Foam::objectRegistry::lookupObjectRef
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:478
fv
labelList fv(nPoints)
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
gamma
const scalar gamma
Definition: EEqn.H:9
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::dimensioned::dimensions
const dimensionSet & dimensions() const
Return const reference to dimensions.
Definition: dimensionedType.C:420
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fv::atmBuoyancyTurbSource::atmBuoyancyTurbSource
atmBuoyancyTurbSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: atmBuoyancyTurbSource.C:129
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:177