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-2020 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_->lookupOrDefault("pressureWorkAlphaLimit", 0.0);
42 
43  if (pressureWorkAlphaLimit > 0)
44  {
45  return
46  (
47  max(alpha - pressureWorkAlphaLimit, scalar(0))
48  /max(alpha - pressureWorkAlphaLimit, pressureWorkAlphaLimit)
49  )*pressureWork;
50  }
51  else
52  {
53  return pressureWork;
54  }
55 }
56 
57 
58 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
59 
60 template<class BasePhaseModel>
62 (
63  const phaseSystem& fluid,
64  const word& phaseName,
65  const bool referencePhase,
66  const label index
67 )
68 :
69  BasePhaseModel(fluid, phaseName, referencePhase, index)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
74 
75 template<class BasePhaseModel>
77 {}
78 
79 
80 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
81 
82 template<class BasePhaseModel>
84 {
86 
87  this->thermo_->correct();
88 }
89 
90 
91 template<class BasePhaseModel>
93 {
94  return false;
95 }
96 
97 
98 template<class BasePhaseModel>
101 {
102  const volScalarField& alpha = *this;
103 
104  const volVectorField U(this->U());
105  const surfaceScalarField alphaPhi(this->alphaPhi());
106  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
107 
108  const volScalarField contErr(this->continuityError());
109  const volScalarField K(this->K());
110 
111  volScalarField& he = this->thermo_->he();
112 
113  tmp<fvScalarMatrix> tEEqn
114  (
115  fvm::ddt(alpha, this->rho(), he)
116  + fvm::div(alphaRhoPhi, he)
117  - fvm::Sp(contErr, he)
118 
119  + fvc::ddt(alpha, this->rho(), K) + fvc::div(alphaRhoPhi, K)
120  - contErr*K
121  + this->divq(he)
122  ==
123  alpha*this->Qdot()
124  );
125 
126  // Add the appropriate pressure-work term
127  if (he.name() == this->thermo_->phasePropertyName("e"))
128  {
129  tEEqn.ref() += filterPressureWork
130  (
131  fvc::div(fvc::absolute(alphaPhi, alpha, U), this->thermo().p())
132  + (fvc::ddt(alpha) - contErr/this->rho())*this->thermo().p()
133  );
134  }
135  else if (this->thermo_->dpdt())
136  {
137  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
138  }
139 
140  return tEEqn;
141 }
142 
143 
144 // ************************************************************************* //
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
virtual bool isothermal() const
Return whether the phase is isothermal.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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 > &)
rhoReactionThermo & thermo
Definition: createFields.H:28
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::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
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
surfaceScalarField alphaPhi(phi.name()+alpha1.name(), fvc::flux(phi, alpha1, alphaScheme))
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.