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-2021 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 #include "fvcLaplacian.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace turbulenceThermophysicalTransportModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class TurbulenceThermophysicalTransportModel>
41 {
42  alphat_ =
43  this->momentumTransport().rho()
44  *this->momentumTransport().nut()/Prt_;
45  alphat_.correctBoundaryConditions();
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 template<class TurbulenceThermophysicalTransportModel>
53 (
54  const momentumTransportModel& momentumTransport,
55  const thermoModel& thermo
56 )
57 :
58  TurbulenceThermophysicalTransportModel
59  (
60  typeName,
61  momentumTransport,
62  thermo
63  ),
64 
65  Prt_("Prt", dimless, this->coeffDict_),
66 
67  alphat_
68  (
69  IOobject
70  (
72  (
73  "alphat",
74  this->momentumTransport().alphaRhoPhi().group()
75  ),
76  momentumTransport.time().timeName(),
77  momentumTransport.mesh(),
80  ),
81  momentumTransport.mesh()
82  )
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
88 template<class TurbulenceThermophysicalTransportModel>
90 {
92  {
93  Prt_.read(this->coeffDict());
94 
95  return true;
96  }
97  else
98  {
99  return false;
100  }
101 }
102 
103 
104 template<class TurbulenceThermophysicalTransportModel>
107 (
108  const volScalarField& Yi
109 ) const
110 {
112  << type() << " supports single component systems only, " << nl
113  << " for multi-component transport select"
114  " nonUnityLewisEddyDiffusivity or unityLewisEddyDiffusivity"
115  << exit(FatalError);
116 
117  return tmp<volScalarField>(nullptr);
118 }
119 
120 
121 template<class TurbulenceThermophysicalTransportModel>
124 (
125  const volScalarField& Yi,
126  const label patchi
127 ) const
128 {
130  << type() << " supports single component systems only, " << nl
131  << " for multi-component transport select"
132  " nonUnityLewisEddyDiffusivity or unityLewisEddyDiffusivity"
133  << exit(FatalError);
134 
135  return tmp<scalarField>(nullptr);
136 }
137 
138 
139 template<class TurbulenceThermophysicalTransportModel>
142 {
144  (
146  (
147  "q",
148  this->momentumTransport().alphaRhoPhi().group()
149  ),
150  -fvc::interpolate(this->alpha()*this->kappaEff())
151  *fvc::snGrad(this->thermo().T())
152  );
153 }
154 
155 
156 template<class TurbulenceThermophysicalTransportModel>
159 (
160  volScalarField& he
161 ) const
162 {
163  // Return heat flux source as an implicit energy correction
164  // to the temperature gradient flux
165  return
166  -correction(fvm::laplacian(this->alpha()*this->alphaEff(), he))
167  -fvc::laplacian(this->alpha()*this->kappaEff(), this->thermo().T());
168 }
169 
170 
171 template<class TurbulenceThermophysicalTransportModel>
174 (
175  const volScalarField& Yi
176 ) const
177 {
179  << type() << " supports single component systems only, " << nl
180  << " for multi-component transport select"
181  " nonUnityLewisEddyDiffusivity or unityLewisEddyDiffusivity"
182  << exit(FatalError);
183 
184  return tmp<surfaceScalarField>(nullptr);
185 }
186 
187 
188 template<class TurbulenceThermophysicalTransportModel>
191 (
192  volScalarField& Yi
193 ) const
194 {
196  << type() << " supports single component systems only, " << nl
197  << " for multi-component transport select"
198  " nonUnityLewisEddyDiffusivity or unityLewisEddyDiffusivity"
199  << exit(FatalError);
200 
201  return tmp<fvScalarMatrix>(nullptr);
202 }
203 
204 
205 template<class TurbulenceThermophysicalTransportModel>
207 {
209  correctAlphat();
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace turbulenceThermophysicalTransportModels
216 } // End namespace Foam
217 
218 // ************************************************************************* //
const char *const group
Group name for atomic constants.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
fluidReactionThermo & thermo
Definition: createFields.H:28
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/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
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
Calculate the matrix for the laplacian of the field.
const dimensionSet dimless
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
static word groupName(Name name, const word &group)
Calculate the laplacian of the given field.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
virtual bool read()
Read thermophysicalTransport dictionary.
virtual void correct()
Correct the eddyDiffusivity viscosity.
static const char nl
Definition: Ostream.H:260
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
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:98
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
TurbulenceThermophysicalTransportModel::thermoModel thermoModel
Namespace for OpenFOAM.