eddyDiffusivity.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "eddyDiffusivity.H"
27 #include "fvmLaplacian.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace turbulenceThermophysicalTransportModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class TurbulenceThermophysicalTransportModel>
40 {
41  alphat_ =
42  this->momentumTransport().rho()
43  *this->momentumTransport().nut()/Prt_;
44  alphat_.correctBoundaryConditions();
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 template<class TurbulenceThermophysicalTransportModel>
52 (
53  const momentumTransportModel& momentumTransport,
54  const thermoModel& thermo
55 )
56 :
58  (
59  typeName,
60  momentumTransport,
61  thermo,
62  false
63  )
64 {}
65 
66 
67 template<class TurbulenceThermophysicalTransportModel>
69 (
70  const word& type,
71  const momentumTransportModel& momentumTransport,
72  const thermoModel& thermo,
73  const bool allowDefaultPrt
74 )
75 :
76  TurbulenceThermophysicalTransportModel
77  (
78  typeName,
79  momentumTransport,
80  thermo
81  ),
82 
83  Prt_
84  (
85  allowDefaultPrt
87  (
88  "Prt",
89  this->coeffDict_,
90  1
91  )
93  (
94  "Prt",
95  dimless,
96  this->coeffDict_
97  )
98  ),
99 
100  alphat_
101  (
102  IOobject
103  (
105  (
106  "alphat",
107  this->momentumTransport().alphaRhoPhi().group()
108  ),
109  momentumTransport.time().timeName(),
110  momentumTransport.mesh(),
113  ),
114  momentumTransport.mesh()
115  )
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
121 template<class TurbulenceThermophysicalTransportModel>
123 {
125  {
126  Prt_.readIfPresent(this->coeffDict());
127 
128  return true;
129  }
130  else
131  {
132  return false;
133  }
134 }
135 
136 
137 template<class TurbulenceThermophysicalTransportModel>
140 {
141  return volVectorField::New
142  (
144  (
145  "q",
146  this->momentumTransport().alphaRhoPhi().group()
147  ),
148  -this->alphaEff()*this->alpha()*fvc::grad(this->thermo().he())
149  );
150 }
151 
152 
153 template<class TurbulenceThermophysicalTransportModel>
156 (
157  volScalarField& he
158 ) const
159 {
160  return -fvm::laplacian(this->alpha()*this->alphaEff(), he);
161 }
162 
163 
164 template<class TurbulenceThermophysicalTransportModel>
167 (
168  const volScalarField& Yi
169 ) const
170 {
171  return volVectorField::New
172  (
174  (
175  "j(" + Yi.name() + ')',
176  this->momentumTransport().alphaRhoPhi().group()
177  ),
178  -this->DEff(Yi)*this->alpha()*fvc::grad(Yi)
179  );
180 }
181 
182 
183 template<class TurbulenceThermophysicalTransportModel>
186 (
187  volScalarField& Yi
188 ) const
189 {
190  return -fvm::laplacian(this->alpha()*this->DEff(Yi), Yi);
191 }
192 
193 
194 template<class TurbulenceThermophysicalTransportModel>
196 {
198  correctAlphat();
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
203 
204 } // End namespace turbulenceThermophysicalTransportModels
205 } // End namespace Foam
206 
207 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
const char *const group
Group name for atomic constants.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
const word & name() const
Return name.
Definition: IOobject.H:303
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
Calculate the matrix for the laplacian of the field.
Eddy-diffusivity based gradient heat flux model for RAS or LES of turbulent flow. Specie fluxes are c...
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
rhoReactionThermo & thermo
Definition: createFields.H:28
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual bool read()
Read thermophysicalTransport dictionary.
virtual void correct()
Correct the eddyDiffusivity viscosity.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
virtual tmp< volVectorField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
eddyDiffusivity(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual tmp< volVectorField > q() const
Return the heat flux.
TurbulenceThermophysicalTransportModel::thermoModel thermoModel
Namespace for OpenFOAM.