AnisothermalPhaseModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2015-2017 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 "AnisothermalPhaseModel.H"
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel>
33 (
34  const phaseSystem& fluid,
35  const word& phaseName,
36  const label index
37 )
38 :
39  BasePhaseModel(fluid, phaseName, index),
40  K_
41  (
42  IOobject
43  (
44  IOobject::groupName("K", this->name()),
45  fluid.mesh().time().timeName(),
46  fluid.mesh()
47  ),
48  fluid.mesh(),
49  dimensionedScalar("K", sqr(dimVelocity), scalar(0))
50  )
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
55 
56 template<class BasePhaseModel>
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
63 template<class BasePhaseModel>
65 {
66  return !this->thermo().incompressible();
67 }
68 
69 
70 template<class BasePhaseModel>
72 {
73  BasePhaseModel::correctKinematics();
74  K_ = 0.5*magSqr(this->U());
75 }
76 
77 
78 template<class BasePhaseModel>
80 {
82 
83  this->thermo_->correct();
84 }
85 
86 
87 template<class BasePhaseModel>
90 (
91  const tmp<volScalarField>& pressureWork
92 ) const
93 {
94  const volScalarField& alpha = *this;
95 
96  scalar pressureWorkAlphaLimit =
97  this->thermo_->lookupOrDefault("pressureWorkAlphaLimit", 0.0);
98 
99  if (pressureWorkAlphaLimit > 0)
100  {
101  return
102  (
103  max(alpha - pressureWorkAlphaLimit, scalar(0))
104  /max(alpha - pressureWorkAlphaLimit, pressureWorkAlphaLimit)
105  )*pressureWork;
106  }
107  else
108  {
109  return pressureWork;
110  }
111 }
112 
113 
114 template<class BasePhaseModel>
117 {
118  const volScalarField& alpha = *this;
119  const volVectorField& U = this->U();
120  const surfaceScalarField& alphaPhi = this->alphaPhi();
121  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi();
122 
123  const volScalarField& contErr(this->continuityError());
124 
125  const volScalarField alphaEff(this->turbulence().alphaEff());
126 
127  volScalarField& he = this->thermo_->he();
128 
129  tmp<fvScalarMatrix> tEEqn
130  (
131  fvm::ddt(alpha, this->rho(), he)
132  + fvm::div(alphaRhoPhi, he)
133  - fvm::Sp(contErr, he)
134 
135  + fvc::ddt(alpha, this->rho(), K_) + fvc::div(alphaRhoPhi, K_)
136  - contErr*K_
137 
139  (
140  fvc::interpolate(alpha)
142  he
143  )
144  ==
145  alpha*this->Qdot()
146  );
147 
148  // Add the appropriate pressure-work term
149  if (he.name() == this->thermo_->phasePropertyName("e"))
150  {
151  tEEqn.ref() += filterPressureWork
152  (
153  fvc::div(fvc::absolute(alphaPhi, alpha, U), this->thermo().p())
154  + this->thermo().p()*fvc::ddt(alpha)
155  );
156  }
157  else if (this->thermo_->dpdt())
158  {
159  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
160  }
161 
162  return tEEqn;
163 }
164 
165 
166 template<class BasePhaseModel>
169 {
170  return K_;
171 }
172 
173 
174 // ************************************************************************* //
virtual void correctThermo()
Correct the thermodynamics.
U
Definition: pEqn.H:83
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
multiphaseSystem & fluid
Definition: createFields.H:11
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Class which represents a phase for which the temperature (strictly energy) varies. Returns the energy equation and corrects the thermodynamic model.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
virtual const volScalarField & K() const
Return the phase kinetic energy.
scalar Qdot
Definition: solveChemistry.H:2
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
virtual void correctKinematics()
Correct the kinematics.
volScalarField & dpdt
word timeName
Definition: getTimeIndex.H:3
dimensioned< scalar > magSqr(const dimensioned< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual bool compressible() const
Return true if the phase is compressible otherwise false.
virtual ~AnisothermalPhaseModel()
Destructor.
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
mixture correctThermo()
const dimensionSet dimVelocity