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-2018 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 label index
66 )
67 :
68  BasePhaseModel(fluid, phaseName, index)
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
73 
74 template<class BasePhaseModel>
76 {}
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
81 template<class BasePhaseModel>
83 {
85 
86  this->thermo_->correct();
87 }
88 
89 
90 template<class BasePhaseModel>
92 {
93  return false;
94 }
95 
96 
97 template<class BasePhaseModel>
100 {
101  const volScalarField& alpha = *this;
102 
103  const volVectorField U(this->U());
104  const surfaceScalarField alphaPhi(this->alphaPhi());
105  const surfaceScalarField alphaRhoPhi(this->alphaRhoPhi());
106 
107  const volScalarField contErr(this->continuityError());
108  const volScalarField K(this->K());
109 
110  volScalarField& he = this->thermo_->he();
111 
112  tmp<fvScalarMatrix> tEEqn
113  (
114  fvm::ddt(alpha, this->rho(), he)
115  + fvm::div(alphaRhoPhi, he)
116  - fvm::Sp(contErr, he)
117 
118  + fvc::ddt(alpha, this->rho(), K) + fvc::div(alphaRhoPhi, K)
119  - contErr*K
120 
122  (
123  fvc::interpolate(alpha)
124  *fvc::interpolate(this->alphaEff()),
125  he
126  )
127  ==
128  alpha*this->Qdot()
129  );
130 
131  // Add the appropriate pressure-work term
132  if (he.name() == this->thermo_->phasePropertyName("e"))
133  {
134  tEEqn.ref() += filterPressureWork
135  (
136  fvc::div(fvc::absolute(alphaPhi, alpha, U), this->thermo().p())
137  + this->thermo().p()*fvc::ddt(alpha)
138  );
139  }
140  else if (this->thermo_->dpdt())
141  {
142  tEEqn.ref() -= filterPressureWork(alpha*this->fluid().dpdt());
143  }
144 
145  return tEEqn;
146 }
147 
148 
149 // ************************************************************************* //
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
multiphaseSystem & fluid
Definition: createFields.H:11
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:55
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:52
volScalarField::Internal contErr((fvc::ddt(rho)+fvc::div(rhoPhi) -(fvOptions(alpha1, mixture.thermo1().rho())&rho1) -(fvOptions(alpha2, mixture.thermo2().rho())&rho2))())
volScalarField & dpdt
mixture correctThermo()
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
AnisothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
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
U
Definition: pEqn.H:72
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual ~AnisothermalPhaseModel()
Destructor.
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.