nonUnityLewisEddyDiffusivity.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 
27 #include "fvcLaplacian.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace turbulenceThermophysicalTransportModels
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class TurbulenceThermophysicalTransportModel>
41 (
42  const momentumTransportModel& momentumTransport,
43  const thermoModel& thermo
44 )
45 :
47  (
48  typeName,
49  momentumTransport,
50  thermo,
51  false
52  ),
53 
54  Sct_
55  (
57  (
58  "Sct",
59  dimless,
60  this->coeffDict_
61  )
62  )
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
68 template<class TurbulenceThermophysicalTransportModel>
69 bool
71 {
73  {
74  Sct_.readIfPresent(this->coeffDict());
75 
76  return true;
77  }
78  else
79  {
80  return false;
81  }
82 }
83 
84 
85 template<class TurbulenceThermophysicalTransportModel>
88 {
89  tmp<volVectorField> tmpq =
91 
92  if (mag(this->Prt_ - Sct_).value() > small)
93  {
94  const basicSpecieMixture& composition = this->thermo().composition();
95 
96  const PtrList<volScalarField>& Y = composition.Y();
97 
98  volScalarField alphatMinusDt
99  (
100  "alphatMinusDt",
101  this->alphat()*(1 - this->Prt_/Sct_)
102  );
103 
104  forAll(Y, i)
105  {
106  tmpq.ref() +=
107  this->alpha()
108  *alphatMinusDt
109  *composition.HE(i, this->thermo().p(), this->thermo().T())
110  *fvc::grad(Y[i]);
111  }
112  }
113 
114  return tmpq;
115 }
116 
117 
118 template<class TurbulenceThermophysicalTransportModel>
121 (
123 ) const
124 {
125  tmp<fvScalarMatrix> tmpDivq =
127 
128  if (mag(this->Prt_ - Sct_).value() > small)
129  {
130  const basicSpecieMixture& composition = this->thermo().composition();
131 
132  const PtrList<volScalarField>& Y = composition.Y();
133 
134  volScalarField alphatMinusDt
135  (
136  "alphatMinusDt",
137  this->alphat()*(1 - this->Prt_/Sct_)
138  );
139 
140  forAll(Y, i)
141  {
142  tmpDivq.ref() +=
144  (
145  this->alpha()
146  *alphatMinusDt
147  *composition.HE(i, this->thermo().p(), this->thermo().T()),
148  Y[i]
149  );
150  }
151  }
152 
153  return tmpDivq;
154 }
155 
156 
157 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158 
159 } // End namespace turbulenceThermophysicalTransportModels
160 } // End namespace Foam
161 
162 // ************************************************************************* //
volScalarField & he
Definition: YEEqn.H:50
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
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
basicSpecieMixture & composition
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Eddy-diffusivity based gradient heat flux model for RAS or LES of turbulent flow. Specie fluxes are c...
rhoReactionThermo & thermo
Definition: createFields.H:28
Specialization of basicMixture for a mixture consisting of a number for molecular species...
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
nonUnityLewisEddyDiffusivity(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
Calculate the laplacian of the given field.
virtual scalar HE(const label speciei, const scalar p, const scalar T) const =0
Enthalpy/Internal energy [J/kg].
TurbulenceThermophysicalTransportModel::momentumTransportModel momentumTransportModel
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
PtrList< volScalarField > & Y
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
virtual tmp< volVectorField > q() const
Return the heat flux.
TurbulenceThermophysicalTransportModel::thermoModel thermoModel
Namespace for OpenFOAM.