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-2022 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 :
47  (
48  typeName,
49  momentumTransport,
50  thermo
51  )
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
57 template<class BasicThermophysicalTransportModel>
58 const dictionary&
60 {
61  return dictionary::null;
62 }
63 
64 
65 template<class BasicThermophysicalTransportModel>
67 {
68  return true;
69 }
70 
71 
72 template<class TurbulenceThermophysicalTransportModel>
75 (
76  const volScalarField& Yi
77 ) const
78 {
80  << type() << " supports single component systems only, " << nl
81  << " for multi-component transport select"
82  " unityLewisFourier"
83  << exit(FatalError);
84 
85  return tmp<volScalarField>(nullptr);
86 }
87 
88 
89 template<class TurbulenceThermophysicalTransportModel>
92 (
93  const volScalarField& Yi,
94  const label patchi
95 ) const
96 {
98  << type() << " supports single component systems only, " << nl
99  << " for multi-component transport select"
100  " unityLewisFourier"
101  << exit(FatalError);
102 
103  return tmp<scalarField>(nullptr);
104 }
105 
106 
107 template<class BasicThermophysicalTransportModel>
109 {
110  const thermoModel& thermo = this->thermo();
111 
113  (
115  (
116  "q",
117  this->momentumTransport().alphaRhoPhi().group()
118  ),
119  -fvc::interpolate(this->alpha()*thermo.kappa())*fvc::snGrad(thermo.T())
120  );
121 }
122 
123 
124 template<class BasicThermophysicalTransportModel>
127 {
128  const thermoModel& thermo = this->thermo();
129 
130  // Return heat flux source as an implicit energy correction
131  // to the temperature gradient flux
132  return
133  -correction(fvm::laplacian(this->alpha()*thermo.alphahe(), he))
134  -fvc::laplacian(this->alpha()*thermo.kappa(), thermo.T());
135 }
136 
137 
138 template<class BasicThermophysicalTransportModel>
140 (
141  const volScalarField& Yi
142 ) const
143 {
145  << type() << " supports single component systems only, " << nl
146  << " for multi-component transport select"
147  " unityLewisFourier"
148  << exit(FatalError);
149 
150  return tmp<surfaceScalarField>(nullptr);
151 }
152 
153 
154 template<class BasicThermophysicalTransportModel>
157 {
159  << type() << " supports single component systems only, " << nl
160  << " for multi-component transport select"
161  " unityLewisFourier"
162  << exit(FatalError);
163 
164  return tmp<fvScalarMatrix>(nullptr);
165 }
166 
167 
168 template<class BasicThermophysicalTransportModel>
170 {
172  <
173  BasicThermophysicalTransportModel
174  >::correct();
175 }
176 
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 } // End namespace laminarThermophysicalTransportModels
181 } // End namespace Foam
182 
183 // ************************************************************************* //
Templated abstract base class for laminar thermophysical transport models.
const char *const group
Group name for atomic constants.
tmp< fvMatrix< Type > > correction(const fvMatrix< Type > &)
Return the correction form of the given matrix.
fluidReactionThermo & thermo
Definition: createFields.H:28
Fourier(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from components.
Definition: Fourier.C:41
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
Definition: Fourier.H:70
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.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.
Definition: Fourier.C:59
BasicThermophysicalTransportModel::thermoModel thermoModel
Definition: Fourier.H:73
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
Definition: Fourier.C:156
static word groupName(Name name, const word &group)
Calculate the laplacian of the given field.
virtual bool read()
Read thermophysicalTransport dictionary.
Definition: Fourier.C:66
static const dictionary null
Null dictionary.
Definition: dictionary.H:242
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:140
virtual void correct()
Correct the Fourier viscosity.
Definition: Fourier.C:169
static const char nl
Definition: Ostream.H:260
thermo he()
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.
label patchi
virtual tmp< volScalarField > DEff(const volScalarField &Yi) const
Effective mass diffusion coefficient.
Definition: Fourier.C:75
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
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Definition: Fourier.C:108
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
Definition: Fourier.C:126
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.