laminarThermophysicalTransportModel.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) 2020-2022 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 "unityLewisFourier.H"
28 
29 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
30 
31 template<class BasicThermophysicalTransportModel>
33 <
34  BasicThermophysicalTransportModel
35 >::printCoeffs(const word& type)
36 {
37  if (printCoeffs_)
38  {
39  Info<< coeffDict_.dictName() << coeffDict_ << endl;
40  }
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
46 template<class BasicThermophysicalTransportModel>
48 <
49  BasicThermophysicalTransportModel
51 (
52  const word& type,
53  const momentumTransportModel& momentumTransport,
54  const thermoModel& thermo
55 )
56 :
57  BasicThermophysicalTransportModel(momentumTransport, thermo),
58  laminarDict_(this->subOrEmptyDict("laminar")),
59  printCoeffs_(laminarDict_.lookupOrDefault<Switch>("printCoeffs", false)),
60  coeffDict_(laminarDict_.optionalSubDict(type + "Coeffs"))
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
65 
66 template<class BasicThermophysicalTransportModel>
68 <
70  <
71  BasicThermophysicalTransportModel
72  >
73 >
75 <
76  BasicThermophysicalTransportModel
77 >::New
78 (
79  const momentumTransportModel& momentumTransport,
80  const thermoModel& thermo
81 )
82 {
84  (
86  (
87  thermophysicalTransportModel::typeName,
88  momentumTransport.alphaRhoPhi().group()
89  ),
90  momentumTransport.time().constant(),
91  momentumTransport.mesh(),
94  false
95  );
96 
97  if (header.headerOk())
98  {
99  IOdictionary modelDict(header);
100 
101  const word modelType(modelDict.subDict("laminar").lookup( "model"));
102 
103  Info<< "Selecting laminar thermophysical transport model "
104  << modelType << endl;
105 
106  typename dictionaryConstructorTable::iterator cstrIter =
107  dictionaryConstructorTablePtr_->find(modelType);
108 
109  if (cstrIter == dictionaryConstructorTablePtr_->end())
110  {
112  << "Unknown laminar thermophysical transport model "
113  << modelType << nl << nl
114  << "Available models:" << endl
115  << dictionaryConstructorTablePtr_->sortedToc()
116  << exit(FatalError);
117  }
118 
119  return autoPtr<laminarThermophysicalTransportModel>
120  (
121  cstrIter()(momentumTransport, thermo)
122  );
123  }
124  else
125  {
126  Info<< "Selecting default laminar thermophysical transport model "
127  << laminarThermophysicalTransportModels::unityLewisFourier<
128  BasicThermophysicalTransportModel>::typeName << endl;
129 
130  return autoPtr<laminarThermophysicalTransportModel>
131  (
132  new laminarThermophysicalTransportModels::unityLewisFourier
133  <
134  BasicThermophysicalTransportModel
135  >(momentumTransport, thermo)
136  );
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
143 template<class BasicThermophysicalTransportModel>
145 <
146  BasicThermophysicalTransportModel
148 {
150  {
151  laminarDict_ <<= this->subDict("laminar");
152 
153  coeffDict_ <<= laminarDict_.optionalSubDict(type() + "Coeffs");
154 
155  return true;
156  }
157  else
158  {
159  return false;
160  }
161 }
162 
163 
164 template<class BasicThermophysicalTransportModel>
166 <
167  BasicThermophysicalTransportModel
168 >::predict()
169 {
170  BasicThermophysicalTransportModel::predict();
171 }
172 
173 
174 template<class BasicThermophysicalTransportModel>
176 <
177  BasicThermophysicalTransportModel
179 {
181 }
182 
183 
184 // ************************************************************************* //
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:134
static word groupName(Name name, const word &group)
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Templated abstract base class for laminar thermophysical transport models.
BasicThermophysicalTransportModel::thermoModel thermoModel
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
Abstract base class for turbulence models (RAS, LES and laminar).
const surfaceScalarField & alphaRhoPhi() const
Access function to phase flux field.
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
bool read(const char *, int32_t &)
Definition: int32IO.C:85
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
error FatalError
static const char nl
Definition: Ostream.H:260
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