regionSizeDistribution.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2013-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "regionSizeDistribution.H"
27 #include "fvcVolumeIntegrate.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34  namespace functionObjects
35  {
37 
39  (
43  );
44  }
45 
46  //- Plus op for FixedList<scalar>
47  template<class T, unsigned Size>
49  {
50  public:
51  void operator()
52  (
54  const FixedList<T, Size>& y
55  ) const
56  {
57  forAll(x, i)
58  {
59  x[i] += y[i];
60  }
61  }
62  };
63 }
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
68 void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
69 (
70  const regionSplit& regions,
71  const Map<label>& patchRegions,
72  const Map<scalar>& regionVolume,
73  const volScalarField& alpha
74 ) const
75 {
76  const scalar maxDropletVol = 1.0/6.0*pow(maxDiam_, 3);
77 
78  // Split alpha field
79  // ~~~~~~~~~~~~~~~~~
80  // Split into
81  // - liquidCore : region connected to inlet patches
82  // - per region a volume : for all other regions
83  // - backgroundAlpha : remaining alpha
84 
85 
86  // Construct field
87  volScalarField liquidCore
88  (
89  IOobject
90  (
91  alphaName_ + "_liquidCore",
92  time_.name(),
93  obr_,
95  ),
96  alpha,
98  );
99 
100  volScalarField backgroundAlpha
101  (
102  IOobject
103  (
104  alphaName_ + "_background",
105  time_.name(),
106  obr_,
108  ),
109  alpha,
111  );
112 
113 
114  // Knock out any cell not in patchRegions
115  forAll(liquidCore, celli)
116  {
117  label regionI = regions[celli];
118  if (patchRegions.found(regionI))
119  {
120  backgroundAlpha[celli] = 0;
121  }
122  else
123  {
124  liquidCore[celli] = 0;
125 
126  scalar regionVol = regionVolume[regionI];
127  if (regionVol < maxDropletVol)
128  {
129  backgroundAlpha[celli] = 0;
130  }
131  }
132  }
133  liquidCore.correctBoundaryConditions();
134  backgroundAlpha.correctBoundaryConditions();
135 
136  Info<< " Volume of liquid-core = "
137  << fvc::domainIntegrate(liquidCore).value()
138  << endl;
139  Info<< " Volume of background = "
140  << fvc::domainIntegrate(backgroundAlpha).value()
141  << endl;
142 
143  Info<< " Writing liquid-core field to " << liquidCore.name() << endl;
144  liquidCore.write();
145  Info<< " Writing background field to " << backgroundAlpha.name() << endl;
146  backgroundAlpha.write();
147 }
148 
149 
151 Foam::functionObjects::regionSizeDistribution::findPatchRegions
152 (
153  const regionSplit& regions
154 ) const
155 {
156  // Mark all regions starting at patches
157  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
158 
159  // Count number of patch faces (just for initial sizing)
160  const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
161 
162  label nPatchFaces = 0;
163  forAllConstIter(labelHashSet, patchIDs, iter)
164  {
165  nPatchFaces += mesh_.boundaryMesh()[iter.key()].size();
166  }
167 
168 
169  Map<label> patchRegions(nPatchFaces);
170  forAllConstIter(labelHashSet, patchIDs, iter)
171  {
172  const polyPatch& pp = mesh_.boundaryMesh()[iter.key()];
173 
174  // Collect all regions on the patch
175  const labelList& faceCells = pp.faceCells();
176 
177  forAll(faceCells, i)
178  {
179  patchRegions.insert
180  (
181  regions[faceCells[i]],
182  Pstream::myProcNo() // dummy value
183  );
184  }
185  }
186 
187 
188  // Make sure all the processors have the same set of regions
189  Pstream::mapCombineGather(patchRegions, minEqOp<label>());
190  Pstream::mapCombineScatter(patchRegions);
191 
192  return patchRegions;
193 }
194 
195 
196 template<class Type>
198 Foam::functionObjects::regionSizeDistribution::divide
199 (
200  const Field<Type>& num,
201  const scalarField& denom
202 )
203 {
204  tmp<Field<Type>> tresult(new Field<Type>(num.size()));
205  Field<Type>& result = tresult.ref();
206 
207  forAll(denom, i)
208  {
209  if (denom[i] != 0)
210  {
211  result[i] = num[i]/denom[i];
212  }
213  else
214  {
215  result[i] = Zero;
216  }
217  }
218 
219  return tresult;
220 }
221 
222 
223 template<class Type>
224 void Foam::functionObjects::regionSizeDistribution::generateFields
225 (
226  const word& fieldName, // name of field
227  const labelList& indices, // index of bin for each region
228  const Field<Type>& sortedField, // per region field data
229  const scalarField& binCount, // per bin number of regions
232 ) const
233 {
234  if (Pstream::master())
235  {
236  // Calculate per-bin sum
237  Field<Type> binSum(nBins_, Zero);
238  forAll(sortedField, i)
239  {
240  binSum[indices[i]] += sortedField[i];
241  }
242 
243  // Calculate per-bin average
244  Field<Type> binAvg(divide(binSum, binCount));
245 
246  // Append
247  fields.setSize(fieldNames.size());
248  fieldNames.append(fieldName + "_sum");
249  fields.append(binSum);
250  fieldNames.append(fieldName + "_avg");
251  fields.append(binAvg);
252  }
253 }
254 
255 
256 void Foam::functionObjects::regionSizeDistribution::generateFields
257 (
258  const word& fieldName, // name of field
259  const labelList& indices, // index of bin for each region
260  const scalarField& sortedField, // per region field data
261  const scalarField& binCount, // per bin number of regions
263  PtrList<scalarField>& fields
264 ) const
265 {
266  if (Pstream::master())
267  {
268  // Calculate per-bin sum
269  scalarField binSum(nBins_, Zero);
270  forAll(sortedField, i)
271  {
272  binSum[indices[i]] += sortedField[i];
273  }
274 
275  // Calculate per-bin average
276  scalarField binAvg(divide(binSum, binCount));
277 
278  // Calculate per-bin deviation
279  scalarField binSqrSum(nBins_, Zero);
280  forAll(sortedField, i)
281  {
282  binSqrSum[indices[i]] += Foam::sqr(sortedField[i]);
283  }
284  scalarField binDev
285  (
286  sqrt(divide(binSqrSum, binCount) - Foam::sqr(binAvg))
287  );
288 
289  // Append
290  fields.setSize(fieldNames.size());
291  fieldNames.append(fieldName + "_sum");
292  fields.append(binSum);
293  fieldNames.append(fieldName + "_avg");
294  fields.append(binAvg);
295  fieldNames.append(fieldName + "_dev");
296  fields.append(binDev);
297  }
298 }
299 
300 
301 template<class Type>
302 void Foam::functionObjects::regionSizeDistribution::generateFields
303 (
304  const word& fieldName, // name of field
305  const Field<Type>& cellField, // per cell field data
306  const regionSplit& regions, // per cell the region(=droplet)
307  const labelList& sortedRegions, // valid regions in sorted order
308  const scalarField& sortedNormalisation,
309  const labelList& indices, // index of bin for each region
310  const scalarField& binCount, // per bin number of regions
312  PtrList<Field<Type>>& fields
313 ) const
314 {
315  // Sum on a per-region basis. Parallel reduced.
316  Map<Type> regionField(regionSum(regions, cellField));
317 
318  // Extract in region order
319  Field<Type> sortedField
320  (
321  sortedNormalisation*extractData(sortedRegions, regionField)
322  );
323 
324  // Generate fields
325  generateFields
326  (
327  fieldName, // name of field
328  indices, // index of bin for each region
329  sortedField, // per region field data
330  binCount, // per bin number of regions
331  fieldNames,
332  fields
333  );
334 }
335 
336 
337 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
338 
340 (
341  const word& name,
342  const Time& runTime,
343  const dictionary& dict
344 )
345 :
346  fvMeshFunctionObject(name, runTime, dict),
347  file_(obr_, name),
348  alphaName_(dict.lookup("alpha")),
349  patchNames_(dict.lookup("patches"))
350 {
351  read(dict);
352 }
353 
354 
355 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
356 
358 {}
359 
360 
361 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
362 
364 {
365  dict.lookup("alpha") >> alphaName_;
366  dict.lookup("patches") >> patchNames_;
367  dict.lookup("threshold") >> threshold_;
368  dict.lookup("maxDiameter") >> maxDiam_;
369  minDiam_ = 0.0;
370  dict.readIfPresent("minDiameter", minDiam_);
371  dict.lookup("nBins") >> nBins_;
372  dict.lookup("fields") >> fields_;
373 
374  formatterPtr_ = setWriter::New(dict.lookup("setFormat"), dict);
375 
376  return true;
377 }
378 
379 
381 {
382  wordList fields(fields_);
383  fields.append(alphaName_);
384  return fields;
385 }
386 
387 
389 {
390  return true;
391 }
392 
393 
395 {
396  Info<< type() << " " << name() << " write:" << nl;
397 
398  autoPtr<volScalarField> alphaPtr;
399  if (obr_.foundObject<volScalarField>(alphaName_))
400  {
401  Info<< " Looking up field " << alphaName_ << endl;
402  }
403  else
404  {
405  Info<< " Reading field " << alphaName_ << endl;
406  alphaPtr.reset
407  (
408  new volScalarField
409  (
410  IOobject
411  (
412  alphaName_,
413  time_.name(),
414  mesh_,
417  ),
418  mesh_
419  )
420  );
421  }
422 
423 
424  const volScalarField& alpha =
425  (
426  alphaPtr.valid()
427  ? alphaPtr()
428  : obr_.lookupObject<volScalarField>(alphaName_)
429  );
430 
431  Info<< " Volume of alpha = "
432  << fvc::domainIntegrate(alpha).value()
433  << endl;
434 
435  const scalar meshVol = gSum(mesh_.V());
436  const scalar maxDropletVol = 1.0/6.0*pow(maxDiam_, 3);
437  const scalar delta = (maxDiam_-minDiam_)/nBins_;
438 
439  Info<< " Mesh volume = " << meshVol << endl;
440  Info<< " Maximum droplet diameter = " << maxDiam_ << endl;
441  Info<< " Maximum droplet volume = " << maxDropletVol << endl;
442 
443 
444  // Determine blocked faces
445  boolList blockedFace(mesh_.nFaces(), false);
446 
447  {
448  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
449  {
450  scalar ownVal = alpha[mesh_.faceOwner()[facei]];
451  scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
452 
453  if
454  (
455  (ownVal < threshold_ && neiVal > threshold_)
456  || (ownVal > threshold_ && neiVal < threshold_)
457  )
458  {
459  blockedFace[facei] = true;
460  }
461  }
462 
463  // Block coupled faces
464  forAll(alpha.boundaryField(), patchi)
465  {
466  const fvPatchScalarField& fvp = alpha.boundaryField()[patchi];
467  if (fvp.coupled())
468  {
469  tmp<scalarField> townFld(fvp.patchInternalField());
470  const scalarField& ownFld = townFld();
471  tmp<scalarField> tnbrFld(fvp.patchNeighbourField());
472  const scalarField& nbrFld = tnbrFld();
473 
474  label start = fvp.patch().patch().start();
475 
476  forAll(ownFld, i)
477  {
478  scalar ownVal = ownFld[i];
479  scalar neiVal = nbrFld[i];
480 
481  if
482  (
483  (ownVal < threshold_ && neiVal > threshold_)
484  || (ownVal > threshold_ && neiVal < threshold_)
485  )
486  {
487  blockedFace[start+i] = true;
488  }
489  }
490  }
491  }
492  }
493 
494 
495  regionSplit regions(mesh_, blockedFace);
496 
497  Info<< " Determined " << regions.nRegions()
498  << " disconnected regions" << endl;
499 
500 
501  if (debug)
502  {
503  volScalarField region
504  (
505  IOobject
506  (
507  "region",
508  time_.name(),
509  mesh_,
512  ),
513  mesh_,
515  );
516  Info<< " Dumping region as volScalarField to " << region.name()
517  << endl;
518 
519  forAll(regions, celli)
520  {
521  region[celli] = regions[celli];
522  }
523  region.correctBoundaryConditions();
524  region.write();
525  }
526 
527 
528  // Determine regions connected to supplied patches
529  Map<label> patchRegions(findPatchRegions(regions));
530 
531 
532  // Sum all regions
533  const scalarField alphaVol(alpha.primitiveField()*mesh_.V());
534  Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
535  Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
536  Map<label> allRegionNumCells
537  (
538  regionSum
539  (
540  regions,
541  labelField(mesh_.nCells(), 1.0)
542  )
543  );
544 
545  if (debug)
546  {
547  Info<< " " << tab << "Region"
548  << tab << "Volume(mesh)"
549  << tab << "Volume(" << alpha.name() << "):"
550  << tab << "nCells"
551  << endl;
552  scalar meshSumVol = 0.0;
553  scalar alphaSumVol = 0.0;
554  label nCells = 0;
555 
556  Map<scalar>::const_iterator vIter = allRegionVolume.begin();
557  Map<scalar>::const_iterator aIter = allRegionAlphaVolume.begin();
558  Map<label>::const_iterator numIter = allRegionNumCells.begin();
559  for
560  (
561  ;
562  vIter != allRegionVolume.end()
563  && aIter != allRegionAlphaVolume.end();
564  ++vIter, ++aIter, ++numIter
565  )
566  {
567  Info<< " " << tab << vIter.key()
568  << tab << vIter()
569  << tab << aIter()
570  << tab << numIter()
571  << endl;
572 
573  meshSumVol += vIter();
574  alphaSumVol += aIter();
575  nCells += numIter();
576  }
577  Info<< " " << tab << "Total:"
578  << tab << meshSumVol
579  << tab << alphaSumVol
580  << tab << nCells
581  << endl;
582  Info<< endl;
583  }
584 
585 
586 
587 
588  {
589  Info<< " Patch connected regions (liquid core):" << endl;
590  Info<< tab << " Region"
591  << tab << "Volume(mesh)"
592  << tab << "Volume(" << alpha.name() << "):"
593  << endl;
594  forAllConstIter(Map<label>, patchRegions, iter)
595  {
596  label regionI = iter.key();
597  Info<< " " << tab << iter.key()
598  << tab << allRegionVolume[regionI]
599  << tab << allRegionAlphaVolume[regionI] << endl;
600 
601  }
602  Info<< endl;
603  }
604 
605  {
606  Info<< " Background regions:" << endl;
607  Info<< " " << tab << "Region"
608  << tab << "Volume(mesh)"
609  << tab << "Volume(" << alpha.name() << "):"
610  << endl;
611  Map<scalar>::const_iterator vIter = allRegionVolume.begin();
612  Map<scalar>::const_iterator aIter = allRegionAlphaVolume.begin();
613 
614  for
615  (
616  ;
617  vIter != allRegionVolume.end()
618  && aIter != allRegionAlphaVolume.end();
619  ++vIter, ++aIter
620  )
621  {
622  if
623  (
624  !patchRegions.found(vIter.key())
625  && vIter() >= maxDropletVol
626  )
627  {
628  Info<< " " << tab << vIter.key()
629  << tab << vIter()
630  << tab << aIter() << endl;
631  }
632  }
633  Info<< endl;
634  }
635 
636 
637 
638  // Split alpha field
639  // ~~~~~~~~~~~~~~~~~
640  // Split into
641  // - liquidCore : region connected to inlet patches
642  // - per region a volume : for all other regions
643  // - backgroundAlpha : remaining alpha
644  writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
645 
646 
647  // Extract droplet-only allRegionVolume, i.e. delete liquid core
648  // (patchRegions) and background regions from maps.
649  // Note that we have to use mesh volume (allRegionVolume) and not
650  // allRegionAlphaVolume since background might not have alpha in it.
651  forAllIter(Map<scalar>, allRegionVolume, vIter)
652  {
653  label regionI = vIter.key();
654  if
655  (
656  patchRegions.found(regionI)
657  || vIter() >= maxDropletVol
658  )
659  {
660  allRegionVolume.erase(vIter);
661  allRegionAlphaVolume.erase(regionI);
662  allRegionNumCells.erase(regionI);
663  }
664  }
665 
666  if (allRegionVolume.size())
667  {
668  // Construct mids of bins for plotting
669  scalarField xBin(nBins_);
670  scalar x = 0.5*delta;
671  forAll(xBin, i)
672  {
673  xBin[i] = x;
674  x += delta;
675  }
676 
677  // Get in region order the alpha*volume and diameter
678  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
679 
680  const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
681 
682  scalarField sortedVols
683  (
684  extractData
685  (
686  sortedRegions,
687  allRegionAlphaVolume
688  )
689  );
690 
691  // Calculate the diameters
692  scalarField sortedDiameters(sortedVols.size());
693  forAll(sortedDiameters, i)
694  {
695  sortedDiameters[i] = Foam::cbrt
696  (
697  sortedVols[i]
699  );
700  }
701 
702  // Determine the bin index for all the diameters
703  labelList indices(sortedDiameters.size());
704  forAll(sortedDiameters, i)
705  {
706  indices[i] = (sortedDiameters[i]-minDiam_)/delta;
707  }
708 
709  // Calculate the counts per diameter bin
710  scalarField binCount(nBins_, 0.0);
711  forAll(sortedDiameters, i)
712  {
713  binCount[indices[i]] += 1.0;
714  }
715 
716  // Write to log
717  {
718  Info<< " Bins:" << endl;
719  Info<< " " << tab << "Bin"
720  << tab << "Min diameter"
721  << tab << "Count:"
722  << endl;
723 
724  scalar diam = 0.0;
725  forAll(binCount, binI)
726  {
727  Info<< " " << tab << binI
728  << tab << diam
729  << tab << binCount[binI] << endl;
730  diam += delta;
731  }
732  Info<< endl;
733  }
734 
735  // Declare fields and field names
737  #define DeclareTypeFields(Type, nullArg) \
738  PtrList<Field<Type>> Type##Fields;
740  #undef DeclareTypeFields
741 
742  // Add the bin count
743  fieldNames.append("binCount");
744  #define TypeFieldsAppend(Type, nullArg) \
745  appendFields(binCount, Type##Fields);
746  #undef TypeFieldsAppend
747 
748  // Add the volumes
749  generateFields
750  (
751  "volume",
752  indices,
753  sortedVols,
754  binCount,
755  fieldNames,
756  scalarFields
757  );
758 
759  // Add other sampled fields
760  forAll(fields_, fieldi)
761  {
762  bool found = false;
763 
764  #define GenerateTypeFields(Type, nullArg) \
765  \
766  if (obr_.foundObject<VolField<Type>>(fields_[fieldi])) \
767  { \
768  found = true; \
769  \
770  const VolField<Type>& field = \
771  obr_.lookupObject<VolField<Type>>(fields_[fieldi]); \
772  \
773  generateFields \
774  ( \
775  fields_[fieldi], \
776  (alphaVol*field)(), \
777  regions, \
778  sortedRegions, \
779  1.0/sortedVols, \
780  indices, \
781  binCount, \
782  fieldNames, \
783  Type##Fields \
784  ); \
785  }
787  #undef GenerateTypeFields
788 
789  if (!found) cannotFindObject(fields_[fieldi]);
790  }
791 
792  // Expand all field lists
793  #define TypeFieldsExpand(Type, nullArg) \
794  Type##Fields.setSize(fieldNames.size());
796  #undef TypeFieldsAppend
797 
798  // Write
799  formatterPtr_().write
800  (
801  file_.baseTimeDir(),
802  typeName,
803  coordSet(true, "diameter", xBin),
804  fieldNames
805  #define TypeFieldsParameter(Type, nullArg) , Type##Fields
807  #undef TypeFieldsParameter
808  );
809  }
810 
811  return true;
812 }
813 
814 
815 // ************************************************************************* //
scalar y
scalar delta
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Pre-declare SubField and related Field type.
Definition: Field.H:82
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:78
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
An STL-conforming const_iterator.
Definition: HashTable.H:484
const Key & key() const
Return the Key corresponding to the iterator.
Definition: HashTableI.H:254
List< Key > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:371
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
Plus op for FixedList<scalar>
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A HashTable to objects of type <T> with a label key.
Definition: Map.H:52
virtual Ostream & write(const char)=0
Write character.
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
Holds list of sampling positions.
Definition: coordSet.H:51
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const word & name() const
Return const reference to name.
Abstract base-class for Time/database functionObjects.
const Time & time_
Reference to time.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const objectRegistry & obr_
Reference to the region objectRegistry.
Creates a size distribution via interrogating a continuous phase fraction field.
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
virtual wordList fields() const
Return the list of fields required.
virtual bool write()
Calculate the regionSizeDistribution and write.
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:327
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
virtual tmp< Field< Type > > patchNeighbourField(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:436
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:355
static const word & calculatedType()
Return the type of the calculated for of fvPatchField.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:139
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
virtual bool write(const bool write=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:118
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:208
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
Definition: setWriter.C:286
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
Volume integrate volField creating a volField.
label patchi
static List< word > fieldNames
Definition: globalFoam.H:46
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
List< word > wordList
A List of words.
Definition: fileName.H:54
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static const char tab
Definition: Ostream.H:259
const dimensionSet dimless
messageStream Info
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define GenerateTypeFields(Type, nullArg)
#define TypeFieldsParameter(Type, nullArg)
#define DeclareTypeFields(Type, nullArg)
#define TypeFieldsExpand(Type, nullArg)
dictionary dict
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112