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-2018 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  {
36  defineTypeNameAndDebug(regionSizeDistribution, 0);
37 
39  (
40  functionObject,
41  regionSizeDistribution,
42  dictionary
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::writeGraph
69 (
70  const coordSet& coords,
71  const word& valueName,
72  const scalarField& values
73 ) const
74 {
75  const wordList valNames(1, valueName);
76 
77  fileName outputPath = file_.baseTimeDir();
78  mkDir(outputPath);
79 
80  OFstream str(outputPath/formatterPtr_().getFileName(coords, valNames));
81 
82  Info<< " Writing distribution of " << valueName << " to " << str.name()
83  << endl;
84 
85  List<const scalarField*> valPtrs(1);
86  valPtrs[0] = &values;
87  formatterPtr_().write(coords, valNames, valPtrs, str);
88 }
89 
90 
91 void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
92 (
93  const regionSplit& regions,
94  const Map<label>& patchRegions,
95  const Map<scalar>& regionVolume,
96  const volScalarField& alpha
97 ) const
98 {
99  const scalar maxDropletVol = 1.0/6.0*pow(maxDiam_, 3);
100 
101  // Split alpha field
102  // ~~~~~~~~~~~~~~~~~
103  // Split into
104  // - liquidCore : region connected to inlet patches
105  // - per region a volume : for all other regions
106  // - backgroundAlpha : remaining alpha
107 
108 
109  // Construct field
110  volScalarField liquidCore
111  (
112  IOobject
113  (
114  alphaName_ + "_liquidCore",
115  obr_.time().timeName(),
116  obr_,
118  ),
119  alpha,
121  );
122 
123  volScalarField backgroundAlpha
124  (
125  IOobject
126  (
127  alphaName_ + "_background",
128  obr_.time().timeName(),
129  obr_,
131  ),
132  alpha,
134  );
135 
136 
137  // Knock out any cell not in patchRegions
138  forAll(liquidCore, celli)
139  {
140  label regionI = regions[celli];
141  if (patchRegions.found(regionI))
142  {
143  backgroundAlpha[celli] = 0;
144  }
145  else
146  {
147  liquidCore[celli] = 0;
148 
149  scalar regionVol = regionVolume[regionI];
150  if (regionVol < maxDropletVol)
151  {
152  backgroundAlpha[celli] = 0;
153  }
154  }
155  }
156  liquidCore.correctBoundaryConditions();
157  backgroundAlpha.correctBoundaryConditions();
158 
159  Info<< " Volume of liquid-core = "
160  << fvc::domainIntegrate(liquidCore).value()
161  << endl;
162  Info<< " Volume of background = "
163  << fvc::domainIntegrate(backgroundAlpha).value()
164  << endl;
165 
166  Info<< " Writing liquid-core field to " << liquidCore.name() << endl;
167  liquidCore.write();
168  Info<< " Writing background field to " << backgroundAlpha.name() << endl;
169  backgroundAlpha.write();
170 }
171 
172 
174 Foam::functionObjects::regionSizeDistribution::findPatchRegions
175 (
176  const regionSplit& regions
177 ) const
178 {
179  // Mark all regions starting at patches
180  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
181 
182  // Count number of patch faces (just for initial sizing)
183  const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
184 
185  label nPatchFaces = 0;
186  forAllConstIter(labelHashSet, patchIDs, iter)
187  {
188  nPatchFaces += mesh_.boundaryMesh()[iter.key()].size();
189  }
190 
191 
192  Map<label> patchRegions(nPatchFaces);
193  forAllConstIter(labelHashSet, patchIDs, iter)
194  {
195  const polyPatch& pp = mesh_.boundaryMesh()[iter.key()];
196 
197  // Collect all regions on the patch
198  const labelList& faceCells = pp.faceCells();
199 
200  forAll(faceCells, i)
201  {
202  patchRegions.insert
203  (
204  regions[faceCells[i]],
205  Pstream::myProcNo() // dummy value
206  );
207  }
208  }
209 
210 
211  // Make sure all the processors have the same set of regions
212  Pstream::mapCombineGather(patchRegions, minEqOp<label>());
213  Pstream::mapCombineScatter(patchRegions);
214 
215  return patchRegions;
216 }
217 
218 
220 Foam::functionObjects::regionSizeDistribution::divide
221 (
222  const scalarField& num,
223  const scalarField& denom
224 )
225 {
226  tmp<scalarField> tresult(new scalarField(num.size()));
227  scalarField& result = tresult.ref();
228 
229  forAll(denom, i)
230  {
231  if (denom[i] != 0)
232  {
233  result[i] = num[i]/denom[i];
234  }
235  else
236  {
237  result[i] = 0.0;
238  }
239  }
240  return tresult;
241 }
242 
243 
244 void Foam::functionObjects::regionSizeDistribution::writeGraphs
245 (
246  const word& fieldName, // name of field
247  const labelList& indices, // index of bin for each region
248  const scalarField& sortedField, // per region field data
249  const scalarField& binCount, // per bin number of regions
250  const coordSet& coords // graph data for bins
251 ) const
252 {
253  if (Pstream::master())
254  {
255  // Calculate per-bin average
256  scalarField binSum(nBins_, 0.0);
257  forAll(sortedField, i)
258  {
259  binSum[indices[i]] += sortedField[i];
260  }
261 
262  scalarField binAvg(divide(binSum, binCount));
263 
264  // Per bin deviation
265  scalarField binSqrSum(nBins_, 0.0);
266  forAll(sortedField, i)
267  {
268  binSqrSum[indices[i]] += Foam::sqr(sortedField[i]);
269  }
270  scalarField binDev
271  (
272  sqrt(divide(binSqrSum, binCount) - Foam::sqr(binAvg))
273  );
274 
275  // Write average
276  writeGraph(coords, fieldName + "_sum", binSum);
277  // Write average
278  writeGraph(coords, fieldName + "_avg", binAvg);
279  // Write deviation
280  writeGraph(coords, fieldName + "_dev", binDev);
281  }
282 }
283 
284 
285 void Foam::functionObjects::regionSizeDistribution::writeGraphs
286 (
287  const word& fieldName, // name of field
288  const scalarField& cellField, // per cell field data
289  const regionSplit& regions, // per cell the region(=droplet)
290  const labelList& sortedRegions, // valid regions in sorted order
291  const scalarField& sortedNormalisation,
292 
293  const labelList& indices, // per region index of bin
294  const scalarField& binCount, // per bin number of regions
295  const coordSet& coords // graph data for bins
296 ) const
297 {
298  // Sum on a per-region basis. Parallel reduced.
299  Map<scalar> regionField(regionSum(regions, cellField));
300 
301  // Extract in region order
302  scalarField sortedField
303  (
304  sortedNormalisation
305  * extractData
306  (
307  sortedRegions,
308  regionField
309  )
310  );
311 
312  writeGraphs
313  (
314  fieldName, // name of field
315  indices, // index of bin for each region
316  sortedField, // per region field data
317  binCount, // per bin number of regions
318  coords // graph data for bins
319  );
320 }
321 
322 
323 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
324 
326 (
327  const word& name,
328  const Time& runTime,
329  const dictionary& dict
330 )
331 :
332  fvMeshFunctionObject(name, runTime, dict),
333  file_(obr_, name),
334  alphaName_(dict.lookup("field")),
335  patchNames_(dict.lookup("patches"))
336 {
337  read(dict);
338 }
339 
340 
341 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
342 
344 {}
345 
346 
347 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
348 
350 {
351  dict.lookup("field") >> alphaName_;
352  dict.lookup("patches") >> patchNames_;
353  dict.lookup("threshold") >> threshold_;
354  dict.lookup("maxDiameter") >> maxDiam_;
355  minDiam_ = 0.0;
356  dict.readIfPresent("minDiameter", minDiam_);
357  dict.lookup("nBins") >> nBins_;
358  dict.lookup("fields") >> fields_;
359 
360  word format(dict.lookup("setFormat"));
361  formatterPtr_ = writer<scalar>::New(format);
362 
363  if (dict.found("coordinateSystem"))
364  {
365  coordSysPtr_.reset(new coordinateSystem(obr_, dict));
366 
367  Info<< "Transforming all vectorFields with coordinate system "
368  << coordSysPtr_().name() << endl;
369  }
370 
371  return true;
372 }
373 
374 
376 {
377  return true;
378 }
379 
380 
382 {
383  Info<< type() << " " << name() << " write:" << nl;
384 
385  autoPtr<volScalarField> alphaPtr;
386  if (obr_.foundObject<volScalarField>(alphaName_))
387  {
388  Info<< " Looking up field " << alphaName_ << endl;
389  }
390  else
391  {
392  Info<< " Reading field " << alphaName_ << endl;
393  alphaPtr.reset
394  (
395  new volScalarField
396  (
397  IOobject
398  (
399  alphaName_,
400  mesh_.time().timeName(),
401  mesh_,
404  ),
405  mesh_
406  )
407  );
408  }
409 
410 
411  const volScalarField& alpha =
412  (
413  alphaPtr.valid()
414  ? alphaPtr()
415  : obr_.lookupObject<volScalarField>(alphaName_)
416  );
417 
418  Info<< " Volume of alpha = "
419  << fvc::domainIntegrate(alpha).value()
420  << endl;
421 
422  const scalar meshVol = gSum(mesh_.V());
423  const scalar maxDropletVol = 1.0/6.0*pow(maxDiam_, 3);
424  const scalar delta = (maxDiam_-minDiam_)/nBins_;
425 
426  Info<< " Mesh volume = " << meshVol << endl;
427  Info<< " Maximum droplet diameter = " << maxDiam_ << endl;
428  Info<< " Maximum droplet volume = " << maxDropletVol << endl;
429 
430 
431  // Determine blocked faces
432  boolList blockedFace(mesh_.nFaces(), false);
433  label nBlocked = 0;
434 
435  {
436  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
437  {
438  scalar ownVal = alpha[mesh_.faceOwner()[facei]];
439  scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
440 
441  if
442  (
443  (ownVal < threshold_ && neiVal > threshold_)
444  || (ownVal > threshold_ && neiVal < threshold_)
445  )
446  {
447  blockedFace[facei] = true;
448  nBlocked++;
449  }
450  }
451 
452  // Block coupled faces
453  forAll(alpha.boundaryField(), patchi)
454  {
455  const fvPatchScalarField& fvp = alpha.boundaryField()[patchi];
456  if (fvp.coupled())
457  {
458  tmp<scalarField> townFld(fvp.patchInternalField());
459  const scalarField& ownFld = townFld();
460  tmp<scalarField> tnbrFld(fvp.patchNeighbourField());
461  const scalarField& nbrFld = tnbrFld();
462 
463  label start = fvp.patch().patch().start();
464 
465  forAll(ownFld, i)
466  {
467  scalar ownVal = ownFld[i];
468  scalar neiVal = nbrFld[i];
469 
470  if
471  (
472  (ownVal < threshold_ && neiVal > threshold_)
473  || (ownVal > threshold_ && neiVal < threshold_)
474  )
475  {
476  blockedFace[start+i] = true;
477  nBlocked++;
478  }
479  }
480  }
481  }
482  }
483 
484 
485  regionSplit regions(mesh_, blockedFace);
486 
487  Info<< " Determined " << regions.nRegions()
488  << " disconnected regions" << endl;
489 
490 
491  if (debug)
492  {
493  volScalarField region
494  (
495  IOobject
496  (
497  "region",
498  mesh_.time().timeName(),
499  mesh_,
502  ),
503  mesh_,
505  );
506  Info<< " Dumping region as volScalarField to " << region.name()
507  << endl;
508 
509  forAll(regions, celli)
510  {
511  region[celli] = regions[celli];
512  }
513  region.correctBoundaryConditions();
514  region.write();
515  }
516 
517 
518  // Determine regions connected to supplied patches
519  Map<label> patchRegions(findPatchRegions(regions));
520 
521 
522  // Sum all regions
523  const scalarField alphaVol(alpha.primitiveField()*mesh_.V());
524  Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
525  Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
526  Map<label> allRegionNumCells
527  (
528  regionSum
529  (
530  regions,
531  labelField(mesh_.nCells(), 1.0)
532  )
533  );
534 
535  if (debug)
536  {
537  Info<< " " << tab << "Region"
538  << tab << "Volume(mesh)"
539  << tab << "Volume(" << alpha.name() << "):"
540  << tab << "nCells"
541  << endl;
542  scalar meshSumVol = 0.0;
543  scalar alphaSumVol = 0.0;
544  label nCells = 0;
545 
546  Map<scalar>::const_iterator vIter = allRegionVolume.begin();
547  Map<scalar>::const_iterator aIter = allRegionAlphaVolume.begin();
548  Map<label>::const_iterator numIter = allRegionNumCells.begin();
549  for
550  (
551  ;
552  vIter != allRegionVolume.end()
553  && aIter != allRegionAlphaVolume.end();
554  ++vIter, ++aIter, ++numIter
555  )
556  {
557  Info<< " " << tab << vIter.key()
558  << tab << vIter()
559  << tab << aIter()
560  << tab << numIter()
561  << endl;
562 
563  meshSumVol += vIter();
564  alphaSumVol += aIter();
565  nCells += numIter();
566  }
567  Info<< " " << tab << "Total:"
568  << tab << meshSumVol
569  << tab << alphaSumVol
570  << tab << nCells
571  << endl;
572  Info<< endl;
573  }
574 
575 
576 
577 
578  {
579  Info<< " Patch connected regions (liquid core):" << endl;
580  Info<< tab << " Region"
581  << tab << "Volume(mesh)"
582  << tab << "Volume(" << alpha.name() << "):"
583  << endl;
584  forAllConstIter(Map<label>, patchRegions, iter)
585  {
586  label regionI = iter.key();
587  Info<< " " << tab << iter.key()
588  << tab << allRegionVolume[regionI]
589  << tab << allRegionAlphaVolume[regionI] << endl;
590 
591  }
592  Info<< endl;
593  }
594 
595  {
596  Info<< " Background regions:" << endl;
597  Info<< " " << tab << "Region"
598  << tab << "Volume(mesh)"
599  << tab << "Volume(" << alpha.name() << "):"
600  << endl;
601  Map<scalar>::const_iterator vIter = allRegionVolume.begin();
602  Map<scalar>::const_iterator aIter = allRegionAlphaVolume.begin();
603 
604  for
605  (
606  ;
607  vIter != allRegionVolume.end()
608  && aIter != allRegionAlphaVolume.end();
609  ++vIter, ++aIter
610  )
611  {
612  if
613  (
614  !patchRegions.found(vIter.key())
615  && vIter() >= maxDropletVol
616  )
617  {
618  Info<< " " << tab << vIter.key()
619  << tab << vIter()
620  << tab << aIter() << endl;
621  }
622  }
623  Info<< endl;
624  }
625 
626 
627 
628  // Split alpha field
629  // ~~~~~~~~~~~~~~~~~
630  // Split into
631  // - liquidCore : region connected to inlet patches
632  // - per region a volume : for all other regions
633  // - backgroundAlpha : remaining alpha
634  writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
635 
636 
637  // Extract droplet-only allRegionVolume, i.e. delete liquid core
638  // (patchRegions) and background regions from maps.
639  // Note that we have to use mesh volume (allRegionVolume) and not
640  // allRegionAlphaVolume since background might not have alpha in it.
641  forAllIter(Map<scalar>, allRegionVolume, vIter)
642  {
643  label regionI = vIter.key();
644  if
645  (
646  patchRegions.found(regionI)
647  || vIter() >= maxDropletVol
648  )
649  {
650  allRegionVolume.erase(vIter);
651  allRegionAlphaVolume.erase(regionI);
652  allRegionNumCells.erase(regionI);
653  }
654  }
655 
656  if (allRegionVolume.size())
657  {
658  // Construct mids of bins for plotting
659  pointField xBin(nBins_);
660 
661  scalar x = 0.5*delta;
662  forAll(xBin, i)
663  {
664  xBin[i] = point(x, 0, 0);
665  x += delta;
666  }
667 
668  const coordSet coords("diameter", "x", xBin, mag(xBin));
669 
670 
671  // Get in region order the alpha*volume and diameter
672  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
673 
674  const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
675 
676  scalarField sortedVols
677  (
678  extractData
679  (
680  sortedRegions,
681  allRegionAlphaVolume
682  )
683  );
684 
685  // Calculate the diameters
686  scalarField sortedDiameters(sortedVols.size());
687  forAll(sortedDiameters, i)
688  {
689  sortedDiameters[i] = Foam::cbrt
690  (
691  sortedVols[i]
693  );
694  }
695 
696  // Determine the bin index for all the diameters
697  labelList indices(sortedDiameters.size());
698  forAll(sortedDiameters, i)
699  {
700  indices[i] = (sortedDiameters[i]-minDiam_)/delta;
701  }
702 
703  // Calculate the counts per diameter bin
704  scalarField binCount(nBins_, 0.0);
705  forAll(sortedDiameters, i)
706  {
707  binCount[indices[i]] += 1.0;
708  }
709 
710  // Write counts
711  if (Pstream::master())
712  {
713  writeGraph(coords, "count", binCount);
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 
736  // Write average and deviation of droplet volume.
737  writeGraphs
738  (
739  "volume", // name of field
740  indices, // per region the bin index
741  sortedVols, // per region field data
742  binCount, // per bin number of regions
743  coords // graph data for bins
744  );
745 
746  // Collect some more field
747  {
748  wordList scalarNames(obr_.names(volScalarField::typeName));
749  labelList selected = findStrings(fields_, scalarNames);
750 
751  forAll(selected, i)
752  {
753  const word& fldName = scalarNames[selected[i]];
754  Info<< " Scalar field " << fldName << endl;
755 
756  const scalarField& fld = obr_.lookupObject
757  <
759  >(fldName).primitiveField();
760 
761  writeGraphs
762  (
763  fldName, // name of field
764  alphaVol*fld, // per cell field data
765 
766  regions, // per cell the region(=droplet)
767  sortedRegions, // valid regions in sorted order
768  1.0/sortedVols, // per region normalisation
769 
770  indices, // index of bin for each region
771  binCount, // per bin number of regions
772  coords // graph data for bins
773  );
774  }
775  }
776  {
777  wordList vectorNames(obr_.names(volVectorField::typeName));
778  labelList selected = findStrings(fields_, vectorNames);
779 
780  forAll(selected, i)
781  {
782  const word& fldName = vectorNames[selected[i]];
783  Info<< " Vector field " << fldName << endl;
784 
785  vectorField fld = obr_.lookupObject
786  <
788  >(fldName).primitiveField();
789 
790  if (coordSysPtr_.valid())
791  {
792  Info<< "Transforming vector field " << fldName
793  << " with coordinate system "
794  << coordSysPtr_().name()
795  << endl;
796 
797  fld = coordSysPtr_().localVector(fld);
798  }
799 
800 
801  // Components
802 
803  for (direction cmp = 0; cmp < vector::nComponents; cmp++)
804  {
805  writeGraphs
806  (
807  fldName + vector::componentNames[cmp],
808  alphaVol*fld.component(cmp),// per cell field data
809 
810  regions, // per cell the region(=droplet)
811  sortedRegions, // valid regions in sorted order
812  1.0/sortedVols, // per region normalisation
813 
814  indices, // index of bin for each region
815  binCount, // per bin number of regions
816  coords // graph data for bins
817  );
818  }
819 
820  // Magnitude
821  writeGraphs
822  (
823  fldName + "mag", // name of field
824  alphaVol*mag(fld), // per cell field data
825 
826  regions, // per cell the region(=droplet)
827  sortedRegions, // valid regions in sorted order
828  1.0/sortedVols, // per region normalisation
829 
830  indices, // index of bin for each region
831  binCount, // per bin number of regions
832  coords // graph data for bins
833  );
834  }
835  }
836  }
837 
838  return true;
839 }
840 
841 
842 // ************************************************************************* //
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Base class for other coordinate system specifications.
scalar delta
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:119
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual bool write()
Calculate the regionSizeDistribution and write.
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
const word & name() const
Return name.
Definition: IOobject.H:295
A class for handling file names.
Definition: fileName.H:79
An STL-conforming const_iterator.
Definition: HashTable.H:481
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:54
static const char tab
Definition: Ostream.H:264
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
static const char *const typeName
Definition: Field.H:105
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
Output to file stream.
Definition: OFstream.H:82
uint8_t direction
Definition: direction.H:45
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Plus op for FixedList<scalar>
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
word format(conversionProperties.lookup("format"))
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
Macros for easy insertion into run-time selection tables.
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
Definition: HashTable.C:371
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
scalar y
bool read(const char *, int32_t &)
Definition: int32IO.C:85
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
const Key & key() const
Return the Key corresponding to the iterator.
Definition: HashTableI.H:254
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Holds list of sampling positions.
Definition: coordSet.H:49
Type gSum(const FieldField< Field, Type > &f)
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:412
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:137
dimensionedScalar cbrt(const dimensionedScalar &ds)
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:313
Volume integrate volField creating a volField.
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:325
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:460
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:35
static const char nl
Definition: Ostream.H:265
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:194
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
vector point
Point is a vector.
Definition: point.H:41
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< label > sortedToc() const
Return the table of contents as a sorted list.
Definition: HashTable.C:217
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:201
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
virtual Ostream & write(const token &)=0
Write next token to stream.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: PtrList.H:53
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583