Fourier.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 "Fourier.H"
27 #include "fvmLaplacian.H"
28 #include "fvcLaplacian.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace laminarThermophysicalTransportModels
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class BasicThermophysicalTransportModel>
41 (
42  const momentumTransportModel& momentumTransport,
43  const thermoModel& thermo
44 )
45 :
46  laminarThermophysicalTransportModel<BasicThermophysicalTransportModel>
47  (
48  typeName,
49  momentumTransport,
50  thermo
51  )
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
57 template<class BasicThermophysicalTransportModel>
59 {
60  return true;
61 }
62 
63 
64 template<class TurbulenceThermophysicalTransportModel>
67 (
68  const volScalarField& Yi
69 ) const
70 {
72  << type() << " supports single component systems only, " << nl
73  << " for multi-component transport select"
74  " unityLewisFourier"
75  << exit(FatalError);
76 
77  return tmp<volScalarField>(nullptr);
78 }
79 
80 
81 template<class TurbulenceThermophysicalTransportModel>
84 (
85  const volScalarField& Yi,
86  const label patchi
87 ) const
88 {
90  << type() << " supports single component systems only, " << nl
91  << " for multi-component transport select"
92  " unityLewisFourier"
93  << exit(FatalError);
94 
95  return tmp<scalarField>(nullptr);
96 }
97 
98 
99 template<class BasicThermophysicalTransportModel>
101 {
102  const thermoModel& thermo = this->thermo();
103 
105  (
107  (
108  "q",
109  this->momentumTransport().alphaRhoPhi().group()
110  ),
111  -fvc::interpolate(this->alpha()*thermo.kappa())*fvc::snGrad(thermo.T())
112  );
113 }
114 
115 
116 template<class BasicThermophysicalTransportModel>
118 (
119  const label patchi
120 ) const
121 {
122  return
123  - (
124  this->alpha().boundaryField()[patchi]
125  *this->thermo().kappa().boundaryField()[patchi]
126  *this->thermo().T().boundaryField()[patchi].snGrad()
127  );
128 }
129 
130 
131 template<class BasicThermophysicalTransportModel>
134 {
135  const thermoModel& thermo = this->thermo();
136 
137  // Return heat flux source as an implicit energy correction
138  // to the temperature gradient flux
139  return
140  -fvc::laplacian(this->alpha()*thermo.kappa(), thermo.T())
142  (
143  this->alpha()*thermo.kappa()/thermo.Cpv(),
144  he
145  );
146 }
147 
148 
149 template<class BasicThermophysicalTransportModel>
151 (
152  const volScalarField& Yi
153 ) const
154 {
156  << type() << " supports single component systems only, " << nl
157  << " for multi-component transport select"
158  " unityLewisFourier"
159  << exit(FatalError);
160 
161  return tmp<surfaceScalarField>(nullptr);
162 }
163 
164 
165 template<class BasicThermophysicalTransportModel>
167 (
168  const volScalarField& Yi,
169  const label patchi
170 ) const
171 {
173  << type() << " supports single component systems only, " << nl
174  << " for multi-component transport select"
175  " unityLewisFourier"
176  << exit(FatalError);
177 
178  return tmp<scalarField>(nullptr);
179 }
180 
181 
182 template<class BasicThermophysicalTransportModel>
185 {
187  << type() << " supports single component systems only, " << nl
188  << " for multi-component transport select"
189  " unityLewisFourier"
190  << exit(FatalError);
191 
192  return tmp<fvScalarMatrix>(nullptr);
193 }
194 
195 
196 template<class BasicThermophysicalTransportModel>
198 {
200  <
201  BasicThermophysicalTransportModel
202  >::predict();
203 }
204 
205 
206 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
207 
208 } // End namespace laminarThermophysicalTransportModels
209 } // End namespace Foam
210 
211 // ************************************************************************* //
Generic GeometricField class.
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,.
static word groupName(Name name, const word &group)
Templated abstract base class for laminar thermophysical transport models.
Fourier(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from components.
Definition: Fourier.C:41
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fourier.C:67
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fourier.C:133
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fourier.C:184
virtual void predict()
Correct the Fourier viscosity.
Definition: Fourier.C:197
virtual tmp< surfaceScalarField > j(const volScalarField &Yi) const
Return the specie flux for the given specie mass-fraction [kg/m^2/s].
Definition: Fourier.C:151
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fourier.C:100
BasicThermophysicalTransportModel::thermoModel thermoModel
Definition: Fourier.H:73
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
Definition: Fourier.H:70
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: Fourier.C:58
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the laplacian of the given field.
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.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > laplacianCorrection(const VolField< scalar > &gamma, const VolField< Type > &vf)
Definition: fvmLaplacian.C:340
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
error FatalError
static const char nl
Definition: Ostream.H:267
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