solverTemplate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2015 OpenFOAM Foundation
9  Copyright (C) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "solverTemplate.H"
30 #include "Time.H"
31 #include "IOPtrList.H"
32 #include "polyMesh.H"
33 #include "regionProperties.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 const Foam::Enum
38 <
40 >
42 ({
43  { solverType::stCompressible, "compressible" },
44  { solverType::stIncompressible, "incompressible" },
45  { solverType::stBuoyant, "buoyant" },
46  { solverType::stUnknown, "unknown" },
47 });
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 Foam::word Foam::solverTemplate::readFromDict
53 (
54  IOobject& dictHeader,
55  const word& entryName
56 ) const
57 {
58  if (!dictHeader.typeHeaderOk<IOdictionary>(true))
59  {
61  << "Unable to open file "
62  << dictHeader.objectPath()
63  << exit(FatalError);
64  }
65 
66  IOdictionary dict(dictHeader);
67  return dict.get<word>(entryName);
68 }
69 
70 
71 Foam::dictionary Foam::solverTemplate::readFluidFieldTemplates
72 (
73  const word& regionName,
74  const fileName& baseDir,
75  const dictionary& solverDict,
76  const Time& runTime
77 ) const
78 {
79  Info<< " Reading fluid field templates";
80 
81  if (!regionName.empty())
82  {
83  Info<< " for region " << regionName;
84  }
85  Info<< endl;
86 
87  dictionary fieldTemplates = solverDict.subDict("fluidFields");
88 
89  const fileName turbModelDir(baseDir/"models"/"turbulence");
90 
91  const dictionary fieldModels(solverDict.subDict("fluidModels"));
92 
93  word turbModel("laminar"); // default to laminar
94  word turbType("none");
95 
96  if (fieldModels.readIfPresent("turbulenceModel", turbType))
97  {
98  if (turbType == "turbulenceModel")
99  {
100  IOdictionary turbPropDict
101  (
102  IOobject
103  (
104  // turbulenceModel::propertiesName
105  "turbulenceProperties",
106  runTime.constant(),
107  regionName,
108  runTime,
112  )
113  );
114 
115  const word modelType(turbPropDict.get<word>("simulationType"));
116 
117  if (modelType == "laminar")
118  {
119  // Leave as laminar
120  }
121  else if (modelType == "RAS")
122  {
123  // "RASModel" for v2006 and earlier
124  turbPropDict.subDict(modelType)
125  .readCompat("model", {{"RASModel", -2006}}, turbModel);
126  }
127  else if (modelType == "LES")
128  {
129  // "LESModel" for v2006 and earlier
130  turbPropDict.subDict(modelType)
131  .readCompat("model", {{"LESModel", -2006}}, turbModel);
132  }
133  else
134  {
136  << "Unhandled turbulence model option " << modelType
137  << ". Valid options are laminar, RAS, LES"
138  << exit(FatalError);
139  }
140  }
141  else
142  {
144  << "Unhandled turbulence model option " << turbType
145  << ". Valid options are turbulenceModel"
146  << exit(FatalError);
147  }
148  }
149 
150  Info<< " Selecting " << turbType << ": " << turbModel << endl;
151 
152  IOdictionary turbModelDict
153  (
154  IOobject
155  (
156  fileName(turbModelDir/turbModel),
157  runTime,
159  )
160  );
161 
162  // Merge common fluid fields
163  fieldTemplates.merge(turbModelDict.subDict("fluidFields"));
164 
165  // Merge specific compressible or incompressible fluid fields
166  switch (solverType_)
167  {
168  case stIncompressible:
169  {
170  fieldTemplates.merge(turbModelDict.subDict("incompressibleFields"));
171  break;
172  }
173  case stCompressible:
174  case stBuoyant:
175  {
176  fieldTemplates.merge(turbModelDict.subDict("compressibleFields"));
177  break;
178  }
179  default:
180  {
181  // do nothing
182  }
183  }
184 
185  return fieldTemplates;
186 }
187 
188 
189 Foam::dictionary Foam::solverTemplate::readSolidFieldTemplates
190 (
191  const word& regionName,
192  const dictionary& solverDict
193 ) const
194 {
195  Info<< " Reading solid field templates for region " << regionName
196  << endl;
197 
198  // nothing more to do for solids
199  return solverDict.subDict("solidFields");
200 }
201 
202 
203 void Foam::solverTemplate::setRegionProperties
204 (
205  const dictionary& fieldDict,
206  const word& regionType,
207  const word& regionName,
208  const label regionI
209 )
210 {
211  regionTypes_[regionI] = regionType,
212  regionNames_[regionI] = regionName,
213  fieldNames_[regionI] = fieldDict.toc();
214  fieldTypes_[regionI].setSize(fieldNames_[regionI].size());
215  fieldDimensions_[regionI].setSize(fieldNames_[regionI].size());
216 
217  forAll(fieldNames_[regionI], i)
218  {
219  const word& fieldName = fieldNames_[regionI][i];
220  const dictionary& dict = fieldDict.subDict(fieldName);
221 
222  dict.readEntry("type", fieldTypes_[regionI][i]);
223  fieldDimensions_[regionI].set
224  (
225  i,
226  new dimensionSet("dimensions", dict)
227  );
228  }
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
233 
235 (
236  const fileName& baseDir,
237  const Time& runTime,
238  const word& solverName
239 )
240 :
241  solverType_(stUnknown),
242  multiRegion_(false),
243  regionTypes_(),
244  fieldNames_(),
245  fieldTypes_(),
246  fieldDimensions_()
247 {
248  IOdictionary solverDict
249  (
250  IOobject
251  (
252  fileName(baseDir/"solvers"/solverName),
253  runTime,
255  )
256  );
257 
258  Info<< "Selecting " << solverName << ": ";
259 
260  solverType_ = solverTypeNames_.get("solverType", solverDict);
261  multiRegion_ = solverDict.get<bool>("multiRegion");
262 
263  Info<< solverTypeNames_[solverType_];
264  if (multiRegion_)
265  {
266  Info<< ", multi-region";
267  }
268  else
269  {
270  Info<< ", single-region";
271  }
272 
273  Info<< " case" << endl;
274 
275  if (multiRegion_)
276  {
277  // read regionProperties
278  regionProperties rp(runTime);
279 
280  const wordList fluidNames(rp["fluid"]);
281  const wordList solidNames(rp["solid"]);
282 
283  const label nRegion = fluidNames.size() + solidNames.size();
284 
285  regionTypes_.setSize(nRegion);
286  regionNames_.setSize(nRegion);
287  fieldNames_.setSize(nRegion);
288  fieldTypes_.setSize(nRegion);
289  fieldDimensions_.setSize(nRegion);
290 
291  // read templates for solver lists of available
292  // - fields and dimensions required for the solver
293  // - models
294  label regionI = 0;
295  forAll(fluidNames, i)
296  {
297  const dictionary fieldDict
298  (
299  readFluidFieldTemplates
300  (
301  fluidNames[i],
302  baseDir,
303  solverDict,
304  runTime
305  )
306  );
307 
308  setRegionProperties(fieldDict, "fluid", fluidNames[i], regionI++);
309  }
310 
311  forAll(solidNames, i)
312  {
313  const dictionary fieldDict
314  (
315  readSolidFieldTemplates
316  (
317  solidNames[i],
318  solverDict
319  )
320  );
321  setRegionProperties(fieldDict, "solid", solidNames[i], regionI++);
322  }
323  }
324  else
325  {
326  regionTypes_.setSize(1);
327  regionNames_.setSize(1);
328  fieldNames_.setSize(1);
329  fieldTypes_.setSize(1);
330  fieldDimensions_.setSize(1);
331 
332  // read templates for solver lists of available
333  // - fields and dimensions required for the solver
334  // - models
335  const dictionary fieldDict
336  (
337  readFluidFieldTemplates
338  (
339  word::null, // use region = null for single-region cases
340  baseDir,
341  solverDict,
342  runTime
343  )
344  );
345 
346  setRegionProperties(fieldDict, "fluid", word::null, 0);
347  }
348 }
349 
350 
351 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
352 
354 {
355  return solverTypeNames_[solverType_];
356 }
357 
358 
360 {
361  return multiRegion_;
362 }
363 
364 
365 Foam::label Foam::solverTemplate::nRegion() const
366 {
367  return regionTypes_.size();
368 }
369 
370 
371 const Foam::word& Foam::solverTemplate::regionType(const label regionI) const
372 {
373  return regionTypes_[regionI];
374 }
375 
376 
377 const Foam::word& Foam::solverTemplate::regionName(const label regionI) const
378 {
379  return regionNames_[regionI];
380 }
381 
382 
384 (
385  const label regionI
386 ) const
387 {
388  return fieldNames_[regionI];
389 }
390 
391 
393 (
394  const label regionI
395 ) const
396 {
397  return fieldTypes_[regionI];
398 }
399 
400 
402 (
403  const label regionI
404 ) const
405 {
406  return fieldDimensions_[regionI];
407 }
408 
409 
410 // ************************************************************************* //
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
solverTemplate(const fileName &baseDir, const Time &runTime, const word &regionName)
Constructor.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
word type() const
Solver type name.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const word & regionType(const label regionI) const
Return the region type.
const wordList solidNames(rp["solid"])
Ignore writing from objectRegistry::writeObject()
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:453
label nRegion() const
Return the number of regions.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
const wordList & fieldNames(const label regionI) const
Return the field names.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
Foam::word regionName(Foam::polyMesh::defaultRegion)
A class for handling words, derived from Foam::string.
Definition: word.H:63
static const Enum< solverType > solverTypeNames_
Solver type names.
regionProperties rp(runTime)
const PtrList< dimensionSet > & fieldDimensions(const label regionI) const
Return the field dimensions.
bool multiRegion() const
Return the multi-region flag.
static const word null
An empty word.
Definition: word.H:84
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:89
List< word > wordList
List of word.
Definition: fileName.H:58
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
solverType
Solver type.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const wordList & fieldTypes(const label regionI) const
Return the field types.
const wordList fluidNames(rp["fluid"])
Do not request registration (bool: false)
const word & regionName(const label regionI) const
Return the region name.