faSchemes.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) 2016-2017 Wikki Ltd
9  Copyright (C) 2019 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 "faSchemes.H"
30 #include "Time.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
35 
36 
37 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
38 
39 void Foam::faSchemes::clear()
40 {
41  ddtSchemes_.clear();
42  defaultDdtScheme_.clear();
43  d2dt2Schemes_.clear();
44  defaultD2dt2Scheme_.clear();
45  interpolationSchemes_.clear();
46  defaultInterpolationScheme_.clear();
47  divSchemes_.clear(); // optional
48  defaultDivScheme_.clear();
49  gradSchemes_.clear(); // optional
50  defaultGradScheme_.clear();
51  lnGradSchemes_.clear();
52  defaultLnGradScheme_.clear();
53  laplacianSchemes_.clear(); // optional
54  defaultLaplacianScheme_.clear();
55  fluxRequired_.clear();
56  defaultFluxRequired_ = false;
57 }
58 
59 
60 void Foam::faSchemes::read(const dictionary& dict)
61 {
62  if (dict.found("ddtSchemes"))
63  {
64  ddtSchemes_ = dict.subDict("ddtSchemes");
65  }
66  else
67  {
68  ddtSchemes_.set("default", "none");
69  }
70 
71  if
72  (
73  ddtSchemes_.found("default")
74  && word(ddtSchemes_.lookup("default")) != "none"
75  )
76  {
77  defaultDdtScheme_ = ddtSchemes_.lookup("default");
78 
79  }
80 
81 
82  if (dict.found("d2dt2Schemes"))
83  {
84  d2dt2Schemes_ = dict.subDict("d2dt2Schemes");
85  }
86  else
87  {
88  d2dt2Schemes_.set("default", "none");
89  }
90 
91  if
92  (
93  d2dt2Schemes_.found("default")
94  && word(d2dt2Schemes_.lookup("default")) != "none"
95  )
96  {
97  defaultD2dt2Scheme_ = d2dt2Schemes_.lookup("default");
98  }
99 
100 
101  if (dict.found("interpolationSchemes"))
102  {
103  interpolationSchemes_ = dict.subDict("interpolationSchemes");
104  }
105  else if (!interpolationSchemes_.found("default"))
106  {
107  interpolationSchemes_.add("default", "linear");
108  }
109 
110  if
111  (
112  interpolationSchemes_.found("default")
113  && word(interpolationSchemes_.lookup("default")) != "none"
114  )
115  {
116  defaultInterpolationScheme_ =
117  interpolationSchemes_.lookup("default");
118  }
119 
120 
121  divSchemes_ = dict.subDict("divSchemes");
122 
123  if
124  (
125  divSchemes_.found("default")
126  && word(divSchemes_.lookup("default")) != "none"
127  )
128  {
129  defaultDivScheme_ = divSchemes_.lookup("default");
130  }
131 
132 
133  gradSchemes_ = dict.subDict("gradSchemes");
134 
135  if
136  (
137  gradSchemes_.found("default")
138  && word(gradSchemes_.lookup("default")) != "none"
139  )
140  {
141  defaultGradScheme_ = gradSchemes_.lookup("default");
142  }
143 
144 
145  if (dict.found("lnGradSchemes"))
146  {
147  lnGradSchemes_ = dict.subDict("lnGradSchemes");
148  }
149  else if (!lnGradSchemes_.found("default"))
150  {
151  lnGradSchemes_.add("default", "corrected");
152  }
153 
154  if
155  (
156  lnGradSchemes_.found("default")
157  && word(lnGradSchemes_.lookup("default")) != "none"
158  )
159  {
160  defaultLnGradScheme_ = lnGradSchemes_.lookup("default");
161  }
162 
163 
164  laplacianSchemes_ = dict.subDict("laplacianSchemes");
165 
166  if
167  (
168  laplacianSchemes_.found("default")
169  && word(laplacianSchemes_.lookup("default")) != "none"
170  )
171  {
172  defaultLaplacianScheme_ = laplacianSchemes_.lookup("default");
173  }
174 
175 
176  if (dict.found("fluxRequired"))
177  {
178  fluxRequired_.merge(dict.subDict("fluxRequired"));
179 
180  if
181  (
182  fluxRequired_.found("default")
183  && fluxRequired_.get<word>("default") != "none"
184  )
185  {
186  defaultFluxRequired_ = fluxRequired_.get<bool>("default");
187  }
188  }
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
193 
194 Foam::faSchemes::faSchemes(const objectRegistry& obr)
195 :
197  (
198  IOobject
199  (
200  "faSchemes",
201  obr.time().system(),
202  obr,
203  (
204  obr.readOpt() == IOobject::MUST_READ
205  || obr.readOpt() == IOobject::READ_IF_PRESENT
206  ? IOobject::MUST_READ_IF_MODIFIED
207  : obr.readOpt()
208  ),
209  IOobject::NO_WRITE
210  )
211  ),
212  ddtSchemes_
213  (
214  ITstream
215  (
216  objectPath() + ".ddtSchemes",
217  tokenList()
218  )()
219  ),
220  defaultDdtScheme_
221  (
222  ddtSchemes_.name() + ".default",
223  tokenList()
224  ),
225  d2dt2Schemes_
226  (
227  ITstream
228  (
229  objectPath() + ".d2dt2Schemes",
230  tokenList()
231  )()
232  ),
233  defaultD2dt2Scheme_
234  (
235  d2dt2Schemes_.name() + ".default",
236  tokenList()
237  ),
238  interpolationSchemes_
239  (
240  ITstream
241  (
242  objectPath() + ".interpolationSchemes",
243  tokenList()
244  )()
245  ),
246  defaultInterpolationScheme_
247  (
248  interpolationSchemes_.name() + ".default",
249  tokenList()
250  ),
251  divSchemes_
252  (
253  ITstream
254  (
255  objectPath() + ".divSchemes",
256  tokenList()
257  )()
258  ),
259  defaultDivScheme_
260  (
261  divSchemes_.name() + ".default",
262  tokenList()
263  ),
264  gradSchemes_
265  (
266  ITstream
267  (
268  objectPath() + ".gradSchemes",
269  tokenList()
270  )()
271  ),
272  defaultGradScheme_
273  (
274  gradSchemes_.name() + ".default",
275  tokenList()
276  ),
277  lnGradSchemes_
278  (
279  ITstream
280  (
281  objectPath() + ".snGradSchemes",
282  tokenList()
283  )()
284  ),
285  defaultLnGradScheme_
286  (
287  lnGradSchemes_.name() + ".default",
288  tokenList()
289  ),
290  laplacianSchemes_
291  (
292  ITstream
293  (
294  objectPath() + ".laplacianSchemes",
295  tokenList()
296  )()
297  ),
298  defaultLaplacianScheme_
299  (
300  laplacianSchemes_.name() + ".default",
301  tokenList()
302  ),
303  fluxRequired_
304  (
305  ITstream
306  (
307  objectPath() + ".fluxRequired",
308  tokenList()
309  )()
310  ),
311  defaultFluxRequired_(false)
312 {
313  if
314  (
318  )
319  {
320  read(schemesDict());
321  }
322 }
323 
324 
325 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
326 
328 {
329  if (regIOobject::read())
330  {
331  // Clear current settings except fluxRequired
332  clear();
333 
334  read(schemesDict());
335 
336  return true;
337  }
338 
339  return false;
340 }
341 
342 
344 {
345  if (found("select"))
346  {
347  return subDict(word(lookup("select")));
348  }
349 
350  return *this;
351 }
352 
353 
355 {
356  if (debug)
357  {
358  Info<< "Lookup ddtScheme for " << name << endl;
359  }
360 
361  if (ddtSchemes_.found(name) || defaultDdtScheme_.empty())
362  {
363  return ddtSchemes_.lookup(name);
364  }
365  else
366  {
367  const_cast<ITstream&>(defaultDdtScheme_).rewind();
368  return const_cast<ITstream&>(defaultDdtScheme_);
369  }
370 }
371 
372 
374 {
375  if (debug)
376  {
377  Info<< "Lookup d2dt2Scheme for " << name << endl;
378  }
379 
380  if (d2dt2Schemes_.found(name) || defaultD2dt2Scheme_.empty())
381  {
382  return d2dt2Schemes_.lookup(name);
383  }
384  else
385  {
386  const_cast<ITstream&>(defaultD2dt2Scheme_).rewind();
387  return const_cast<ITstream&>(defaultD2dt2Scheme_);
388  }
389 }
390 
391 
393 {
394  if (debug)
395  {
396  Info<< "Lookup interpolationScheme for " << name << endl;
397  }
398 
399  if
400  (
401  interpolationSchemes_.found(name)
402  || defaultInterpolationScheme_.empty()
403  )
404  {
405  return interpolationSchemes_.lookup(name);
406  }
407  else
408  {
409  const_cast<ITstream&>(defaultInterpolationScheme_).rewind();
410  return const_cast<ITstream&>(defaultInterpolationScheme_);
411  }
412 }
413 
414 
416 {
417  if (debug)
418  {
419  Info<< "Lookup divScheme for " << name << endl;
420  }
421 
422  if (divSchemes_.found(name) || defaultDivScheme_.empty())
423  {
424  return divSchemes_.lookup(name);
425  }
426  else
427  {
428  const_cast<ITstream&>(defaultDivScheme_).rewind();
429  return const_cast<ITstream&>(defaultDivScheme_);
430  }
431 }
432 
433 
435 {
436  if (debug)
437  {
438  Info<< "Lookup gradScheme for " << name << endl;
439  }
440 
441  if (gradSchemes_.found(name) || defaultGradScheme_.empty())
442  {
443  return gradSchemes_.lookup(name);
444  }
445  else
446  {
447  const_cast<ITstream&>(defaultGradScheme_).rewind();
448  return const_cast<ITstream&>(defaultGradScheme_);
449  }
450 }
451 
452 
454 {
455  if (debug)
456  {
457  Info<< "Lookup snGradScheme for " << name << endl;
458  }
459 
460  if (lnGradSchemes_.found(name) || defaultLnGradScheme_.empty())
461  {
462  return lnGradSchemes_.lookup(name);
463  }
464  else
465  {
466  const_cast<ITstream&>(defaultLnGradScheme_).rewind();
467  return const_cast<ITstream&>(defaultLnGradScheme_);
468  }
469 }
470 
471 
473 {
474  if (debug)
475  {
476  Info<< "Lookup laplacianScheme for " << name << endl;
477  }
478 
479  if (laplacianSchemes_.found(name) || defaultLaplacianScheme_.empty())
480  {
481  return laplacianSchemes_.lookup(name);
482  }
483  else
484  {
485  const_cast<ITstream&>(defaultLaplacianScheme_).rewind();
486  return const_cast<ITstream&>(defaultLaplacianScheme_);
487  }
488 }
489 
490 
492 {
493  if (debug)
494  {
495  Info<< "Setting fluxRequired for " << name << endl;
496  }
497 
498  fluxRequired_.add(name, true, true);
499 }
500 
501 
503 {
504  if (debug)
505  {
506  Info<< "Lookup fluxRequired for " << name << endl;
507  }
508 
509  if (fluxRequired_.found(name))
510  {
511  return true;
512  }
513  else
514  {
515  return defaultFluxRequired_;
516  }
517 }
518 
519 
521 {
522  // Write dictionaries
523  os << nl << "ddtSchemes";
524  ddtSchemes_.write(os, true);
525 
526  os << nl << "d2dt2Schemes";
527  d2dt2Schemes_.write(os, true);
528 
529  os << nl << "interpolationSchemes";
530  interpolationSchemes_.write(os, true);
531 
532  os << nl << "divSchemes";
533  divSchemes_.write(os, true);
534 
535  os << nl << "gradSchemes";
536  gradSchemes_.write(os, true);
537 
538  os << nl << "lnGradSchemes";
539  lnGradSchemes_.write(os, true);
540 
541  os << nl << "laplacianSchemes";
542  laplacianSchemes_.write(os, true);
543 
544  os << nl << "fluxRequired";
545  fluxRequired_.write(os, true);
546 
547  return true;
548 }
549 
550 
551 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::faSchemes::schemesDict
const dictionary & schemesDict() const
Definition: faSchemes.C:343
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:364
Foam::debug::debugSwitch
int debugSwitch(const char *name, const int deflt=0)
Lookup debug switch or add default value.
Definition: debug.C:225
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::system
int system(const std::string &command, const bool bg=false)
Execute the specified command via the shell.
Definition: MSwindows.C:1140
Foam::faSchemes::writeData
virtual bool writeData(Ostream &) const
WriteData function required for regIOobject write operation.
Definition: faSchemes.C:520
Foam::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:202
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:847
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::faSchemes::laplacianScheme
ITstream & laplacianScheme(const word &name) const
Definition: faSchemes.C:472
Foam::ITstream
An input stream of tokens.
Definition: ITstream.H:55
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::faSchemes::divScheme
ITstream & divScheme(const word &name) const
Definition: faSchemes.C:415
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::faSchemes::interpolationScheme
ITstream & interpolationScheme(const word &name) const
Definition: faSchemes.C:392
Foam::faSchemes::lnGradScheme
ITstream & lnGradScheme(const word &name) const
Definition: faSchemes.C:453
Foam::radiation::lookup
Lookup type of boundary radiation properties.
Definition: lookup.H:63
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:424
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::faSchemes::debug
static int debug
Debug switch.
Definition: faSchemes.H:101
Foam::faSchemes::ddtScheme
ITstream & ddtScheme(const word &name) const
Definition: faSchemes.C:354
found
bool found
Definition: TABSMDCalcMethod2.H:32
Time.H
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:385
faSchemes.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::IOobject::readOpt
readOption readOpt() const
The read option.
Definition: IOobjectI.H:165
Foam::faSchemes::d2dt2Scheme
ITstream & d2dt2Scheme(const word &name) const
Definition: faSchemes.C:373
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:115
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:121
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:708
Foam::faSchemes::setFluxRequired
void setFluxRequired(const word &name) const
Definition: faSchemes.C:491
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::dictionary::clear
void clear()
Clear the dictionary.
Definition: dictionary.C:924
Foam::faSchemes::gradScheme
ITstream & gradScheme(const word &name) const
Definition: faSchemes.C:434
Foam::regIOobject::headerOk
bool headerOk()
Read and check header info.
Definition: regIOobject.C:453
Foam::faSchemes::fluxRequired
dictionary & fluxRequired()
Return access to flux required.
Definition: faSchemes.H:180
Foam::faSchemes::read
bool read()
Read the faSchemes.
Definition: faSchemes.C:327
Foam::IOobject::MUST_READ
Definition: IOobject.H:120