IsothermalPhaseModel.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 "IsothermalPhaseModel.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 {}
41 
42 
43 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
44 
45 template<class BasePhaseModel>
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
52 template<class BasePhaseModel>
54 {
56 
57  // Correct the thermo, but make sure that the temperature remains the same
58  tmp<volScalarField> TCopy
59  (
61  (
62  this->thermo().T().name() + ":Copy",
63  this->thermo().T()
64  )
65  );
66  this->thermo_->he() = this->thermo().he(this->thermo().p(), TCopy);
67  this->thermo_->correct();
68  this->thermo_->T() = TCopy;
69 }
70 
71 
72 template<class BasePhaseModel>
74 {
75  return true;
76 }
77 
78 
79 template<class BasePhaseModel>
82 {
84  << "Cannot construct an energy equation for an isothermal phase"
85  << exit(FatalError);
86 
87  return tmp<fvScalarMatrix>();
88 }
89 
90 
91 // ************************************************************************* //
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual bool isothermal() const
Return whether the phase is isothermal.
IsothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual void correctThermo()
Correct the thermodynamics.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
mixture correctThermo()
virtual ~IsothermalPhaseModel()
Destructor.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53