atmLengthScaleTurbSource.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 
30 #include "geometricOneField.H"
32 
33 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace fv
38 {
39  defineTypeNameAndDebug(atmLengthScaleTurbSource, 0);
40  addToRunTimeSelectionTable(option, atmLengthScaleTurbSource, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 Foam::tmp<Foam::volScalarField::Internal>Foam::fv::atmLengthScaleTurbSource::
48 calcC1Star
49 (
50  const volScalarField::Internal& k,
51  const volScalarField::Internal& epsilon
52 ) const
53 {
54  // Mixing-length scale estimation (P:Eq. 10.37 & p. 374)
55  tmp<volScalarField::Internal> L_(pow(Cmu_, 0.75)*pow(k, 1.5)/epsilon);
56 
57  // (AC:Eq. 16) wherein the exponentiation "n_" is not present.
58  // "n_" is an ad-hoc implementation.
59  return (C2_ - C1_)*pow(L_/Lmax_, n_);
60 }
61 
62 
63 Foam::tmp<Foam::volScalarField::Internal> Foam::fv::atmLengthScaleTurbSource::
64 calcGammaStar
65 (
67  const volScalarField::Internal& omega,
70 ) const
71 {
72  // (L:Eq. 3.20)
73  tmp<volScalarField::Internal> L_(sqrt(k)/(pow025(Cmu_)*omega));
74 
75  // (L:Eq. 3.34)
76  return (gamma - beta)*pow(L_/Lmax_, n_);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
83 (
84  const word& sourceName,
85  const word& modelType,
86  const dictionary& dict,
87  const fvMesh& mesh
88 )
89 :
90  fv::cellSetOption(sourceName, modelType, dict, mesh),
91  isEpsilon_(false),
92  rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
93  Lmax_
94  (
96  (
97  dimLength,
98  coeffs_.getCheckOrDefault<scalar>
99  (
100  "Lmax",
101  41.575,
102  [&](const scalar Lmax){ return Lmax > SMALL; }
103  )
104  )
105  ),
106  n_
107  (
109  (
110  dimless,
111  coeffs_.getCheckOrDefault<scalar>
112  (
113  "n",
114  3.0,
115  [&](const scalar n){ return n > SMALL; }
116  )
117  )
118  ),
119  Cmu_(Zero),
120  C1_(Zero),
121  C2_(Zero),
122  C3_(Zero)
123 {
124  const auto* turbPtr =
125  mesh_.findObject<turbulenceModel>
126  (
128  );
129 
130  if (!turbPtr)
131  {
133  << "Unable to find a turbulence model."
134  << abort(FatalError);
135  }
136 
137  fieldNames_.resize(1);
138 
139  tmp<volScalarField> tepsilon = turbPtr->epsilon();
140  tmp<volScalarField> tomega = turbPtr->omega();
141 
142  if (tepsilon.is_reference())
143  {
144  isEpsilon_ = true;
145  fieldNames_[0] = tepsilon().name();
146 
147  const dictionary& turbDict = turbPtr->coeffDict();
148  Cmu_.read("Cmu", turbDict);
149  C1_.read("C1", turbDict);
150  C2_.read("C2", turbDict);
151  C3_.read("C3", turbDict);
152  }
153  else if (tomega.is_reference())
154  {
155  isEpsilon_ = false;
156  fieldNames_[0] = tomega().name();
157 
158  const dictionary& turbDict = turbPtr->coeffDict();
159  Cmu_.read("betaStar", turbDict);
160  }
161  else
162  {
164  << "Needs either epsilon or omega field."
165  << abort(FatalError);
166  }
167 
169 
170  Log << " Applying atmLengthScaleTurbSource to: " << fieldNames_[0]
171  << endl;
172 }
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
178 (
179  fvMatrix<scalar>& eqn,
180  const label fieldi
181 )
182 {
183  if (isEpsilon_)
184  {
185  atmLengthScaleTurbSourceEpsilon
186  (
189  eqn,
190  fieldi
191  );
192  }
193  else
194  {
195  atmLengthScaleTurbSourceOmega
196  (
199  eqn,
200  fieldi
201  );
202  }
203 }
204 
205 
207 (
208  const volScalarField& rho,
209  fvMatrix<scalar>& eqn,
210  const label fieldi
211 )
212 {
213  if (isEpsilon_)
214  {
215  atmLengthScaleTurbSourceEpsilon(geometricOneField(), rho, eqn, fieldi);
216  }
217  else
218  {
219  atmLengthScaleTurbSourceOmega(geometricOneField(), rho, eqn, fieldi);
220  }
221 }
222 
223 
225 (
226  const volScalarField& alpha,
227  const volScalarField& rho,
228  fvMatrix<scalar>& eqn,
229  const label fieldi
230 )
231 {
232  if (isEpsilon_)
233  {
234  atmLengthScaleTurbSourceEpsilon(alpha, rho, eqn, fieldi);
235  }
236  else
237  {
238  atmLengthScaleTurbSourceOmega(alpha, rho, eqn, fieldi);
239  }
240 }
241 
242 
243 // ************************************************************************* //
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to epsilon or omega equation for incompressible flow computations.
dictionary dict
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
atmLengthScaleTurbSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
dimensionedScalar pow025(const dimensionedScalar &ds)
label k
Boltzmann constant.
const dimensionSet dimless
Dimensionless.
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
static const word propertiesName
Default name of the turbulence properties dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
errorManip< error > abort(error &err)
Definition: errorManip.H:139
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
scalar epsilon
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
#define Log
Definition: PDRblock.C:28
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
label n
Intermediate abstract class for handling cell-set options for the derived fvOptions.
const scalar gamma
Definition: EEqn.H:9
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127