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-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 "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 :
65  (
66  type,
67  momentumTransport,
68  thermo
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 template<class BasicThermophysicalTransportModel>
76 const dictionary&
78 {
79  return dictionary::null;
80 }
81 
82 
83 template<class BasicThermophysicalTransportModel>
85 {
86  return true;
87 }
88 
89 
90 template<class BasicThermophysicalTransportModel>
93 {
95  (
97  (
98  "q",
99  this->momentumTransport().alphaRhoPhi().group()
100  ),
101  -fvc::interpolate(this->alpha()*this->thermo().alphahe())
102  *fvc::snGrad(this->thermo().he())
103  );
104 }
105 
106 
107 template<class BasicThermophysicalTransportModel>
111 {
112  return -fvm::laplacian(this->alpha()*this->thermo().alphahe(), he);
113 }
114 
115 
116 template<class BasicThermophysicalTransportModel>
118 (
119  const volScalarField& Yi
120 ) const
121 {
123  (
125  (
126  "j(" + Yi.name() + ')',
127  this->momentumTransport().alphaRhoPhi().group()
128  ),
129  -fvc::interpolate(this->alpha()*this->DEff(Yi))
130  *fvc::snGrad(Yi)
131  );
132 }
133 
134 
135 template<class BasicThermophysicalTransportModel>
139 {
140  return -fvm::laplacian(this->alpha()*this->DEff(Yi), Yi);
141 }
142 
143 
144 template<class BasicThermophysicalTransportModel>
146 {
148  <
149  BasicThermophysicalTransportModel
150  >::correct();
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
155 
156 } // End namespace laminarThermophysicalTransportModels
157 } // End namespace Foam
158 
159 // ************************************************************************* //
Templated abstract base class for laminar thermophysical transport models.
const char *const group
Group name for atomic constants.
const word & name() const
Return name.
Definition: IOobject.H:315
fluidReactionThermo & thermo
Definition: createFields.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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,.
virtual tmp< surfaceScalarField > q() const
Return the heat flux [W/m^2].
Calculate the matrix for the laplacian of the field.
Calculate the snGrad of the given volField.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const
Return the source term for the energy equation.
A class for handling words, derived from string.
Definition: word.H:59
unityLewisFourier(const momentumTransportModel &momentumTransport, const thermoModel &thermo)
Construct from a momentum transport model and a thermo model.
static word groupName(Name name, const word &group)
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].
virtual tmp< fvScalarMatrix > divj(volScalarField &Yi) const
Return the source term for the given specie mass-fraction equation.
thermo he()
virtual bool read()
Read thermophysicalTransport dictionary.
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.
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
unityLewisFourier&#39;s energy gradient heat flux model for laminar flow. Specie fluxes are computed assu...
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Correct the unityLewisFourier viscosity.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.
virtual const dictionary & coeffDict() const
Const access to the coefficients dictionary.