SolidThermoPhaseModel.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) 2023 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 "SolidThermoPhaseModel.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  viscosity(),
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>
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
63 template<class BasePhaseModel, class ThermoModel>
64 bool
66 {
67  return thermo_().incompressible();
68 }
69 
70 
71 template<class BasePhaseModel, class ThermoModel>
73 {
74  return thermo_().isochoric();
75 }
76 
77 
78 template<class BasePhaseModel, class ThermoModel>
79 const Foam::rhoThermo&
81 {
82  return thermo_();
83 }
84 
85 
86 template<class BasePhaseModel, class ThermoModel>
89 {
90  return thermo_();
91 }
92 
93 
94 template<class BasePhaseModel, class ThermoModel>
97 {
99  return refCast<const rhoFluidThermo>(thermo_());
100 }
101 
102 
103 template<class BasePhaseModel, class ThermoModel>
106 {
108  return refCast<rhoFluidThermo>(thermo_());
109 }
110 
111 
112 template<class BasePhaseModel, class ThermoModel>
115 {
116  return thermo_->rho();
117 }
118 
119 
120 template<class BasePhaseModel, class ThermoModel>
123 {
124  return thermo_->rho();
125 }
126 
127 
128 template<class BasePhaseModel, class ThermoModel>
131 {
133  return tmp<volScalarField>(nullptr);
134 }
135 
136 
137 template<class BasePhaseModel, class ThermoModel>
140 (
141  const label patchi
142 ) const
143 {
145  return tmp<scalarField>(nullptr);
146 }
147 
148 
149 template<class BasePhaseModel, class ThermoModel>
152 {
154  return tmp<volScalarField>(nullptr);
155 }
156 
157 
158 template<class BasePhaseModel, class ThermoModel>
161 (
162  const label patchi
163 ) const
164 {
166  return tmp<scalarField>(nullptr);
167 }
168 
169 
170 // ************************************************************************* //
static const char *const typeName
Definition: Field.H:106
Generic GeometricField class.
static word groupName(Name name, const word &group)
virtual const rhoThermo & thermo() const
Return the thermophysical model.
virtual const volScalarField & rho() const
Return the density field.
virtual ~SolidThermoPhaseModel()
Destructor.
virtual tmp< volScalarField > mu() const
Return the laminar dynamic viscosity.
virtual tmp< volScalarField > nu() const
Return the laminar kinematic viscosity.
SolidThermoPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual const rhoFluidThermo & fluidThermo() const
Return the thermophysical model.
autoPtr< ThermoModel > thermo_
Thermophysical model.
virtual bool isochoric() const
Return whether the phase is constant density.
virtual bool incompressible() const
Return whether the phase is incompressible.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
Base-class for fluid thermodynamic properties based on density.
Base-class for thermodynamic properties based on density.
Definition: rhoThermo.H:54
A class for managing temporary objects.
Definition: tmp.H:55
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
label patchi
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39