IsothermalPhaseModel.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015-2018 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class BasePhaseModel>
35 (
36  const phaseSystem& fluid,
37  const word& phaseName,
38  const label index
39 )
40 :
41  BasePhaseModel(fluid, phaseName, index)
42 {}
43 
44 
45 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
46 
47 template<class BasePhaseModel>
49 {
51 
52  // Correct the thermo, but make sure that the temperature remains the same
54  (
56  (
57  this->thermo().T().name() + ":Copy",
58  this->thermo().T()
59  )
60  );
61  this->thermo_->he() = this->thermo().he(this->thermo().p(), TCopy);
62  this->thermo_->correct();
63  this->thermo_->T() = TCopy;
64 }
65 
66 
67 template<class BasePhaseModel>
69 {
70  return true;
71 }
72 
73 
74 template<class BasePhaseModel>
77 {
79  << "Cannot construct an energy equation for an isothermal phase"
80  << exit(FatalError);
81 
82  return nullptr;
83 }
84 
85 
86 // ************************************************************************* //
twoPhaseSystem & fluid
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
IsothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
virtual void correctThermo()
Correct the thermodynamics.
psiReactionThermo & thermo
Definition: createFields.H:28
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
A class for handling words, derived from Foam::string.
Definition: word.H:63
const volScalarField & T
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual bool isothermal() const
Return whether the phase is isothermal.
mixture correctThermo()