AnisothermalPhaseModel.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) 2015-2021 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 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel>
34 (
35  const tmp<volScalarField>& pressureWork
36 ) const
37 {
38  const volScalarField& alpha = *this;
39 
40  scalar pressureWorkAlphaLimit =
41  this->thermo_->properties()
42  .lookupOrDefault("pressureWorkAlphaLimit", 0.0);
43 
44  if (pressureWorkAlphaLimit > 0)
45  {
46  return
47  (
48  max(alpha - pressureWorkAlphaLimit, scalar(0))
49  /max(alpha - pressureWorkAlphaLimit, pressureWorkAlphaLimit)
50  )*pressureWork;
51  }
52  else
53  {
54  return pressureWork;
55  }
56 }
57 
58 
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
60 
61 template<class BasePhaseModel>
63 (
64  const phaseSystem& fluid,
65  const word& phaseName,
66  const bool referencePhase,
67  const label index
68 )
69 :
70  BasePhaseModel(fluid, phaseName, referencePhase, index)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
75 
76 template<class BasePhaseModel>
78 {}
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
83 template<class BasePhaseModel>
85 {
87 
88  this->thermo_->correct();
89 }
90 
91 
92 template<class BasePhaseModel>
94 {
95  return false;
96 }
97 
98 
99 template<class BasePhaseModel>
102 {
103  const volScalarField& alpha = *this;
104  const volScalarField& rho = this->rho();
105 
106  const tmp<volVectorField> tU(this->U());
107  const volVectorField& U(tU());
108 
109  const tmp<surfaceScalarField> talphaRhoPhi(this->alphaRhoPhi());
110  const surfaceScalarField& alphaRhoPhi(talphaRhoPhi());
111 
112  const tmp<volScalarField> tcontErr(this->continuityError());
113  const volScalarField& contErr(tcontErr());
114 
115  tmp<volScalarField> tK(this->K());
116  const volScalarField& K(tK());
117 
118  volScalarField& he = this->thermo_->he();
119 
120  tmp<fvScalarMatrix> tEEqn
121  (
122  fvm::ddt(alpha, rho, he)
123  + fvm::div(alphaRhoPhi, he)
124  - fvm::Sp(contErr, he)
125 
126  + fvc::ddt(alpha, rho, K) + fvc::div(alphaRhoPhi, K)
127  - contErr*K
128 
129  + this->divq(he)
130  ==
131  alpha*this->Qdot()
132  );
133 
134  // Add the appropriate pressure-work term
135  if (he.name() == this->thermo_->phasePropertyName("e"))
136  {
137  tEEqn.ref() += filterPressureWork
138  (
139  fvc::div
140  (
141  fvc::absolute(alphaRhoPhi, alpha, rho, U),
142  this->thermo().p()/rho
143  )
144  + (fvc::ddt(alpha) - contErr/rho)*this->thermo().p()
145  );
146  }
147  else if (this->thermo_->dpdt())
148  {
149  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
150  }
151 
152  return tEEqn;
153 }
154 
155 
156 // ************************************************************************* //
virtual void correctThermo()
Correct the thermodynamics.
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
fluidReactionThermo & thermo
Definition: createFields.H:28
virtual bool isothermal() const
Return whether the phase is isothermal.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvModels.source(alpha1, mixture.thermo1().rho())&rho1) -(fvModels.source(alpha2, mixture.thermo2().rho())&rho2))())
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:58
CGAL::Exact_predicates_exact_constructions_kernel K
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:57
volScalarField & dpdt
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
mixture correctThermo()
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phaseSystem & fluid
Definition: createFields.H:11
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
virtual ~AnisothermalPhaseModel()
Destructor.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.