fvSchemes.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) 2011-2015 OpenFOAM Foundation
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 "fvSchemes.H"
30 #include "Time.H"
31 #include "steadyStateDdtScheme.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
36 
37 
38 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
39 
40 void Foam::fvSchemes::clear()
41 {
42  ddtSchemes_.clear();
43  defaultDdtScheme_.clear();
44  d2dt2Schemes_.clear();
45  defaultD2dt2Scheme_.clear();
46  interpolationSchemes_.clear();
47  defaultInterpolationScheme_.clear();
48  divSchemes_.clear();
49  defaultDivScheme_.clear();
50  gradSchemes_.clear();
51  defaultGradScheme_.clear();
52  snGradSchemes_.clear();
53  defaultSnGradScheme_.clear();
54  laplacianSchemes_.clear();
55  defaultLaplacianScheme_.clear();
56  // Do not clear fluxRequired settings
57 }
58 
59 
60 void Foam::fvSchemes::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  steady_ =
79  (
80  word(defaultDdtScheme_)
81  == fv::steadyStateDdtScheme<scalar>::typeName
82  );
83  }
84 
85 
86  if (dict.found("d2dt2Schemes"))
87  {
88  d2dt2Schemes_ = dict.subDict("d2dt2Schemes");
89  }
90  else
91  {
92  d2dt2Schemes_.set("default", "none");
93  }
94 
95  if
96  (
97  d2dt2Schemes_.found("default")
98  && word(d2dt2Schemes_.lookup("default")) != "none"
99  )
100  {
101  defaultD2dt2Scheme_ = d2dt2Schemes_.lookup("default");
102  }
103 
104 
105  if (dict.found("interpolationSchemes"))
106  {
107  interpolationSchemes_ = dict.subDict("interpolationSchemes");
108  }
109  else if (!interpolationSchemes_.found("default"))
110  {
111  interpolationSchemes_.add("default", "linear");
112  }
113 
114  if
115  (
116  interpolationSchemes_.found("default")
117  && word(interpolationSchemes_.lookup("default")) != "none"
118  )
119  {
120  defaultInterpolationScheme_ =
121  interpolationSchemes_.lookup("default");
122  }
123 
124 
125  divSchemes_ = dict.subDict("divSchemes");
126 
127  if
128  (
129  divSchemes_.found("default")
130  && word(divSchemes_.lookup("default")) != "none"
131  )
132  {
133  defaultDivScheme_ = divSchemes_.lookup("default");
134  }
135 
136 
137  gradSchemes_ = dict.subDict("gradSchemes");
138 
139  if
140  (
141  gradSchemes_.found("default")
142  && word(gradSchemes_.lookup("default")) != "none"
143  )
144  {
145  defaultGradScheme_ = gradSchemes_.lookup("default");
146  }
147 
148 
149  if (dict.found("snGradSchemes"))
150  {
151  snGradSchemes_ = dict.subDict("snGradSchemes");
152  }
153  else if (!snGradSchemes_.found("default"))
154  {
155  snGradSchemes_.add("default", "corrected");
156  }
157 
158  if
159  (
160  snGradSchemes_.found("default")
161  && word(snGradSchemes_.lookup("default")) != "none"
162  )
163  {
164  defaultSnGradScheme_ = snGradSchemes_.lookup("default");
165  }
166 
167 
168  laplacianSchemes_ = dict.subDict("laplacianSchemes");
169 
170  if
171  (
172  laplacianSchemes_.found("default")
173  && word(laplacianSchemes_.lookup("default")) != "none"
174  )
175  {
176  defaultLaplacianScheme_ = laplacianSchemes_.lookup("default");
177  }
178 
179 
180  if (dict.found("fluxRequired"))
181  {
182  fluxRequired_.merge(dict.subDict("fluxRequired"));
183 
184  if
185  (
186  fluxRequired_.found("default")
187  && fluxRequired_.get<word>("default") != "none"
188  )
189  {
190  defaultFluxRequired_ = fluxRequired_.get<bool>("default");
191  }
192  }
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
197 
198 Foam::fvSchemes::fvSchemes(const objectRegistry& obr)
199 :
201  (
202  IOobject
203  (
204  "fvSchemes",
205  obr.time().system(),
206  obr,
207  (
208  obr.readOpt() == IOobject::MUST_READ
209  || obr.readOpt() == IOobject::READ_IF_PRESENT
210  ? IOobject::MUST_READ_IF_MODIFIED
211  : obr.readOpt()
212  ),
213  IOobject::NO_WRITE
214  )
215  ),
216  ddtSchemes_
217  (
218  ITstream
219  (
220  objectPath() + ".ddtSchemes",
221  tokenList()
222  )()
223  ),
224  defaultDdtScheme_
225  (
226  ddtSchemes_.name() + ".default",
227  tokenList()
228  ),
229  d2dt2Schemes_
230  (
231  ITstream
232  (
233  objectPath() + ".d2dt2Schemes",
234  tokenList()
235  )()
236  ),
237  defaultD2dt2Scheme_
238  (
239  d2dt2Schemes_.name() + ".default",
240  tokenList()
241  ),
242  interpolationSchemes_
243  (
244  ITstream
245  (
246  objectPath() + ".interpolationSchemes",
247  tokenList()
248  )()
249  ),
250  defaultInterpolationScheme_
251  (
252  interpolationSchemes_.name() + ".default",
253  tokenList()
254  ),
255  divSchemes_
256  (
257  ITstream
258  (
259  objectPath() + ".divSchemes",
260  tokenList()
261  )()
262  ),
263  defaultDivScheme_
264  (
265  divSchemes_.name() + ".default",
266  tokenList()
267  ),
268  gradSchemes_
269  (
270  ITstream
271  (
272  objectPath() + ".gradSchemes",
273  tokenList()
274  )()
275  ),
276  defaultGradScheme_
277  (
278  gradSchemes_.name() + ".default",
279  tokenList()
280  ),
281  snGradSchemes_
282  (
283  ITstream
284  (
285  objectPath() + ".snGradSchemes",
286  tokenList()
287  )()
288  ),
289  defaultSnGradScheme_
290  (
291  snGradSchemes_.name() + ".default",
292  tokenList()
293  ),
294  laplacianSchemes_
295  (
296  ITstream
297  (
298  objectPath() + ".laplacianSchemes",
299  tokenList()
300  )()
301  ),
302  defaultLaplacianScheme_
303  (
304  laplacianSchemes_.name() + ".default",
305  tokenList()
306  ),
307  fluxRequired_
308  (
309  ITstream
310  (
311  objectPath() + ".fluxRequired",
312  tokenList()
313  )()
314  ),
315  defaultFluxRequired_(false),
316  steady_(false)
317 {
318  if
319  (
323  )
324  {
325  read(schemesDict());
326  }
327 }
328 
329 
330 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
331 
333 {
334  if (regIOobject::read())
335  {
336  // Clear current settings except fluxRequired
337  clear();
338 
339  read(schemesDict());
340 
341  return true;
342  }
343 
344  return false;
345 }
346 
347 
349 {
350  if (found("select"))
351  {
352  return subDict(word(lookup("select")));
353  }
354 
355  return *this;
356 }
357 
358 
360 {
361  if (debug)
362  {
363  Info<< "Lookup ddtScheme for " << name << endl;
364  }
365 
366  if (ddtSchemes_.found(name) || defaultDdtScheme_.empty())
367  {
368  return ddtSchemes_.lookup(name);
369  }
370  else
371  {
372  const_cast<ITstream&>(defaultDdtScheme_).rewind();
373  return const_cast<ITstream&>(defaultDdtScheme_);
374  }
375 }
376 
377 
379 {
380  if (debug)
381  {
382  Info<< "Lookup d2dt2Scheme for " << name << endl;
383  }
384 
385  if (d2dt2Schemes_.found(name) || defaultD2dt2Scheme_.empty())
386  {
387  return d2dt2Schemes_.lookup(name);
388  }
389  else
390  {
391  const_cast<ITstream&>(defaultD2dt2Scheme_).rewind();
392  return const_cast<ITstream&>(defaultD2dt2Scheme_);
393  }
394 }
395 
396 
398 {
399  if (debug)
400  {
401  Info<< "Lookup interpolationScheme for " << name << endl;
402  }
403 
404  if
405  (
406  interpolationSchemes_.found(name)
407  || defaultInterpolationScheme_.empty()
408  )
409  {
410  return interpolationSchemes_.lookup(name);
411  }
412  else
413  {
414  const_cast<ITstream&>(defaultInterpolationScheme_).rewind();
415  return const_cast<ITstream&>(defaultInterpolationScheme_);
416  }
417 }
418 
419 
421 {
422  if (debug)
423  {
424  Info<< "Lookup divScheme for " << name << endl;
425  }
426 
427  if (divSchemes_.found(name) || defaultDivScheme_.empty())
428  {
429  return divSchemes_.lookup(name);
430  }
431  else
432  {
433  const_cast<ITstream&>(defaultDivScheme_).rewind();
434  return const_cast<ITstream&>(defaultDivScheme_);
435  }
436 }
437 
438 
440 {
441  if (debug)
442  {
443  Info<< "Lookup gradScheme for " << name << endl;
444  }
445 
446  if (gradSchemes_.found(name) || defaultGradScheme_.empty())
447  {
448  return gradSchemes_.lookup(name);
449  }
450  else
451  {
452  const_cast<ITstream&>(defaultGradScheme_).rewind();
453  return const_cast<ITstream&>(defaultGradScheme_);
454  }
455 }
456 
457 
459 {
460  if (debug)
461  {
462  Info<< "Lookup snGradScheme for " << name << endl;
463  }
464 
465  if (snGradSchemes_.found(name) || defaultSnGradScheme_.empty())
466  {
467  return snGradSchemes_.lookup(name);
468  }
469  else
470  {
471  const_cast<ITstream&>(defaultSnGradScheme_).rewind();
472  return const_cast<ITstream&>(defaultSnGradScheme_);
473  }
474 }
475 
476 
478 {
479  if (debug)
480  {
481  Info<< "Lookup laplacianScheme for " << name << endl;
482  }
483 
484  if (laplacianSchemes_.found(name) || defaultLaplacianScheme_.empty())
485  {
486  return laplacianSchemes_.lookup(name);
487  }
488  else
489  {
490  const_cast<ITstream&>(defaultLaplacianScheme_).rewind();
491  return const_cast<ITstream&>(defaultLaplacianScheme_);
492  }
493 }
494 
495 
497 {
498  if (debug)
499  {
500  Info<< "Setting fluxRequired for " << name << endl;
501  }
502 
503  fluxRequired_.add(name, true, true);
504 }
505 
506 
508 {
509  if (debug)
510  {
511  Info<< "Lookup fluxRequired for " << name << endl;
512  }
513 
514  if (fluxRequired_.found(name))
515  {
516  return true;
517  }
518 
519  return defaultFluxRequired_;
520 }
521 
522 
523 // ************************************************************************* //
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::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
steadyStateDdtScheme.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fvSchemes::ddtScheme
ITstream & ddtScheme(const word &name) const
Definition: fvSchemes.C:359
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:359
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::regIOobject::read
virtual bool read()
Read object.
Definition: regIOobjectRead.C:202
Foam::fvSchemes::divScheme
ITstream & divScheme(const word &name) const
Definition: fvSchemes.C:420
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::set
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:842
Foam::fvSchemes::read
bool read()
Read the fvSchemes.
Definition: fvSchemes.C:332
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::fvSchemes::debug
static int debug
Debug switch.
Definition: fvSchemes.H:105
Foam::fvSchemes::setFluxRequired
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:496
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::cellModeller::lookup
const cellModel * lookup(const word &modelName)
Deprecated(2017-11) equivalent to cellModel::ptr static method.
Definition: cellModeller.H:46
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:523
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:122
Foam::fvSchemes::gradScheme
ITstream & gradScheme(const word &name) const
Definition: fvSchemes.C:439
Foam::fvSchemes::schemesDict
const dictionary & schemesDict() const
Definition: fvSchemes.C:348
fvSchemes.H
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:419
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::fvSchemes::snGradScheme
ITstream & snGradScheme(const word &name) const
Definition: fvSchemes.C:458
Foam::fvSchemes::fluxRequired
bool fluxRequired(const word &name) const
Definition: fvSchemes.C:507
found
bool found
Definition: TABSMDCalcMethod2.H:32
Time.H
clear
patchWriters clear()
Foam::fvSchemes::interpolationScheme
ITstream & interpolationScheme(const word &name) const
Definition: fvSchemes.C:397
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:102
Foam::IOobject::readOpt
readOption readOpt() const
The read option.
Definition: IOobjectI.H:141
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::fvSchemes::d2dt2Scheme
ITstream & d2dt2Scheme(const word &name) const
Definition: fvSchemes.C:378
Foam::dictionary::add
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition: dictionary.C:703
Foam::fvSchemes::laplacianScheme
ITstream & laplacianScheme(const word &name) const
Definition: fvSchemes.C:477
Foam::dictionary::clear
void clear()
Clear the dictionary.
Definition: dictionary.C:919
Foam::regIOobject::headerOk
bool headerOk()
Read and check header info.
Definition: regIOobject.C:449
Foam::IOobject::MUST_READ
Definition: IOobject.H:120