Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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\*---------------------------------------------------------------------------*/
29#include "phaseSystem.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<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
47template<class BasePhaseModel>
49{
50 BasePhaseModel::correctThermo();
51
52 // Correct the thermo, but make sure that the temperature remains the same
53 auto TCopy = volScalarField::New
54 (
55 IOobject::scopedName(this->thermo().T().name(), "Copy"),
57 this->thermo().T()
58 );
59 this->thermo_->he() = this->thermo().he(this->thermo().p(), TCopy());
60 this->thermo_->correct();
61 this->thermo_->T() = TCopy;
62}
63
64
65template<class BasePhaseModel>
68 return true;
69}
70
71
72template<class BasePhaseModel>
75{
77 << "Cannot construct an energy equation for an isothermal phase"
78 << exit(FatalError);
79
80 return nullptr;
81}
82
83
84// ************************************************************************* //
twoPhaseSystem & fluid
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field).
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
IsothermalPhaseModel(const phaseSystem &fluid, const word &phaseName, const label index)
virtual void correctThermo()
Correct the thermodynamics.
virtual bool isothermal() const
Return whether the phase is isothermal.
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
Class to represent a system of phases and model interfacial transfers between them.
Definition phaseSystem.H:72
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const volScalarField & T
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
psiReactionThermo & thermo