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-2023 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_,
72  ),
73  mesh_,
76  );
77 
78  mesh_.objectRegistry::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 = tmp<volScalarField>::New
325  (
326  IOobject
327  (
328  "magGradCC",
329  time_.timeName(),
330  mesh_,
333  ),
334  mesh_,
337  );
338  auto& magGradCC = tmagGradCC.ref();
339 
340  for (direction i=0; i<vector::nComponents; i++)
341  {
342  // Create field with zero grad
343  volScalarField cci
344  (
345  IOobject
346  (
347  "cc" + word(vector::componentNames[i]),
348  mesh_.time().timeName(),
349  mesh_,
352  ),
353  mesh_,
356  );
357  cci = mesh_.C().component(i);
358  cci.correctBoundaryConditions();
359  magGradCC += mag(fvc::grad(cci)).ref();
360  }
361 
362  if (first)
363  {
364  Log << " Max magGradCc : " << max(magGradCC.ref()).value()
365  << endl;
366  }
367 
368  indicator =
369  max
370  (
371  indicator,
372  min
373  (
374  max
375  (
376  scalar(0),
377  (magGradCC - maxGradCc_)
378  / (minGradCc_ - maxGradCc_)
379  ),
380  scalar(1)
381  )
382  );
383  }
384 
385 
386  const auto* UNamePtr = mesh_.findObject<volVectorField>(UName_);
387 
388  if (Co_)
389  {
390  if (Co1_ >= Co2_)
391  {
393  << "Co2 should be larger than Co1."
394  << exit(FatalError);
395  }
396 
397  auto CoPtr = tmp<volScalarField>::New
398  (
399  IOobject
400  (
401  "Co",
402  time_.timeName(),
403  mesh_,
406  ),
407  mesh_,
410  );
411 
412  auto& Co = CoPtr.ref();
413 
414  Co.primitiveFieldRef() =
415  mesh_.time().deltaT()*mag(*UNamePtr)/cbrt(mesh_.V());
416 
417  indicator =
418  max
419  (
420  indicator,
421  clamp((Co - Co1_)/(Co2_ - Co1_), zero_one{})
422  );
423 
424  if (first)
425  {
426  Log << " Max Co : " << max(Co).value()
427  << endl;
428  }
429  }
430 
431  if (first)
432  {
433  Log << nl;
434  }
435 
436  if (log)
437  {
438  label nCellsScheme1 = 0;
439  label nCellsScheme2 = 0;
440  label nCellsBlended = 0;
441 
442  calcStats(nCellsScheme1, nCellsScheme2, nCellsBlended);
443 
444  Log << nl << name() << " execute :" << nl
445  << " scheme 1 cells : " << nCellsScheme1 << nl
446  << " scheme 2 cells : " << nCellsScheme2 << nl
447  << " blended cells : " << nCellsBlended << nl
448  << endl;
449  }
450 
451  indicator.correctBoundaryConditions();
452  indicator.clamp_range(zero_one{});
453 
454  // Update the blended surface field
455  auto& surBlended = mesh_.lookupObjectRef<surfaceScalarField>(resultName_);
456 
457  surBlended = fvc::interpolate(indicator);
458 
459  return true;
460 }
461 
462 
463 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
464 
466 (
467  const word& name,
468  const Time& runTime,
469  const dictionary& dict
470 )
471 :
473  writeFile(obr_, name, typeName, dict),
474  nonOrthogonality_(dict.getOrDefault<Switch>("switchNonOrtho", false)),
475  gradCc_(dict.getOrDefault<Switch>("switchGradCc", false)),
476  residuals_(dict.getOrDefault<Switch>("switchResiduals", false)),
477  faceWeight_(dict.getOrDefault<Switch>("switchFaceWeight", false)),
478  skewness_(dict.getOrDefault<Switch>("switchSkewness", false)),
479  Co_(dict.getOrDefault<Switch>("switchCo", false)),
480 
481  maxNonOrthogonality_
482  (
483  dict.getOrDefault<scalar>("maxNonOrthogonality", 20.0)
484  ),
485  minNonOrthogonality_
486  (
487  dict.getOrDefault<scalar>("minNonOrthogonality", 60.0)
488  ),
489  maxGradCc_(dict.getOrDefault<scalar>("maxGradCc", 3.0)),
490  minGradCc_(dict.getOrDefault<scalar>("minGradCc", 4.0)),
491  maxResidual_(dict.getOrDefault<scalar>("maxResidual", 10.0)),
492  minFaceWeight_(dict.getOrDefault<scalar>("minFaceWeight", 0.3)),
493  maxFaceWeight_(dict.getOrDefault<scalar>("maxFaceWeight", 0.2)),
494  maxSkewness_(dict.getOrDefault<scalar>("maxSkewness", 2.0)),
495  minSkewness_(dict.getOrDefault<scalar>("minSkewness", 3.0)),
496  Co1_(dict.getOrDefault<scalar>("Co1", 1.0)),
497  Co2_(dict.getOrDefault<scalar>("Co2", 10.0)),
498 
499  nonOrthogonalityName_
500  (
501  dict.getOrDefault<word>("nonOrthogonality", "nonOrthoAngle")
502  ),
503  faceWeightName_
504  (
505  dict.getOrDefault<word>("faceWeight", "faceWeight")
506  ),
507  skewnessName_
508  (
509  dict.getOrDefault<word>("skewness", "skewness")
510  ),
511  residualName_
512  (
513  dict.getOrDefault<word>
514  (
515  "residual",
516  IOobject::scopedName("initialResidual", "p")
517  )
518  ),
519  UName_
520  (
521  dict.getOrDefault<word>("U", "U")
522  ),
523 
524  tolerance_(0.001),
525  error_(mesh_.nCells(), Zero),
526  errorIntegral_(mesh_.nCells(), Zero),
527  oldError_(mesh_.nCells(), Zero),
528  oldErrorIntegral_(mesh_.nCells(), Zero),
529  P_(dict.getOrDefault<scalar>("P", 3)),
530  I_(dict.getOrDefault<scalar>("I", 0.0)),
531  D_(dict.getOrDefault<scalar>("D", 0.25))
532 {
533  read(dict);
534  setResultName(typeName, "");
535 
536  auto faceBlendedPtr = tmp<surfaceScalarField>::New
537  (
538  IOobject
539  (
540  resultName_,
541  time_.timeName(),
542  mesh_,
546  ),
547  mesh_,
549  );
550  store(resultName_, faceBlendedPtr);
551 
552  const auto* nonOrthPtr =
553  mesh_.findObject<volScalarField>(nonOrthogonalityName_);
554 
555  if (nonOrthogonality_ && !nonOrthPtr)
556  {
557  IOobject fieldHeader
558  (
559  nonOrthogonalityName_,
560  mesh_.time().constant(),
561  mesh_,
565  );
566 
567  if (fieldHeader.typeHeaderOk<volScalarField>(true, true, false))
568  {
569  auto* vfPtr = new volScalarField(fieldHeader, mesh_);
570  mesh_.objectRegistry::store(vfPtr);
571  }
572  else
573  {
575  << "Field : " << nonOrthogonalityName_ << " not found."
576  << " The function object will not be used"
577  << exit(FatalError);
578  }
579  }
580 
581 
582  const auto* faceWeightsPtr =
583  mesh_.findObject<volScalarField>(faceWeightName_);
584 
585  if (faceWeight_ && !faceWeightsPtr)
586  {
587  IOobject fieldHeader
588  (
589  faceWeightName_,
590  mesh_.time().constant(),
591  mesh_,
595  );
596 
597  if (fieldHeader.typeHeaderOk<volScalarField>(true, true, false))
598  {
599  auto* vfPtr = new volScalarField(fieldHeader, mesh_);
600  mesh_.objectRegistry::store(vfPtr);
601  }
602  else
603  {
605  << "Field : " << faceWeightName_ << " not found."
606  << " The function object will not be used"
607  << exit(FatalError);
608  }
609  }
610 
611  const auto* skewnessPtr = mesh_.findObject<volScalarField>(skewnessName_);
612 
613  if (skewness_ && !skewnessPtr)
614  {
615  IOobject fieldHeader
616  (
617  skewnessName_,
618  mesh_.time().constant(),
619  mesh_,
623  );
624 
625  if (fieldHeader.typeHeaderOk<volScalarField>(true, true, false))
626  {
627  auto* vfPtr = new volScalarField(fieldHeader, mesh_);
628  mesh_.objectRegistry::store(vfPtr);
629  }
630  else
631  {
633  << "Field : " << skewnessName_ << " not found."
634  << " The function object will not be used"
635  << exit(FatalError);
636  }
637  }
638 
639  if (log)
640  {
641  indicator().writeOpt(IOobject::AUTO_WRITE);
642  }
643 
644  if (writeToFile_)
645  {
647  }
648 
649  init(true);
650 }
651 
652 
653 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
654 
656 (
657  const dictionary& dict
658 )
659 {
661  {
662  dict.readEntry("switchNonOrtho", nonOrthogonality_);
663  dict.readEntry("switchGradCc", gradCc_);
664  dict.readEntry("switchResiduals", residuals_);
665  dict.readEntry("switchFaceWeight", faceWeight_);
666  dict.readEntry("switchSkewness", skewness_);
667  dict.readEntry("switchCo", Co_);
668 
669  dict.readIfPresent("maxNonOrthogonality", maxNonOrthogonality_);
670  dict.readIfPresent("maxGradCc", maxGradCc_);
671  dict.readIfPresent("maxResidual", maxResidual_);
672  dict.readIfPresent("maxSkewness", maxSkewness_);
673  dict.readIfPresent("maxFaceWeight", maxFaceWeight_);
674  dict.readIfPresent("Co2", Co2_);
675 
676  dict.readIfPresent("minFaceWeight", minFaceWeight_);
677  dict.readIfPresent("minNonOrthogonality", minNonOrthogonality_);
678  dict.readIfPresent("minGradCc", minGradCc_);
679  dict.readIfPresent("minSkewness", minSkewness_);
680  dict.readIfPresent("Co1", Co1_);
681 
682 
683  dict.readIfPresent("P", P_);
684  dict.readIfPresent("I", I_);
685  dict.readIfPresent("D", D_);
686 
687  tolerance_ = 0.001;
688  if
689  (
690  dict.readIfPresent("tolerance", tolerance_)
691  && (tolerance_ < 0 || tolerance_ > 1)
692  )
693  {
695  << "tolerance must be in the range 0 to 1. Supplied value: "
696  << tolerance_ << exit(FatalError);
697  }
698 
699  Info<< type() << " " << name() << ":" << nl;
700  if (nonOrthogonality_)
701  {
702  Info<< " Including nonOrthogonality between: "
703  << minNonOrthogonality_ << " and " << maxNonOrthogonality_
704  << endl;
705  }
706  if (gradCc_)
707  {
708  Info<< " Including gradient between: "
709  << minGradCc_ << " and " << maxGradCc_ << endl;
710  }
711  if (residuals_)
712  {
713  Info<< " Including residuals" << endl;
714  }
715  if (faceWeight_)
716  {
717  Info<< " Including faceWeight between: "
718  << minFaceWeight_ << " and " << maxFaceWeight_ << endl;
719  }
720  if (skewness_)
721  {
722  Info<< " Including skewness between: "
723  << minSkewness_ << " and " << maxSkewness_ << endl;
724  }
725  if (Co_)
726  {
727  Info<< " Including Co between: "
728  << Co2_ << " and " << Co1_ << endl;
729  }
730 
731  return true;
732  }
733 
734  return false;
735 }
736 
737 
739 {
740  label nCellsScheme1 = 0;
741  label nCellsScheme2 = 0;
742  label nCellsBlended = 0;
743 
744  calcStats(nCellsScheme1, nCellsScheme2, nCellsBlended);
745 
746  if (writeToFile_)
747  {
748  writeCurrentTime(file());
749 
750  file()
751  << tab << nCellsScheme1
752  << tab << nCellsScheme2
753  << tab << nCellsBlended
754  << endl;
755  }
756 
757  return true;
758 }
759 
760 
761 // ************************************************************************* //
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:598
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
tmp< GeometricField< cmptType, PatchField, GeoMesh > > component(const direction) const
Return a component of the field.
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.
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:82
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.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
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.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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
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
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
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:172
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
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.