IsothermalSolidPhaseModel.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 
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel>
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 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
46 template<class BasePhaseModel>
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
53 template<class BasePhaseModel>
55 {
56  BasePhaseModel::correctThermo();
57 }
58 
59 
60 template<class BasePhaseModel>
62 {
63  return true;
64 }
65 
66 
67 template<class BasePhaseModel>
70 (
71  const label patchi
72 ) const
73 {
75  return this->thermo().kappa().boundaryField()[patchi];
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 // ************************************************************************* //
const Boundary & boundaryField() const
Return const-reference to the boundary field.
virtual void correctThermo()
Correct the thermodynamics.
virtual bool isothermal() const
Return whether the phase is isothermal.
virtual tmp< fvScalarMatrix > heEqn()
Return the enthalpy equation.
IsothermalSolidPhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
virtual tmp< scalarField > kappaEff(const label patchi) const
Return the effective thermal conductivity on a patch.
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:73
A class for managing temporary objects.
Definition: tmp.H:55
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
error FatalError
fluidMulticomponentThermo & thermo
Definition: createFields.H:31