ThermoPhaseModel.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-2021 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 "ThermoPhaseModel.H"
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel, class ThermoModel>
33 (
34  const phaseSystem& fluid,
35  const word& phaseName,
36  const bool referencePhase,
37  const label index
38 )
39 :
40  BasePhaseModel(fluid, phaseName, referencePhase, index),
41  dynamicTransportModel(),
42  thermo_(ThermoModel::New(fluid.mesh(), this->name()))
43 {
44  thermo_->validate
45  (
47  "h",
48  "e"
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
54 
55 template<class BasePhaseModel, class ThermoModel>
57 {}
58 
59 
60 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
61 
62 template<class BasePhaseModel, class ThermoModel>
64 {
65  return thermo_().incompressible();
66 }
67 
68 
69 template<class BasePhaseModel, class ThermoModel>
71 {
72  return thermo_().isochoric();
73 }
74 
75 
76 template<class BasePhaseModel, class ThermoModel>
77 const Foam::rhoThermo&
79 {
80  return thermo_();
81 }
82 
83 
84 template<class BasePhaseModel, class ThermoModel>
87 {
88  return thermo_();
89 }
90 
91 
92 template<class BasePhaseModel, class ThermoModel>
95 {
96  return thermo_->rho();
97 }
98 
99 
100 template<class BasePhaseModel, class ThermoModel>
103 {
104  return thermo_->mu();
105 }
106 
107 
108 template<class BasePhaseModel, class ThermoModel>
111 (
112  const label patchi
113 ) const
114 {
115  return thermo_->mu(patchi);
116 }
117 
118 
119 template<class BasePhaseModel, class ThermoModel>
122 {
123  return thermo_->nu();
124 }
125 
126 
127 template<class BasePhaseModel, class ThermoModel>
130 (
131  const label patchi
132 ) const
133 {
134  return thermo_->nu(patchi);
135 }
136 
137 
138 // ************************************************************************* //
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
static const char *const typeName
Definition: Field.H:105
virtual bool incompressible() const
Return whether the phase is incompressible.
virtual bool isochoric() const
Return whether the phase is constant density.
virtual rhoThermo & thermoRef()
Access the thermophysical model.
virtual const rhoThermo & thermo() const
Return the thermophysical model.
virtual tmp< volScalarField > mu() const
Return the laminar dynamic viscosity.
dynamicFvMesh & mesh
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
static word groupName(Name name, const word &group)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual tmp< volScalarField > rho() const
Return the density field.
Base-class for fluid thermodynamic properties based on density.
Definition: rhoThermo.H:52
ThermoPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
A class for managing temporary objects.
Definition: PtrList.H:53
virtual ~ThermoPhaseModel()
Destructor.
virtual tmp< volScalarField > nu() const
Return the laminar kinematic viscosity.