multiphaseStabilizedTurbulence.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) 2019 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 
29 #include "fvMatrices.H"
31 #include "gravityMeshObject.H"
33 
34 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace fv
39 {
40  defineTypeNameAndDebug(multiphaseStabilizedTurbulence, 0);
41 
43  (
44  option,
45  multiphaseStabilizedTurbulence,
46  dictionary
47  );
48 }
49 }
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
53 Foam::fv::multiphaseStabilizedTurbulence::multiphaseStabilizedTurbulence
54 (
55  const word& sourceName,
56  const word& modelType,
57  const dictionary& dict,
58  const fvMesh& mesh
59 )
60 :
61  option(sourceName, modelType, dict, mesh),
62  rhoName_(coeffs_.getOrDefault<word>("rho", "rho")),
63  Cmu_
64  (
65  dimensionedScalar::lookupOrAddToDict
66  (
67  "Cmu",
68  coeffs_,
69  0.09
70  )
71  ),
72  C_
73  (
74  dimensionedScalar::lookupOrAddToDict
75  (
76  "C",
77  coeffs_,
78  1.51
79  )
80  ),
81  lambda2_
82  (
83  dimensionedScalar::lookupOrAddToDict
84  (
85  "lambda2",
86  coeffs_,
87  0.1
88  )
89  ),
90  alpha_
91  (
92  dimensionedScalar::lookupOrAddToDict
93  (
94  "alpha",
95  coeffs_,
96  1.36
97  )
98  )
99 {
100  fieldNames_.setSize(2, "undefined");
101 
102  // Note: incompressible only
103  const auto* turbPtr =
104  mesh_.findObject<incompressible::turbulenceModel>
105  (
106  turbulenceModel::propertiesName
107  );
108 
109  if (turbPtr)
110  {
111  const tmp<volScalarField>& tk = turbPtr->k();
112  fieldNames_[0] = tk().name();
113 
114  const tmp<volScalarField>& tnut = turbPtr->nut();
115  fieldNames_[1] = tnut().name();
116 
117  Log << " Applying model to " << fieldNames_[0]
118  << " and " << fieldNames_[1] << endl;
119  }
120  else
121  {
123  << "Unable to find incompressible turbulence model"
124  << exit(FatalError);
125  }
126 
127  applied_.setSize(fieldNames_.size(), false);
128 }
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
134 (
135  const volScalarField& rho,
136  fvMatrix<scalar>& eqn,
137  const label fieldi
138 )
139 {
140  // Not applicable to compressible cases
142 }
143 
144 
146 (
147  fvMatrix<scalar>& eqn,
148  const label fieldi
149 )
150 {
151  if (fieldi != 0)
152  {
153  return;
154  }
155 
156  Log << this->name() << ": applying bouyancy production term to "
157  << eqn.psi().name() << endl;
158 
159  // Buoyancy production in k eqn
160 
161  const auto* turbPtr =
162  mesh_.findObject<incompressible::turbulenceModel>
163  (
164  turbulenceModel::propertiesName
165  );
166 
167  if (!turbPtr)
168  {
170  << "Unable to find incompressible turbulence model"
171  << exit(FatalError);
172  }
173 
174 
175  tmp<volScalarField> tepsilon = turbPtr->epsilon();
176  const volScalarField& epsilon = tepsilon();
177  const volScalarField& k = eqn.psi();
178 
179  // Note: using solver density field for incompressible multiphase cases
180  const auto& rho = mesh_.lookupObject<volScalarField>(rhoName_);
181 
182  const auto& g = meshObjects::gravity::New(mesh_.time());
183 
184  const dimensionedScalar eps0("eps0", epsilon.dimensions(), SMALL);
185 
186  // Note: differing from reference by replacing nut/k by Cmu*k/epsilon
187  const volScalarField GbyK
188  (
189  "GbyK",
190  Cmu_*k/(epsilon + eps0)*alpha_*(g & fvc::grad(rho))/rho
191  );
192 
193  return eqn -= fvm::SuSp(GbyK, k);
194 }
195 
196 
198 {
199  if (field.name() != fieldNames_[1])
200  {
201  return;
202  }
203 
204  Log << this->name() << ": correcting " << field.name() << endl;
205 
206  const auto* turbPtr =
208  (
210  );
211 
212  // nut correction
213  const auto& U = turbPtr->U();
214  tmp<volScalarField> tepsilon = turbPtr->epsilon();
215  const auto& epsilon = tepsilon();
216  tmp<volScalarField> tk = turbPtr->k();
217  const auto& k = tk();
218 
219  tmp<volTensorField> tgradU = fvc::grad(U);
220  const auto& gradU = tgradU();
221  const dimensionedScalar pSmall("pSmall", dimless/sqr(dimTime), SMALL);
222  const volScalarField pRat
223  (
224  magSqr(symm(gradU))/(magSqr(skew(gradU)) + pSmall)
225  );
226 
227  const volScalarField epsilonTilde
228  (
229  max
230  (
231  epsilon,
232  lambda2_*C_*pRat*epsilon
233  )
234  );
235 
236  field = Cmu_*sqr(k)/epsilonTilde;
237  field.correctBoundaryConditions();
238 }
239 
240 
241 // ************************************************************************* //
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:84
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
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:62
Foam::skew
dimensionedTensor skew(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:137
Foam::tmp< volScalarField >
Foam::fv::option::name
const word & name() const
Return const access to the source name.
Definition: fvOptionI.H:30
turbulentTransportModel.H
gravityMeshObject.H
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:100
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
rho
rho
Definition: readInitialConditions.H:96
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:285
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:82
Foam::fam::SuSp
tmp< faMatrix< Type > > SuSp(const areaScalarField &sp, const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition: famSup.C:151
Foam::magSqr
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:419
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::fv::option
Finite volume options abstract base class. Provides a base set of controls, e.g.:
Definition: fvOption.H:69
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fv::option::fieldNames_
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:94
field
rDeltaTY field()
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
dict
dictionary dict
Definition: searchingEngine.H:14
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:121
Foam::fv::multiphaseStabilizedTurbulence::correct
virtual void correct(volScalarField &field)
Correct the turbulence viscosity.
Definition: multiphaseStabilizedTurbulence.C:197
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:84
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:24
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
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
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
U
U
Definition: pEqn.H:72
Foam::fv::multiphaseStabilizedTurbulence::addSup
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible k equation.
Definition: multiphaseStabilizedTurbulence.C:134
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::objectRegistry::findObject
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
Definition: objectRegistryTemplates.C:401
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
Foam::GeometricField< scalar, fvPatchField, volMesh >
multiphaseStabilizedTurbulence.H
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332