solidThermophysicalTransportModel.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) 2022-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 "isotropic.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 }
35 
36 
37 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
38 
40 (const word& type)
41 {
42  if (printCoeffs_)
43  {
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& type,
54  const alphaField& alpha,
55  const solidThermo& thermo
56 )
57 :
59  alpha_(alpha),
60  thermo_(thermo),
61  printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
62  coeffDict_(optionalSubDict(type + "Coeffs"))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
67 
70 {
72  (
73  solidThermophysicalTransportModel::typeName,
74  thermo.mesh().time().constant(),
75  thermo.mesh(),
78  false
79  );
80 
81  if (header.headerOk())
82  {
83  const word modelType(IOdictionary(header).lookup("model"));
84 
85  Info<< "Selecting solid thermophysical transport model "
86  << modelType << endl;
87 
88  typename dictionaryConstructorTable::iterator cstrIter =
89  dictionaryConstructorTablePtr_->find(modelType);
90 
91  if (cstrIter == dictionaryConstructorTablePtr_->end())
92  {
94  << "Unknown solid thermophysical transport model "
95  << modelType << nl << nl
96  << "Available models:" << endl
97  << dictionaryConstructorTablePtr_->sortedToc()
98  << exit(FatalError);
99  }
100 
102  (
103  cstrIter()(geometricOneField(), thermo)
104  );
105  }
106  else
107  {
108  Info<< "Selecting default solid thermophysical transport model "
109  << solidThermophysicalTransportModels::
110  isotropic<solidThermophysicalTransportModel>::typeName
111  << endl;
112 
114  (
115  new solidThermophysicalTransportModels::
116  isotropic<solidThermophysicalTransportModel>
117  (
119  thermo
120  )
121  );
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
127 
130 {
131  return thermo().kappa();
132 }
133 
134 
137 (
138  const label patchi
139 ) const
140 {
141  return thermo().kappa().boundaryField()[patchi];
142 }
143 
144 
146 {
147  if (regIOobject::read())
148  {
149  coeffDict_ <<= optionalSubDict(type() + "Coeffs");
150  return true;
151  }
152  else
153  {
154  return false;
155  }
156 }
157 
158 
160 {}
161 
162 
164 {}
165 
166 
167 // ************************************************************************* //
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const volScalarField & kappa() const =0
Thermal conductivity of mixture [W/m/K].
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
virtual bool read()
Read object.
Base-class for solid thermodynamic properties.
Definition: solidThermo.H:59
Abstract base class for solid thermophysical transport models.
dictionary coeffDict_
Model coefficients dictionary.
virtual void printCoeffs(const word &type)
Print model coefficients.
virtual void correct()
Solve the thermophysical transport model equations.
virtual bool read()=0
Read model coefficients if they have changed.
static autoPtr< solidThermophysicalTransportModel > New(const solidThermo &thermo)
Return a reference to the selected thermophysical transport model.
virtual void predict()=0
Predict the thermophysical transport coefficients if possible.
virtual tmp< volScalarField > kappa() const
Thermal conductivity [W/m/K].
solidThermophysicalTransportModel(const word &type, const alphaField &alpha, const solidThermo &thermo)
Construct from solid thermophysical properties.
Switch printCoeffs_
Flag to print the model coeffs at run-time.
Abstract base class for all fluid and solid thermophysical transport models.
A class for managing temporary objects.
Definition: tmp.H:55
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
error FatalError
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
fluidMulticomponentThermo & thermo
Definition: createFields.H:31