RASThermophysicalTransportModel.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-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 \*---------------------------------------------------------------------------*/
25 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasicThermophysicalTransportModel>
33 <
34  BasicThermophysicalTransportModel
36 (
37  const word& type,
38  const momentumTransportModel& momentumTransport,
39  const thermoModel& thermo
40 )
41 :
42  BasicThermophysicalTransportModel(momentumTransport, thermo)
43 {}
44 
45 
46 // * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * * //
47 
48 template<class BasicThermophysicalTransportModel>
50 <
52  <
53  BasicThermophysicalTransportModel
54  >
55 >
57 <
58  BasicThermophysicalTransportModel
59 >::New
60 (
61  const momentumTransportModel& momentumTransport,
62  const thermoModel& thermo
63 )
64 {
66  (
68  (
69  thermophysicalTransportModel::typeName,
70  momentumTransport.alphaRhoPhi().group()
71  ),
72  momentumTransport.time().constant(),
73  momentumTransport.mesh(),
76  false
77  );
78 
79  if (header.headerOk())
80  {
81  IOdictionary modelDict(header);
82 
83  const word modelType(modelDict.subDict("RAS").lookup( "model"));
84 
85  Info<< indent
86  << "Selecting RAS thermophysical transport model "
87  << modelType << endl;
88 
89  typename dictionaryConstructorTable::iterator cstrIter =
90  dictionaryConstructorTablePtr_->find(modelType);
91 
92  if (cstrIter == dictionaryConstructorTablePtr_->end())
93  {
95  << "Unknown RAS thermophysical transport model "
96  << modelType << nl << nl
97  << "Available models:" << endl
98  << dictionaryConstructorTablePtr_->sortedToc()
99  << exit(FatalError);
100  }
101 
102  Info<< incrIndent;
103 
104  autoPtr<RASThermophysicalTransportModel> modelPtr
105  (
106  cstrIter()(momentumTransport, thermo)
107  );
108 
109  Info<< decrIndent;
110 
111  return modelPtr;
112  }
113  else
114  {
115  typedef
116  turbulenceThermophysicalTransportModels::unityLewisEddyDiffusivity
117  <
118  RASThermophysicalTransportModel
119  <
120  BasicThermophysicalTransportModel
121  >
122  > RASunityLewisEddyDiffusivity;
123 
124  Info<< indent
125  << "Selecting default RAS thermophysical transport model "
126  << RASunityLewisEddyDiffusivity::typeName << endl;
127 
128  Info<< incrIndent;
129 
130  autoPtr<RASThermophysicalTransportModel> modelPtr
131  (
132  new RASunityLewisEddyDiffusivity
133  (
134  RASunityLewisEddyDiffusivity::typeName,
135  momentumTransport,
136  thermo,
137  true
138  )
139  );
140 
141  Info<< decrIndent;
142 
143  return modelPtr;
144  }
145 }
146 
147 
148 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
149 
150 template<class BasicThermophysicalTransportModel>
152 <
153  BasicThermophysicalTransportModel
154 >::coeffDict() const
155 {
156  return this->subOrEmptyDict("RAS").optionalSubDict(type() + "Coeffs");
157 }
158 
159 
160 template<class BasicThermophysicalTransportModel>
162 <
163  BasicThermophysicalTransportModel
165 {
167 }
168 
169 
170 template<class BasicThermophysicalTransportModel>
172 predict()
173 {
174  BasicThermophysicalTransportModel::predict();
175 }
176 
177 
178 template<class BasicThermophysicalTransportModel>
180 correct()
181 {
183 }
184 
185 
186 // ************************************************************************* //
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:131
static word groupName(Name name, const word &group)
Templated abstract base class for RAS thermophysical transport models.
virtual void correct()
Solve the thermophysical transport model equations.
virtual void predict()
Predict the RAS transport coefficients if possible.
BasicThermophysicalTransportModel::thermoModel thermoModel
BasicThermophysicalTransportModel::momentumTransportModel momentumTransportModel
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
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 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:530
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
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
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 & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:242
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:235
error FatalError
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:228
tmp< DimensionedField< TypeR, GeoMesh, Field > > New(const tmp< DimensionedField< TypeR, GeoMesh, Field >> &tdf1, const word &name, const dimensionSet &dimensions)
static const char nl
Definition: Ostream.H:267
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