buoyancyTurbSource.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-2023 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "buoyancyTurbSource.H"
29 #include "gravityMeshObject.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
38  defineTypeNameAndDebug(buoyancyTurbSource, 0);
39  addToRunTimeSelectionTable(option, buoyancyTurbSource, dictionary);
40 }
41 }
42 
43 
44 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
45 
46 Foam::tmp<Foam::volScalarField::Internal> Foam::fv::buoyancyTurbSource::
47 B() const
48 {
49  const auto& alphat = mesh_.lookupObjectRef<volScalarField>(alphatName_);
50  const auto& T = mesh_.lookupObjectRef<volScalarField>(Tname_);
51 
52  // (BMA:Eq. 8)
53  return beta_*alphat()*(fvc::grad(T) & g_)();
54 }
55 
56 
57 void Foam::fv::buoyancyTurbSource::buoyancyTurbSourceEpsilon
58 (
59  fvMatrix<scalar>& eqn
60 ) const
61 {
62  const auto* turbPtr =
63  mesh_.findObject<turbulenceModel>
64  (
66  );
67  const dictionary& turbDict = turbPtr->coeffDict();
68  const dimensionedScalar C1(turbDict.getOrDefault<scalar>("C1", 1.44));
69  const dimensionedScalar Cmu(turbDict.getOrDefault<scalar>("Cmu", 0.09));
70 
71  const volScalarField& epsilon = eqn.psi();
72  const volScalarField::Internal& k = turbPtr->k()();
73  const volVectorField::Internal& U = turbPtr->U()();
74  const dimensionedScalar k0(k.dimensions(), SMALL);
75 
76  // (BMA:Eq. 9)
77  const vector gHat(g_.value()/mag(g_.value()));
78 
79  volScalarField::Internal v(gHat & U);
81  (
82  mag(U - gHat*v)
84  );
85 
86  // (BMA:Eq. 6)
87  eqn -= fvm::SuSp(C1*tanh(mag(u/v))*B()/(k + k0), epsilon);
88 }
89 
90 
91 void Foam::fv::buoyancyTurbSource::buoyancyTurbSourceOmega
92 (
93  fvMatrix<scalar>& eqn
94 ) const
95 {
96  const auto* turbPtr =
97  mesh_.findObject<turbulenceModel>
98  (
100  );
101 
102  const volScalarField::Internal& nut = turbPtr->nut()();
103  const auto& gamma =
104  mesh_.lookupObjectRef<volScalarField::Internal>
105  (
106  IOobject::scopedName(turbPtr->type(), "gamma")
107  );
108 
109  // (Heuristic approximation to BMA:Eq. 6)
110  eqn -= gamma/(nut + dimensionedScalar(nut.dimensions(), SMALL))*B();
111 }
112 
113 
114 void Foam::fv::buoyancyTurbSource::buoyancyTurbSourceK
115 (
116  fvMatrix<scalar>& eqn
117 ) const
118 {
119  const volScalarField& k = eqn.psi();
120  const dimensionedScalar k0(k.dimensions(), SMALL);
121 
122  // (BMA: Eq. 5)
123  eqn -= fvm::Sp(B()/(k() + k0), k);
124 }
125 
126 
127 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
128 
130 (
131  const word& sourceName,
132  const word& modelType,
133  const dictionary& dict,
134  const fvMesh& mesh
135 )
136 :
137  fv::cellSetOption(sourceName, modelType, dict, mesh),
138  isEpsilon_(false),
139  rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
140  alphatName_(coeffs_.getOrDefault<word>("alphat", "alphat")),
141  Tname_(coeffs_.getOrDefault<word>("T", "T")),
142  beta_
143  (
145  (
147  coeffs_.getCheckOrDefault<scalar>
148  (
149  "beta",
150  3.3e-3,
151  [=](const scalar x){ return x > SMALL; }
152  )
153  )
154  ),
155  g_
156  (
158  meshObjects::gravity::New(mesh_.time()).value()
159  )
160 {
161  if (mag(g_.value()) < SMALL)
162  {
164  << "Gravitational field cannot be equal to or less than zero"
165  << exit(FatalError);
166  }
167 
168  const auto* turbPtr =
169  mesh_.findObject<turbulenceModel>
170  (
172  );
173 
174  if (!turbPtr)
175  {
177  << "Unable to find a turbulence model."
178  << exit(FatalError);
179  }
180 
181  fieldNames_.resize(2);
182 
183  tmp<volScalarField> tepsilon = turbPtr->epsilon();
184  tmp<volScalarField> tomega = turbPtr->omega();
185 
186  if (tepsilon.is_reference())
187  {
188  isEpsilon_ = true;
189  fieldNames_[0] = tepsilon().name();
190  }
191  else if (tomega.is_reference())
192  {
193  isEpsilon_ = false;
194  fieldNames_[0] = tomega().name();
195  }
196  else
197  {
199  << "Unable to find an omega or epsilon field." << nl
200  << "buoyancyTurbSource needs an omega- or epsilon-based model."
201  << exit(FatalError);
202  }
203 
204  fieldNames_[1] = turbPtr->k()().name();
205 
207 
208  Log << " Applying buoyancyTurbSource to: "
209  << fieldNames_[0] << " and " << fieldNames_[1]
210  << endl;
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
217 (
218  fvMatrix<scalar>& eqn,
219  const label fieldi
220 )
221 {
222  if (fieldi == 1)
223  {
224  buoyancyTurbSourceK(eqn);
225  return;
226  }
227 
228  if (isEpsilon_)
229  {
230  buoyancyTurbSourceEpsilon(eqn);
231  }
232  else
233  {
234  buoyancyTurbSourceOmega(eqn);
235  }
236 }
237 
238 
240 (
241  const volScalarField& rho,
242  fvMatrix<scalar>& eqn,
243  const label fieldi
244 )
245 {
246  if (fieldi == 1)
247  {
248  buoyancyTurbSourceK(geometricOneField(), rho, eqn, fieldi);
249  return;
250  }
251 }
252 
253 
255 (
256  const volScalarField& alpha,
257  const volScalarField& rho,
258  fvMatrix<scalar>& eqn,
259  const label fieldi
260 )
261 {
262  if (fieldi == 1)
263  {
264  buoyancyTurbSourceK(alpha, rho, eqn, fieldi);
265  return;
266  }
267 }
268 
269 
270 // ************************************************************************* //
dimensionedScalar tanh(const dimensionedScalar &ds)
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Lookup and return non-const reference to the object of the given Type. Fatal if not found or the wron...
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
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to k and epsilon/omega equation for incompressible flow computations...
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
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
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
buoyancyTurbSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Vector< scalar > vector
Definition: vector.H:57
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
U
Definition: pEqn.H:72
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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: [].
static const Foam::dimensionedScalar B("", Foam::dimless, 18.678)
scalar nut
Namespace for OpenFOAM.
const dimensionSet dimVelocity