phaseModel.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) 2013-2016 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 
28 #include "phaseModel.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 (
34  const word& phaseName,
35  const volScalarField& p,
36  const volScalarField& T
37 )
38 :
40  (
41  IOobject
42  (
43  IOobject::groupName("alpha", phaseName),
44  p.mesh().time().timeName(),
45  p.mesh(),
46  IOobject::MUST_READ,
47  IOobject::AUTO_WRITE
48  ),
49  p.mesh()
50  ),
51  name_(phaseName),
52  p_(p),
53  T_(T),
54  thermo_(nullptr),
55  dgdt_
56  (
57  IOobject
58  (
59  IOobject::groupName("dgdt", phaseName),
60  p.mesh().time().timeName(),
61  p.mesh(),
62  IOobject::READ_IF_PRESENT,
63  IOobject::AUTO_WRITE
64  ),
65  p.mesh(),
67  )
68 {
69  {
70  volScalarField Tp(IOobject::groupName("T", phaseName), T);
71  Tp.write();
72  }
73 
74  thermo_ = rhoThermo::New(p.mesh(), phaseName);
75  thermo_->validate(phaseName, "e");
76 
77  correct();
78 }
79 
80 
81 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
82 
84 {
86  return nullptr;
87 }
88 
89 
91 {
92  thermo_->he() = thermo_->he(p_, T_);
93  thermo_->correct();
94 }
95 
96 
97 // ************************************************************************* //
autoPtr< phaseModel > clone() const
Return clone.
Definition: phaseModel.C:201
static autoPtr< rhoThermo > New(const fvMesh &, const word &phaseName=word::null)
Selector.
Definition: rhoThermo.C:190
const dimensionSet dimless
Dimensionless.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
void correct()
Correct the phase properties.
Definition: phaseModel.C:209
tmp< GeometricField< Type, PatchField, GeoMesh > > T() const
Return transpose (only if it is a tensor field)
const volScalarField & T
phaseModel(const word &phaseName, const dictionary &phaseDict, const fvMesh &mesh)
Definition: phaseModel.C:33
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
volScalarField & p
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127