unityLewisFourier.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-2024 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 "unityLewisFourier.H"
27 #include "fvmLaplacian.H"
28 #include "fvcSnGrad.H"
29 #include "surfaceInterpolate.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace laminarThermophysicalTransportModels
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasicThermophysicalTransportModel>
42 (
43  const momentumTransportModel& momentumTransport,
44  const thermoModel& thermo
45 )
46 :
48  (
49  typeName,
50  momentumTransport,
51  thermo
52  )
53 {}
54 
55 
56 template<class BasicThermophysicalTransportModel>
58 (
59  const word& type,
60  const momentumTransportModel& momentumTransport,
61  const thermoModel& thermo
62 )
63 :
64  laminarThermophysicalTransportModel<BasicThermophysicalTransportModel>
65  (
66  type,
67  momentumTransport,
68  thermo
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class BasicThermophysicalTransportModel>
77 {
78  return true;
79 }
80 
81 
82 template<class BasicThermophysicalTransportModel>
85 {
87  (
89  (
90  "q",
91  this->momentumTransport().alphaRhoPhi().group()
92  ),
94  (
95  this->alpha()*this->thermo().kappa()/this->thermo().Cpv()
96  )
97  *fvc::snGrad(this->thermo().he())
98  );
99 }
100 
101 
102 template<class BasicThermophysicalTransportModel>
104 (
105  const label patchi
106 )
107 const
108 {
109  return
110  - (
111  this->alpha().boundaryField()[patchi]
112  *this->thermo().kappa().boundaryField()[patchi]
113  /this->thermo().Cpv()[patchi]
114  *this->thermo().he().boundaryField()[patchi].snGrad()
115  );
116 }
117 
118 
119 template<class BasicThermophysicalTransportModel>
122 divq(volScalarField& he) const
123 {
124  volScalarField alphahe
125  (
127  (
128  "alphahe",
129  this->thermo().kappa()/this->thermo().Cpv()
130  )
131  );
132 
133  return -fvm::laplacian(this->alpha()*alphahe, he);
134 }
135 
136 
137 template<class BasicThermophysicalTransportModel>
139 (
140  const volScalarField& Yi
141 ) const
142 {
144  (
146  (
147  "j(" + Yi.name() + ')',
148  this->momentumTransport().alphaRhoPhi().group()
149  ),
150  -fvc::interpolate(this->alpha()*this->DEff(Yi))
151  *fvc::snGrad(Yi)
152  );
153 }
154 
155 
156 template<class BasicThermophysicalTransportModel>
158 (
159  const volScalarField& Yi,
160  const label patchi
161 ) const
162 {
163  return
164  - (
165  this->alpha().boundaryField()[patchi]
166  *this->DEff(Yi, patchi)
167  *Yi.boundaryField()[patchi].snGrad()
168  );
169 }
170 
171 
172 template<class BasicThermophysicalTransportModel>
175 divj(volScalarField& Yi) const
176 {
177  return -fvm::laplacian(this->alpha()*this->DEff(Yi), Yi);
178 }
179 
180 
181 template<class BasicThermophysicalTransportModel>
183 {
185  <
186  BasicThermophysicalTransportModel
187  >::predict();
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace laminarThermophysicalTransportModels
194 } // End namespace Foam
195 
196 // ************************************************************************* //
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
const word & name() const
Return name.
Definition: IOobject.H:307
static word groupName(Name name, const word &group)
Templated abstract base class for laminar thermophysical transport models.
BasicThermophysicalTransportModel::thermoModel thermoModel
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
unityLewisFourier's energy gradient heat flux model for laminar flow. Specie fluxes are computed assu...
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
virtual void predict()
Correct the unityLewisFourier viscosity.
unityLewisFourier(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
virtual bool read()
Read thermophysicalTransport dictionary.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the snGrad of the given volField.
Calculate the matrix for the laplacian of the field.
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const char *const group
Group name for atomic constants.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
Namespace for OpenFOAM.
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:59
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()
fluidMulticomponentThermo & thermo
Definition: createFields.H:31