phaseSolidThermophysicalTransportModel.H
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-2024 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 Class
25  Foam::phaseSolidThermophysicalTransportModel
26 
27 Description
28  Abstract base class for solid thermophysical transport models
29 
30 SourceFiles
31  phaseSolidThermophysicalTransportModel.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef phaseSolidThermophysicalTransportModel_H
36 #define phaseSolidThermophysicalTransportModel_H
37 
39 #include "solidThermo.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class phaseSolidThermophysicalTransportModel Declaration
48 \*---------------------------------------------------------------------------*/
49 
51 :
53 {
54 public:
55 
56  typedef volScalarField alphaField;
57 
58 
59 protected:
60 
61  // Protected data
62 
63  const alphaField& alpha_;
64 
65  //- Reference to the solid thermophysical properties
66  const solidThermo& thermo_;
67 
68 
69  // Protected member functions
70 
71  //- Const access to the coefficients dictionary
72  const dictionary& coeffDict() const;
73 
74 
75 public:
76 
77  // Declare run-time constructor selection table
78 
80  (
81  autoPtr,
83  dictionary,
84  (
85  const alphaField& alpha,
86  const solidThermo& thermo
87  ),
88  (alpha, thermo)
89  );
90 
91 
92  // Constructors
93 
94  //- Construct from solid thermophysical properties
96  (
97  const word& type,
98  const alphaField& alpha,
99  const solidThermo& thermo
100  );
101 
102 
103  // Selectors
104 
105  //- Return a reference to the selected thermophysical transport model
107  (
108  const alphaField& alpha,
109  const solidThermo& thermo
110  );
111 
112 
113  //- Destructor
115  {}
116 
117 
118  // Member Functions
119 
120  //- Read model coefficients if they have changed
121  virtual bool read() = 0;
122 
123  //- Return the phase fraction field
124  const alphaField& alpha() const
125  {
126  return alpha_;
127  }
128 
129  //- Access function to solid thermophysical properties
130  virtual const solidThermo& thermo() const
131  {
132  return thermo_;
133  }
134 
135  //- Thermal conductivity [W/m/K]
136  virtual tmp<volScalarField> kappa() const;
137 
138  //- Thermal conductivity for patch [W/m/K]
139  virtual tmp<scalarField> kappa(const label patchi) const;
140 
141  //- Effective thermal conductivity
142  // of mixture [W/m/K]
143  virtual tmp<volScalarField> kappaEff() const
144  {
145  return kappa();
146  }
147 
148  //- Effective thermal conductivity
149  // of mixture for patch [W/m/K]
150  virtual tmp<scalarField> kappaEff(const label patchi) const
151  {
152  return kappa(patchi);
153  }
154 
155  //- Return the heat flux [W/m^2]
156  virtual tmp<surfaceScalarField> q() const = 0;
157 
158  //- Return the patch heat flux correction [W/m^2]
159  // For isotropic or patch-aligned thermal conductivity qCorr is null
160  virtual tmp<scalarField> qCorr(const label patchi) const = 0;
161 
162  //- Return the source term for the energy equation
163  virtual tmp<fvScalarMatrix> divq(volScalarField& he) const = 0;
164 
165  //- Predict the thermophysical transport coefficients if possible
166  // without solving thermophysical transport model equations
167  virtual void predict() = 0;
168 
169  //- Solve the thermophysical transport model equations
170  // and correct the thermophysical transport coefficients
171  virtual void correct();
172 };
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace Foam
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 #endif
182 
183 // ************************************************************************* //
Generic GeometricField class.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for solid thermophysical transport models.
virtual tmp< scalarField > qCorr(const label patchi) const =0
Return the patch heat flux correction [W/m^2].
declareRunTimeSelectionTable(autoPtr, phaseSolidThermophysicalTransportModel, dictionary,(const alphaField &alpha, const solidThermo &thermo),(alpha, thermo))
virtual void correct()
Solve the thermophysical transport model equations.
static autoPtr< phaseSolidThermophysicalTransportModel > New(const alphaField &alpha, const solidThermo &thermo)
Return a reference to the selected thermophysical transport model.
virtual tmp< fvScalarMatrix > divq(volScalarField &he) const =0
Return the source term for the energy equation.
virtual bool read()=0
Read model coefficients if they have changed.
phaseSolidThermophysicalTransportModel(const word &type, const alphaField &alpha, const solidThermo &thermo)
Construct from solid thermophysical properties.
const solidThermo & thermo_
Reference to the solid thermophysical properties.
const dictionary & coeffDict() const
Const access to the coefficients dictionary.
const alphaField & alpha() const
Return the phase fraction field.
virtual void predict()=0
Predict the thermophysical transport coefficients if possible.
virtual tmp< volScalarField > kappa() const
Thermal conductivity [W/m/K].
virtual tmp< surfaceScalarField > q() const =0
Return the heat flux [W/m^2].
virtual tmp< volScalarField > kappaEff() const
Effective thermal conductivity.
virtual const solidThermo & thermo() const
Access function to solid thermophysical properties.
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
Abstract base class for all fluid and solid thermophysical transport models.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
label patchi
Namespace for OpenFOAM.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
thermo he()