stabilityBlendingFactor.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) 2018-2024 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
31 #include "fvcGrad.H"
32 #include "surfaceInterpolate.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38 namespace functionObjects
39 {
40  defineTypeNameAndDebug(stabilityBlendingFactor, 0);
42  (
43  functionObject,
44  stabilityBlendingFactor,
45  dictionary
46  );
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52 
54 Foam::functionObjects::stabilityBlendingFactor::indicator()
55 {
56  const word fldName = "blendedIndicator" + fieldName_;
57 
58  auto* fldPtr = mesh_.getObjectPtr<volScalarField>(fldName);
59 
60  if (!fldPtr)
61  {
62  fldPtr = new volScalarField
63  (
64  IOobject
65  (
66  "blendedIndicator" + fieldName_,
67  time_.timeName(),
68  mesh_.thisDb(),
72  ),
73  mesh_,
76  );
77 
78  regIOobject::store(fldPtr);
79  }
80 
81  return *fldPtr;
82 }
83 
84 
85 void Foam::functionObjects::stabilityBlendingFactor::calcStats
86 (
87  label& nCellsScheme1,
88  label& nCellsScheme2,
89  label& nCellsBlended
90 ) const
91 {
92  const auto* indicatorPtr =
93  mesh_.cfindObject<volScalarField>("blendedIndicator" + fieldName_);
94 
95  if (!indicatorPtr)
96  {
98  << "Indicator field not set"
99  << abort(FatalError);
100  }
101 
102  const auto& indicator = *indicatorPtr;
103 
104  for (const scalar i : indicator)
105  {
106  if (i < tolerance_)
107  {
108  ++nCellsScheme2;
109  }
110  else if (i > (1 - tolerance_))
111  {
112  ++nCellsScheme1;
113  }
114  else
115  {
116  ++nCellsBlended;
117  }
118  }
119 
120  reduce(nCellsScheme1, sumOp<label>());
121  reduce(nCellsScheme2, sumOp<label>());
122  reduce(nCellsBlended, sumOp<label>());
123 }
124 
125 
127 (
128  Ostream& os
129 ) const
130 {
131  writeHeader(os, "Stabilization blending factor");
132  writeCommented(os, "Time");
133  writeTabbed(os, "Scheme1");
134  writeTabbed(os, "Scheme2");
135  writeTabbed(os, "Blended");
136  os << endl;
137 }
138 
139 
140 bool Foam::functionObjects::stabilityBlendingFactor::calc()
141 {
142  init(false);
143  return true;
144 }
145 
146 
147 bool Foam::functionObjects::stabilityBlendingFactor::init(bool first)
148 {
149  const auto* residualPtr = mesh_.findObject<IOField<scalar>>(residualName_);
150 
151  auto& indicator = this->indicator();
152 
153  if (residuals_)
154  {
155  if (!residualPtr)
156  {
158  << "Could not find residual field : " << residualName_
159  << ". The residual mode will not be " << nl
160  << " considered for the blended field in the stability "
161  << "blending factor. " << nl
162  << " Please add the corresponding 'solverInfo'"
163  << " function object with 'writeResidualFields true'." << nl
164  << " If the solverInfo function object is already enabled "
165  << "you might need to wait " << nl
166  << " for the first iteration."
167  << nl << endl;
168  }
169  else
170  {
171  scalar meanRes = gAverage(mag(*residualPtr)) + VSMALL;
172 
173  Log << nl << name() << " : " << nl
174  << " Average(mag(residuals)) : " << meanRes << endl;
175 
176  oldError_ = error_;
177  oldErrorIntegral_ = errorIntegral_;
178  error_ = mag(meanRes - mag(*residualPtr));
179  errorIntegral_ = oldErrorIntegral_ + 0.5*(error_ + oldError_);
180  const scalarField errorDifferential(error_ - oldError_);
181 
182  const scalarField factorList
183  (
184  + P_*error_
185  + I_*errorIntegral_
186  + D_*errorDifferential
187  );
188 
189  const scalarField indicatorResidual
190  (
191  max
192  (
193  min
194  (
195  mag(factorList - meanRes)/(maxResidual_*meanRes),
196  scalar(1)
197  ),
198  scalar(0)
199  )
200  );
201 
202  forAll(indicator, i)
203  {
204  indicator[i] = indicatorResidual[i];
205  }
206  }
207  }
208 
209  const volScalarField* nonOrthPtr =
210  mesh_.findObject<volScalarField>(nonOrthogonalityName_);
211 
212  if (nonOrthogonality_)
213  {
214  if (maxNonOrthogonality_ >= minNonOrthogonality_)
215  {
217  << "minNonOrthogonality should be larger than "
218  << "maxNonOrthogonality."
219  << exit(FatalError);
220  }
221 
222  indicator =
223  max
224  (
225  indicator,
226  min
227  (
228  max
229  (
230  scalar(0),
231  (*nonOrthPtr - maxNonOrthogonality_)
232  / (minNonOrthogonality_ - maxNonOrthogonality_)
233  ),
234  scalar(1)
235  )
236  );
237 
238  if (first)
239  {
240  Log << " Max non-orthogonality : " << max(*nonOrthPtr).value()
241  << endl;
242  }
243  }
244 
245  const auto* skewnessPtr = mesh_.cfindObject<volScalarField>(skewnessName_);
246 
247  if (skewness_)
248  {
249  if (maxSkewness_ >= minSkewness_)
250  {
252  << "minSkewness should be larger than maxSkewness."
253  << exit(FatalError);
254  }
255 
256  indicator =
257  max
258  (
259  indicator,
260  min
261  (
262  max
263  (
264  scalar(0),
265  (*skewnessPtr - maxSkewness_)
266  / (minSkewness_ - maxSkewness_)
267  ),
268  scalar(1)
269  )
270  );
271 
272  if (first)
273  {
274  Log << " Max skewness : " << max(*skewnessPtr).value()
275  << endl;
276  }
277  }
278 
279  const auto* faceWeightsPtr =
280  mesh_.cfindObject<volScalarField>(faceWeightName_);
281 
282  if (faceWeight_)
283  {
284  if (maxFaceWeight_ >= minFaceWeight_)
285  {
287  << "minFaceWeight should be larger than maxFaceWeight."
288  << exit(FatalError);
289  }
290 
291  indicator =
292  max
293  (
294  indicator,
295  min
296  (
297  max
298  (
299  scalar(0),
300  (minFaceWeight_ - *faceWeightsPtr)
301  / (minFaceWeight_ - maxFaceWeight_)
302  ),
303  scalar(1)
304  )
305  );
306 
307  if (first)
308  {
309  Log << " Min face weight: " << min(*faceWeightsPtr).value()
310  << endl;
311  }
312  }
313 
314 
315  if (gradCc_)
316  {
317  if (maxGradCc_ >= minGradCc_)
318  {
320  << "minGradCc should be larger than maxGradCc."
321  << exit(FatalError);
322  }
323 
324  auto tmagGradCC = volScalarField::New
325  (
326  "magGradCC",
328  mesh_,
331  );
332  auto& magGradCC = tmagGradCC.ref();
333 
334  for (direction i=0; i<vector::nComponents; i++)
335  {
336  // Create field with zero grad
337  auto tcci = volScalarField::New
338  (
339  "cc" + word(vector::componentNames[i]),
341  mesh_,
344  );
345  auto& cci = tcci.ref();
346 
347  cci = mesh_.C().component(i);
348  cci.correctBoundaryConditions();
349  magGradCC += mag(fvc::grad(cci)).ref();
350  }
351 
352  if (first)
353  {
354  Log << " Max magGradCc : " << max(magGradCC.ref()).value()
355  << endl;
356  }
357 
358  indicator =
359  max
360  (
361  indicator,
362  min
363  (
364  max
365  (
366  scalar(0),
367  (magGradCC - maxGradCc_)
368  / (minGradCc_ - maxGradCc_)
369  ),
370  scalar(1)
371  )
372  );
373  }
374 
375 
376  const auto* UNamePtr = mesh_.findObject<volVectorField>(UName_);
377 
378  if (Co_)
379  {
380  if (Co1_ >= Co2_)
381  {
383  << "Co2 should be larger than Co1."
384  << exit(FatalError);
385  }
386 
387  auto CoPtr = volScalarField::New
388  (
389  "Co",
391  mesh_,
394  );
395  auto& Co = CoPtr.ref();
396 
397  Co.primitiveFieldRef() =
398  mesh_.time().deltaT()*mag(*UNamePtr)/cbrt(mesh_.V());
399 
400  indicator =
401  max
402  (
403  indicator,
404  clamp((Co - Co1_)/(Co2_ - Co1_), zero_one{})
405  );
406 
407  if (first)
408  {
409  Log << " Max Co : " << max(Co).value()
410  << endl;
411  }
412  }
413 
414  if (first)
415  {
416  Log << nl;
417  }
418 
419  if (log)
420  {
421  label nCellsScheme1 = 0;
422  label nCellsScheme2 = 0;
423  label nCellsBlended = 0;
424 
425  calcStats(nCellsScheme1, nCellsScheme2, nCellsBlended);
426 
427  Log << nl << name() << " execute :" << nl
428  << " scheme 1 cells : " << nCellsScheme1 << nl
429  << " scheme 2 cells : " << nCellsScheme2 << nl
430  << " blended cells : " << nCellsBlended << nl
431  << endl;
432  }
433 
434  indicator.correctBoundaryConditions();
435  indicator.clamp_range(zero_one{});
436 
437  // Update the blended surface field
438  auto& surBlended = mesh_.lookupObjectRef<surfaceScalarField>(resultName_);
439 
440  surBlended = fvc::interpolate(indicator);
441 
442  return true;
443 }
444 
445 
446 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
447 
449 (
450  const word& name,
451  const Time& runTime,
452  const dictionary& dict
453 )
454 :
456  writeFile(obr_, name, typeName, dict),
457  nonOrthogonality_(dict.getOrDefault<Switch>("switchNonOrtho", false)),
458  gradCc_(dict.getOrDefault<Switch>("switchGradCc", false)),
459  residuals_(dict.getOrDefault<Switch>("switchResiduals", false)),
460  faceWeight_(dict.getOrDefault<Switch>("switchFaceWeight", false)),
461  skewness_(dict.getOrDefault<Switch>("switchSkewness", false)),
462  Co_(dict.getOrDefault<Switch>("switchCo", false)),
463 
464  maxNonOrthogonality_
465  (
466  dict.getOrDefault<scalar>("maxNonOrthogonality", 20.0)
467  ),
468  minNonOrthogonality_
469  (
470  dict.getOrDefault<scalar>("minNonOrthogonality", 60.0)
471  ),
472  maxGradCc_(dict.getOrDefault<scalar>("maxGradCc", 3.0)),
473  minGradCc_(dict.getOrDefault<scalar>("minGradCc", 4.0)),
474  maxResidual_(dict.getOrDefault<scalar>("maxResidual", 10.0)),
475  minFaceWeight_(dict.getOrDefault<scalar>("minFaceWeight", 0.3)),
476  maxFaceWeight_(dict.getOrDefault<scalar>("maxFaceWeight", 0.2)),
477  maxSkewness_(dict.getOrDefault<scalar>("maxSkewness", 2.0)),
478  minSkewness_(dict.getOrDefault<scalar>("minSkewness", 3.0)),
479  Co1_(dict.getOrDefault<scalar>("Co1", 1.0)),
480  Co2_(dict.getOrDefault<scalar>("Co2", 10.0)),
481 
482  nonOrthogonalityName_
483  (
484  dict.getOrDefault<word>("nonOrthogonality", "nonOrthoAngle")
485  ),
486  faceWeightName_
487  (
488  dict.getOrDefault<word>("faceWeight", "faceWeight")
489  ),
490  skewnessName_
491  (
492  dict.getOrDefault<word>("skewness", "skewness")
493  ),
494  residualName_
495  (
496  dict.getOrDefault<word>
497  (
498  "residual",
499  IOobject::scopedName("initialResidual", "p")
500  )
501  ),
502  UName_
503  (
504  dict.getOrDefault<word>("U", "U")
505  ),
506 
507  tolerance_(0.001),
508  error_(mesh_.nCells(), Zero),
509  errorIntegral_(mesh_.nCells(), Zero),
510  oldError_(mesh_.nCells(), Zero),
511  oldErrorIntegral_(mesh_.nCells(), Zero),
512  P_(dict.getOrDefault<scalar>("P", 3)),
513  I_(dict.getOrDefault<scalar>("I", 0.0)),
514  D_(dict.getOrDefault<scalar>("D", 0.25))
515 {
516  read(dict);
517  setResultName(typeName, "");
518 
519  auto faceBlendedPtr = surfaceScalarField::New
520  (
521  resultName_,
523  mesh_,
525  );
526  store(resultName_, faceBlendedPtr);
527 
528  const auto* nonOrthPtr =
529  mesh_.findObject<volScalarField>(nonOrthogonalityName_);
530 
531  if (nonOrthogonality_ && !nonOrthPtr)
532  {
533  IOobject fieldHeader
534  (
535  nonOrthogonalityName_,
536  mesh_.time().constant(),
537  mesh_.thisDb(),
541  );
542 
543  if (fieldHeader.typeHeaderOk<volScalarField>(true, true, false))
544  {
545  auto* vfPtr = new volScalarField(fieldHeader, mesh_);
546  regIOobject::store(vfPtr);
547  }
548  else
549  {
551  << "Field : " << nonOrthogonalityName_ << " not found."
552  << " The function object will not be used"
553  << exit(FatalError);
554  }
555  }
556 
557 
558  const auto* faceWeightsPtr =
559  mesh_.findObject<volScalarField>(faceWeightName_);
560 
561  if (faceWeight_ && !faceWeightsPtr)
562  {
563  IOobject fieldHeader
564  (
565  faceWeightName_,
566  mesh_.time().constant(),
567  mesh_.thisDb(),
571  );
572 
573  if (fieldHeader.typeHeaderOk<volScalarField>(true, true, false))
574  {
575  auto* vfPtr = new volScalarField(fieldHeader, mesh_);
576  regIOobject::store(vfPtr);
577  }
578  else
579  {
581  << "Field : " << faceWeightName_ << " not found."
582  << " The function object will not be used"
583  << exit(FatalError);
584  }
585  }
586 
587  const auto* skewnessPtr = mesh_.findObject<volScalarField>(skewnessName_);
588 
589  if (skewness_ && !skewnessPtr)
590  {
591  IOobject fieldHeader
592  (
593  skewnessName_,
594  mesh_.time().constant(),
595  mesh_.thisDb(),
599  );
600 
601  if (fieldHeader.typeHeaderOk<volScalarField>(true, true, false))
602  {
603  auto* vfPtr = new volScalarField(fieldHeader, mesh_);
604  regIOobject::store(vfPtr);
605  }
606  else
607  {
609  << "Field : " << skewnessName_ << " not found."
610  << " The function object will not be used"
611  << exit(FatalError);
612  }
613  }
614 
615  if (log)
616  {
617  indicator().writeOpt(IOobject::AUTO_WRITE);
618  }
619 
620  if (writeToFile_)
621  {
623  }
624 
625  init(true);
626 }
627 
628 
629 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
630 
632 (
633  const dictionary& dict
634 )
635 {
637  {
638  dict.readEntry("switchNonOrtho", nonOrthogonality_);
639  dict.readEntry("switchGradCc", gradCc_);
640  dict.readEntry("switchResiduals", residuals_);
641  dict.readEntry("switchFaceWeight", faceWeight_);
642  dict.readEntry("switchSkewness", skewness_);
643  dict.readEntry("switchCo", Co_);
644 
645  dict.readIfPresent("maxNonOrthogonality", maxNonOrthogonality_);
646  dict.readIfPresent("maxGradCc", maxGradCc_);
647  dict.readIfPresent("maxResidual", maxResidual_);
648  dict.readIfPresent("maxSkewness", maxSkewness_);
649  dict.readIfPresent("maxFaceWeight", maxFaceWeight_);
650  dict.readIfPresent("Co2", Co2_);
651 
652  dict.readIfPresent("minFaceWeight", minFaceWeight_);
653  dict.readIfPresent("minNonOrthogonality", minNonOrthogonality_);
654  dict.readIfPresent("minGradCc", minGradCc_);
655  dict.readIfPresent("minSkewness", minSkewness_);
656  dict.readIfPresent("Co1", Co1_);
657 
658 
659  dict.readIfPresent("P", P_);
660  dict.readIfPresent("I", I_);
661  dict.readIfPresent("D", D_);
662 
663  tolerance_ = 0.001;
664  if
665  (
666  dict.readIfPresent("tolerance", tolerance_)
667  && (tolerance_ < 0 || tolerance_ > 1)
668  )
669  {
671  << "tolerance must be in the range 0 to 1. Supplied value: "
672  << tolerance_ << exit(FatalError);
673  }
674 
675  Info<< type() << " " << name() << ":" << nl;
676  if (nonOrthogonality_)
677  {
678  Info<< " Including nonOrthogonality between: "
679  << minNonOrthogonality_ << " and " << maxNonOrthogonality_
680  << endl;
681  }
682  if (gradCc_)
683  {
684  Info<< " Including gradient between: "
685  << minGradCc_ << " and " << maxGradCc_ << endl;
686  }
687  if (residuals_)
688  {
689  Info<< " Including residuals" << endl;
690  }
691  if (faceWeight_)
692  {
693  Info<< " Including faceWeight between: "
694  << minFaceWeight_ << " and " << maxFaceWeight_ << endl;
695  }
696  if (skewness_)
697  {
698  Info<< " Including skewness between: "
699  << minSkewness_ << " and " << maxSkewness_ << endl;
700  }
701  if (Co_)
702  {
703  Info<< " Including Co between: "
704  << Co2_ << " and " << Co1_ << endl;
705  }
706 
707  return true;
708  }
709 
710  return false;
711 }
712 
713 
715 {
716  label nCellsScheme1 = 0;
717  label nCellsScheme2 = 0;
718  label nCellsBlended = 0;
719 
720  calcStats(nCellsScheme1, nCellsScheme2, nCellsBlended);
721 
722  if (writeToFile_)
723  {
724  writeCurrentTime(file());
725 
726  file()
727  << tab << nCellsScheme1
728  << tab << nCellsScheme2
729  << tab << nCellsBlended
730  << endl;
731  }
732 
733  return true;
734 }
735 
736 
737 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
defineTypeNameAndDebug(ObukhovLength, 0)
uint8_t direction
Definition: direction.H:46
word resultName_
Name of result field.
static void writeHeader(Ostream &os, const word &fieldName)
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
word fieldName_
Name of field to process.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool read(const dictionary &dict)
Read the fieldExpression data.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
virtual void writeFileHeader(Ostream &os) const
Write the file header.
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition: fvMesh.H:376
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
void setResultName(const word &typeName, const word &defaultArg)
Set the name of result field.
virtual bool read(const dictionary &)
Read the stabilityBlendingFactor data.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar cbrt(const dimensionedScalar &ds)
bool store(word &fieldName, const tmp< ObjectType > &tfield, bool cacheable=false)
Store the field in the (sub) objectRegistry under the given name.
stabilityBlendingFactor(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
bool log
Flag to write log into Info.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
bool writeToFile_
Flag to enable/disable writing to file.
Definition: writeFile.H:146
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
Intermediate class for handling field expression function objects (e.g. blendingFactor etc...
Nothing to be read.
Automatically write from objectRegistry::writeObject()
#define Log
Definition: PDRblock.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool write()
Write the stabilityBlendingFactor.
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Base class for writing single files from the function objects.
Definition: writeFile.H:112
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Do not request registration (bool: false)
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const Time & time_
Reference to the time database.