regionSizeDistribution.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 | Copyright (C) 2013-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "volFields.H"
28 #include "regionSplit.H"
29 #include "fvcVolumeIntegrate.H"
30 #include "mathematicalConstants.H"
31 #include "stringListOps.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37  defineTypeNameAndDebug(regionSizeDistribution, 0);
38 
39  //- Plus op for FixedList<scalar>
40  template<class T, unsigned Size>
42  {
43  public:
44  void operator()
45  (
47  const FixedList<T, Size>& y
48  ) const
49  {
50  forAll(x, i)
51  {
52  x[i] += y[i];
53  }
54  }
55  };
56 }
57 
58 
59 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
60 
62 (
63  const coordSet& coords,
64  const word& valueName,
65  const scalarField& values
66 ) const
67 {
68  const wordList valNames(1, valueName);
69 
70  fileName outputPath = baseTimeDir();
71  mkDir(outputPath);
72 
73  OFstream str(outputPath/formatterPtr_().getFileName(coords, valNames));
74 
75  if (log_) Info
76  << "Writing distribution of " << valueName << " to " << str.name()
77  << endl;
78 
79  List<const scalarField*> valPtrs(1);
80  valPtrs[0] = &values;
81  formatterPtr_().write(coords, valNames, valPtrs, str);
82 }
83 
84 
86 (
87  const regionSplit& regions,
88  const Map<label>& patchRegions,
89  const Map<scalar>& regionVolume,
90  const volScalarField& alpha
91 ) const
92 {
93  const scalar maxDropletVol = 1.0/6.0*pow(maxDiam_, 3);
94 
95  // Split alpha field
96  // ~~~~~~~~~~~~~~~~~
97  // Split into
98  // - liquidCore : region connected to inlet patches
99  // - per region a volume : for all other regions
100  // - backgroundAlpha : remaining alpha
101 
102 
103  // Construct field
104  volScalarField liquidCore
105  (
106  IOobject
107  (
108  alphaName_ + "_liquidCore",
109  obr_.time().timeName(),
110  obr_,
111  IOobject::NO_READ
112  ),
113  alpha,
115  );
116 
117  volScalarField backgroundAlpha
118  (
119  IOobject
120  (
121  alphaName_ + "_background",
122  obr_.time().timeName(),
123  obr_,
124  IOobject::NO_READ
125  ),
126  alpha,
128  );
129 
130 
131  // Knock out any cell not in patchRegions
132  forAll(liquidCore, cellI)
133  {
134  label regionI = regions[cellI];
135  if (patchRegions.found(regionI))
136  {
137  backgroundAlpha[cellI] = 0;
138  }
139  else
140  {
141  liquidCore[cellI] = 0;
142 
143  scalar regionVol = regionVolume[regionI];
144  if (regionVol < maxDropletVol)
145  {
146  backgroundAlpha[cellI] = 0;
147  }
148  }
149  }
150  liquidCore.correctBoundaryConditions();
151  backgroundAlpha.correctBoundaryConditions();
152 
153  if (log_)
154  {
155  Info<< " Volume of liquid-core = "
156  << fvc::domainIntegrate(liquidCore).value()
157  << endl;
158  Info<< " Volume of background = "
159  << fvc::domainIntegrate(backgroundAlpha).value()
160  << endl;
161  }
162 
163  if (log_) Info
164  << "Writing liquid-core field to " << liquidCore.name() << endl;
165  liquidCore.write();
166 
167  if (log_) Info
168  << "Writing background field to " << backgroundAlpha.name() << endl;
169  backgroundAlpha.write();
170 }
171 
172 
174 (
175  const polyMesh& mesh,
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 (
221  const scalarField& num,
222  const scalarField& denom
223 )
224 {
225  tmp<scalarField> tresult(new scalarField(num.size()));
226  scalarField& result = tresult();
227 
228  forAll(denom, i)
229  {
230  if (denom[i] != 0)
231  {
232  result[i] = num[i]/denom[i];
233  }
234  else
235  {
236  result[i] = 0.0;
237  }
238  }
239  return tresult;
240 }
241 
242 
244 (
245  const word& fieldName, // name of field
246  const labelList& indices, // index of bin for each region
247  const scalarField& sortedField, // per region field data
248  const scalarField& binCount, // per bin number of regions
249  const coordSet& coords // graph data for bins
250 ) const
251 {
252  if (Pstream::master())
253  {
254  // Calculate per-bin average
255  scalarField binSum(nBins_, 0.0);
256  forAll(sortedField, i)
257  {
258  binSum[indices[i]] += sortedField[i];
259  }
260 
261  scalarField binAvg(divide(binSum, binCount));
262 
263  // Per bin deviation
264  scalarField binSqrSum(nBins_, 0.0);
265  forAll(sortedField, i)
266  {
267  binSqrSum[indices[i]] += Foam::sqr(sortedField[i]);
268  }
269  scalarField binDev
270  (
271  sqrt(divide(binSqrSum, binCount) - Foam::sqr(binAvg))
272  );
273 
274  // Write average
275  writeGraph(coords, fieldName + "_sum", binSum);
276  // Write average
277  writeGraph(coords, fieldName + "_avg", binAvg);
278  // Write deviation
279  writeGraph(coords, fieldName + "_dev", binDev);
280  }
281 }
282 
283 
285 (
286  const word& fieldName, // name of field
287  const scalarField& cellField, // per cell field data
288  const regionSplit& regions, // per cell the region(=droplet)
289  const labelList& sortedRegions, // valid regions in sorted order
290  const scalarField& sortedNormalisation,
291 
292  const labelList& indices, // per region index of bin
293  const scalarField& binCount, // per bin number of regions
294  const coordSet& coords // graph data for bins
295 ) const
296 {
297  // Sum on a per-region basis. Parallel reduced.
298  Map<scalar> regionField(regionSum(regions, cellField));
299 
300  // Extract in region order
301  scalarField sortedField
302  (
303  sortedNormalisation
304  * extractData
305  (
306  sortedRegions,
307  regionField
308  )
309  );
310 
311  writeGraphs
312  (
313  fieldName, // name of field
314  indices, // index of bin for each region
315  sortedField, // per region field data
316  binCount, // per bin number of regions
317  coords // graph data for bins
318  );
319 }
320 
321 
322 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
323 
325 (
326  const word& name,
327  const objectRegistry& obr,
328  const dictionary& dict,
329  const bool loadFromFiles
330 )
331 :
332  functionObjectFile(obr, name),
333  name_(name),
334  obr_(obr),
335  active_(true),
336  alphaName_(dict.lookup("field")),
337  patchNames_(dict.lookup("patches")),
338  log_(true)
339 {
340  // Check if the available mesh is an fvMesh, otherwise deactivate
341  if (isA<fvMesh>(obr_))
342  {
343  read(dict);
344  }
345  else
346  {
347  active_ = false;
349  << "No fvMesh available, deactivating " << name_ << nl
350  << endl;
351  }
352 }
353 
354 
355 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
356 
358 {}
359 
360 
361 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
362 
364 {
365  if (active_)
366  {
368 
369  log_.readIfPresent("log", dict);
370 
371  dict.lookup("field") >> alphaName_;
372  dict.lookup("patches") >> patchNames_;
373  dict.lookup("threshold") >> threshold_;
374  dict.lookup("maxDiameter") >> maxDiam_;
375  minDiam_ = 0.0;
376  dict.readIfPresent("minDiameter", minDiam_);
377  dict.lookup("nBins") >> nBins_;
378  dict.lookup("fields") >> fields_;
379 
380  word format(dict.lookup("setFormat"));
381  formatterPtr_ = writer<scalar>::New(format);
382 
383  if (dict.found("coordinateSystem"))
384  {
385  coordSysPtr_.reset(new coordinateSystem(obr_, dict));
386 
387  if (log_) Info
388  << "Transforming all vectorFields with coordinate system "
389  << coordSysPtr_().name() << endl;
390  }
391  }
392 }
393 
394 
396 {
397  // Do nothing - only valid on write
398 }
399 
400 
402 {
403  // Do nothing - only valid on write
404 }
405 
406 
408 {
409  // Do nothing - only valid on write
410 }
411 
412 
414 {
415  if (active_)
416  {
417  if (log_) Info << type() << " " << name_ << " output:" << nl;
418 
419  const fvMesh& mesh = refCast<const fvMesh>(obr_);
420 
422  if (obr_.foundObject<volScalarField>(alphaName_))
423  {
424  if (log_) Info << " Looking up field " << alphaName_ << endl;
425  }
426  else
427  {
428  if (log_) Info << " Reading field " << alphaName_ << endl;
429  alphaPtr.reset
430  (
431  new volScalarField
432  (
433  IOobject
434  (
435  alphaName_,
436  mesh.time().timeName(),
437  mesh,
440  ),
441  mesh
442  )
443  );
444  }
445 
446 
447  const volScalarField& alpha =
448  (
449  alphaPtr.valid()
450  ? alphaPtr()
451  : obr_.lookupObject<volScalarField>(alphaName_)
452  );
453 
454  if (log_) Info
455  << " Volume of alpha = "
457  << endl;
458 
459  const scalar meshVol = gSum(mesh.V());
460  const scalar maxDropletVol = 1.0/6.0*pow(maxDiam_, 3);
461  const scalar delta = (maxDiam_-minDiam_)/nBins_;
462 
463  if (log_) Info
464  << " Mesh volume = " << meshVol << nl
465  << " Maximum droplet diameter = " << maxDiam_ << nl
466  << " Maximum droplet volume = " << maxDropletVol
467  << endl;
468 
469 
470  // Determine blocked faces
471  boolList blockedFace(mesh.nFaces(), false);
472  label nBlocked = 0;
473 
474  {
475  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
476  {
477  scalar ownVal = alpha[mesh.faceOwner()[faceI]];
478  scalar neiVal = alpha[mesh.faceNeighbour()[faceI]];
479 
480  if
481  (
482  (ownVal < threshold_ && neiVal > threshold_)
483  || (ownVal > threshold_ && neiVal < threshold_)
484  )
485  {
486  blockedFace[faceI] = true;
487  nBlocked++;
488  }
489  }
490 
491  // Block coupled faces
492  forAll(alpha.boundaryField(), patchI)
493  {
494  const fvPatchScalarField& fvp = alpha.boundaryField()[patchI];
495  if (fvp.coupled())
496  {
497  tmp<scalarField> townFld(fvp.patchInternalField());
498  const scalarField& ownFld = townFld();
499  tmp<scalarField> tnbrFld(fvp.patchNeighbourField());
500  const scalarField& nbrFld = tnbrFld();
501 
502  label start = fvp.patch().patch().start();
503 
504  forAll(ownFld, i)
505  {
506  scalar ownVal = ownFld[i];
507  scalar neiVal = nbrFld[i];
508 
509  if
510  (
511  (ownVal < threshold_ && neiVal > threshold_)
512  || (ownVal > threshold_ && neiVal < threshold_)
513  )
514  {
515  blockedFace[start+i] = true;
516  nBlocked++;
517  }
518  }
519  }
520  }
521  }
522 
523 
524  regionSplit regions(mesh, blockedFace);
525 
526  if (log_) Info
527  << " Determined " << regions.nRegions()
528  << " disconnected regions" << endl;
529 
530 
531  if (debug)
532  {
533  volScalarField region
534  (
535  IOobject
536  (
537  "region",
538  mesh.time().timeName(),
539  mesh,
542  ),
543  mesh,
544  dimensionedScalar("zero", dimless, 0)
545  );
546 
547  if (log_) Info
548  << " Dumping region as " << volScalarField::typeName
549  << " to " << region.name() << endl;
550 
551  forAll(regions, cellI)
552  {
553  region[cellI] = regions[cellI];
554  }
555  region.correctBoundaryConditions();
556  region.write();
557  }
558 
559 
560  // Determine regions connected to supplied patches
561  Map<label> patchRegions(findPatchRegions(mesh, regions));
562 
563 
564 
565  // Sum all regions
566  const scalarField alphaVol(alpha.internalField()*mesh.V());
567  Map<scalar> allRegionVolume(regionSum(regions, mesh.V()));
568  Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
569  Map<label> allRegionNumCells
570  (
571  regionSum
572  (
573  regions,
574  labelField(mesh.nCells(), 1.0)
575  )
576  );
577 
578  if (debug)
579  {
580  if (log_) Info
581  << " " << token::TAB << "Region"
582  << token::TAB << "Volume(mesh)"
583  << token::TAB << "Volume(" << alpha.name() << "):"
584  << token::TAB << "nCells"
585  << endl;
586 
587  scalar meshSumVol = 0.0;
588  scalar alphaSumVol = 0.0;
589  label nCells = 0;
590 
591  Map<scalar>::const_iterator vIter = allRegionVolume.begin();
592  Map<scalar>::const_iterator aIter = allRegionAlphaVolume.begin();
593  Map<label>::const_iterator numIter = allRegionNumCells.begin();
594  for
595  (
596  ;
597  vIter != allRegionVolume.end()
598  && aIter != allRegionAlphaVolume.end();
599  ++vIter, ++aIter, ++numIter
600  )
601  {
602  if (log_) Info
603  << " " << token::TAB << vIter.key()
604  << token::TAB << vIter()
605  << token::TAB << aIter()
606  << token::TAB << numIter()
607  << endl;
608 
609  meshSumVol += vIter();
610  alphaSumVol += aIter();
611  nCells += numIter();
612  }
613 
614  if (log_) Info
615  << " " << token::TAB << "Total:"
616  << token::TAB << meshSumVol
617  << token::TAB << alphaSumVol
618  << token::TAB << nCells
619  << nl << endl;
620  }
621 
622 
623 
624  if (log_)
625  {
626  Info<< " Patch connected regions (liquid core):" << nl
627  << token::TAB << " Region"
628  << token::TAB << "Volume(mesh)"
629  << token::TAB << "Volume(" << alpha.name() << "):"
630  << endl;
631 
632  forAllConstIter(Map<label>, patchRegions, iter)
633  {
634  label regionI = iter.key();
635 
636  Info<< " " << token::TAB << iter.key()
637  << token::TAB << allRegionVolume[regionI]
638  << token::TAB << allRegionAlphaVolume[regionI] << endl;
639 
640  }
641 
642  Info<< endl;
643  }
644 
645  if (log_)
646  {
647  Info<< " Background regions:" << nl
648  << " " << token::TAB << "Region"
649  << token::TAB << "Volume(mesh)"
650  << token::TAB << "Volume(" << alpha.name() << "):"
651  << endl;
652 
653  Map<scalar>::const_iterator vIter = allRegionVolume.begin();
654  Map<scalar>::const_iterator aIter = allRegionAlphaVolume.begin();
655 
656  for
657  (
658  ;
659  vIter != allRegionVolume.end()
660  && aIter != allRegionAlphaVolume.end();
661  ++vIter, ++aIter
662  )
663  {
664  if
665  (
666  !patchRegions.found(vIter.key())
667  && vIter() >= maxDropletVol
668  )
669  {
670  Info<< " " << token::TAB << vIter.key()
671  << token::TAB << vIter()
672  << token::TAB << aIter() << endl;
673  }
674  }
675 
676  Info<< endl;
677  }
678 
679 
680 
681  // Split alpha field
682  // ~~~~~~~~~~~~~~~~~
683  // Split into
684  // - liquidCore : region connected to inlet patches
685  // - per region a volume : for all other regions
686  // - backgroundAlpha : remaining alpha
687  writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
688 
689 
690  // Extract droplet-only allRegionVolume, i.e. delete liquid core
691  // (patchRegions) and background regions from maps.
692  // Note that we have to use mesh volume (allRegionVolume) and not
693  // allRegionAlphaVolume since background might not have alpha in it.
694  forAllIter(Map<scalar>, allRegionVolume, vIter)
695  {
696  label regionI = vIter.key();
697  if
698  (
699  patchRegions.found(regionI)
700  || vIter() >= maxDropletVol
701  )
702  {
703  allRegionVolume.erase(vIter);
704  allRegionAlphaVolume.erase(regionI);
705  allRegionNumCells.erase(regionI);
706  }
707  }
708 
709  if (allRegionVolume.size())
710  {
711  // Construct mids of bins for plotting
712  pointField xBin(nBins_);
713 
714  scalar x = 0.5*delta;
715  forAll(xBin, i)
716  {
717  xBin[i] = point(x, 0, 0);
718  x += delta;
719  }
720 
721  const coordSet coords("diameter", "x", xBin, mag(xBin));
722 
723 
724  // Get in region order the alpha*volume and diameter
725  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
726 
727  const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
728 
729  scalarField sortedVols
730  (
731  extractData
732  (
733  sortedRegions,
734  allRegionAlphaVolume
735  )
736  );
737 
738  // Calculate the diameters
739  scalarField sortedDiameters(sortedVols.size());
740  forAll(sortedDiameters, i)
741  {
742  sortedDiameters[i] = Foam::cbrt
743  (
744  sortedVols[i]
746  );
747  }
748 
749  // Determine the bin index for all the diameters
750  labelList indices(sortedDiameters.size());
751  forAll(sortedDiameters, i)
752  {
753  indices[i] = (sortedDiameters[i]-minDiam_)/delta;
754  }
755 
756  // Calculate the counts per diameter bin
757  scalarField binCount(nBins_, 0.0);
758  forAll(sortedDiameters, i)
759  {
760  binCount[indices[i]] += 1.0;
761  }
762 
763  // Write counts
764  if (Pstream::master())
765  {
766  writeGraph(coords, "count", binCount);
767  }
768 
769  // Write to screen
770  if (log_)
771  {
772  Info<< " Bins:" << nl
773  << " " << token::TAB << "Bin"
774  << token::TAB << "Min diameter"
775  << token::TAB << "Count:"
776  << endl;
777 
778  scalar diam = 0.0;
779  forAll(binCount, binI)
780  {
781  Info<< " " << token::TAB << binI
782  << token::TAB << diam
783  << token::TAB << binCount[binI] << endl;
784 
785  diam += delta;
786  }
787 
788  Info<< endl;
789  }
790 
791 
792  // Write average and deviation of droplet volume.
793  writeGraphs
794  (
795  "volume", // name of field
796  indices, // per region the bin index
797  sortedVols, // per region field data
798  binCount, // per bin number of regions
799  coords // graph data for bins
800  );
801 
802  // Collect some more field
803  {
804  wordList scalarNames(obr_.names(volScalarField::typeName));
805  labelList selected = findStrings(fields_, scalarNames);
806 
807  forAll(selected, i)
808  {
809  const word& fldName = scalarNames[selected[i]];
810  if (log_) Info << " Scalar field " << fldName << endl;
811 
812  const scalarField& fld = obr_.lookupObject
813  <
815  >(fldName).internalField();
816 
817  writeGraphs
818  (
819  fldName, // name of field
820  alphaVol*fld, // per cell field data
821 
822  regions, // per cell the region(=droplet)
823  sortedRegions, // valid regions in sorted order
824  1.0/sortedVols, // per region normalisation
825 
826  indices, // index of bin for each region
827  binCount, // per bin number of regions
828  coords // graph data for bins
829  );
830  }
831  }
832  {
833  wordList vectorNames(obr_.names(volVectorField::typeName));
834  labelList selected = findStrings(fields_, vectorNames);
835 
836  forAll(selected, i)
837  {
838  const word& fldName = vectorNames[selected[i]];
839  if (log_) Info << " Vector field " << fldName << endl;
840 
841  vectorField fld = obr_.lookupObject
842  <
844  >(fldName).internalField();
845 
846  if (coordSysPtr_.valid())
847  {
848  if (log_) Info
849  << "Transforming vector field " << fldName
850  << " with coordinate system "
851  << coordSysPtr_().name() << endl;
852 
853  fld = coordSysPtr_().localVector(fld);
854  }
855 
856 
857  // Components
858 
859  for (direction cmp = 0; cmp < vector::nComponents; cmp++)
860  {
861  writeGraphs
862  (
863  fldName + vector::componentNames[cmp],
864  alphaVol*fld.component(cmp),// per cell field data
865 
866  regions, // per cell the region(=droplet)
867  sortedRegions, // valid regions in sorted order
868  1.0/sortedVols, // per region normalisation
869 
870  indices, // index of bin for each region
871  binCount, // per bin number of regions
872  coords // graph data for bins
873  );
874  }
875 
876  // Magnitude
877  writeGraphs
878  (
879  fldName + "mag", // name of field
880  alphaVol*mag(fld), // per cell field data
881 
882  regions, // per cell the region(=droplet)
883  sortedRegions, // valid regions in sorted order
884  1.0/sortedVols, // per region normalisation
885 
886  indices, // index of bin for each region
887  binCount, // per bin number of regions
888  coords // graph data for bins
889  );
890  }
891  }
892  }
893  }
894 }
895 
896 
897 // ************************************************************************* //
Foam::token::TAB
@ TAB
Definition: token.H:96
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
regionSizeDistribution.H
volFields.H
Foam::regionSizeDistribution::findPatchRegions
Map< label > findPatchRegions(const polyMesh &, const regionSplit &) const
Mark all regions starting at patches.
Definition: regionSizeDistribution.C:174
Foam::regionSizeDistribution::regionSizeDistribution
regionSizeDistribution(const regionSizeDistribution &)
Disallow default bitwise copy construct.
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
mathematicalConstants.H
format
word format(conversionProperties.lookup("format"))
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::regionSizeDistribution::divide
static tmp< scalarField > divide(const scalarField &, const scalarField &)
Helper: divide if denom != 0.
Definition: regionSizeDistribution.C:220
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dictionary::readIfPresent
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Definition: dictionaryTemplates.C:94
Foam::ListPlusEqOp
Plus op for FixedList<scalar>
Definition: regionSizeDistribution.C:41
Foam::fvc::domainIntegrate
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcVolumeIntegrate.C:86
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
alphaPtr
autoPtr< volScalarField > alphaPtr
Definition: readThermalProperties.H:162
Foam::HashTable::const_iterator
An STL-conforming const_iterator.
Definition: HashTable.H:470
Foam::regionSizeDistribution::writeGraphs
void writeGraphs(const word &fieldName, const labelList &indices, const scalarField &sortedField, const scalarField &binCount, const coordSet &coords) const
Given per-region data calculate per-bin average/deviation and graph.
Definition: regionSizeDistribution.C:244
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::dimensioned::name
const word & name() const
Return const reference to name.
Definition: dimensionedType.C:235
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::nComponents
@ nComponents
Number of components in this vector space.
Definition: VectorSpace.H:88
Foam::Map< label >
Foam::functionObjectFile::read
void read(const dictionary &dict)
Read.
Definition: functionObjectFile.C:178
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::regionSizeDistribution::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: regionSizeDistribution.C:407
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::Vector< scalar >::componentNames
static const char * componentNames[]
Definition: Vector.H:79
Foam::HashSet< label, Hash< label > >
Foam::regionSizeDistribution::~regionSizeDistribution
virtual ~regionSizeDistribution()
Definition: regionSizeDistribution.C:357
Foam::regionSizeDistribution::writeAlphaFields
void writeAlphaFields(const regionSplit &regions, const Map< label > &keepRegions, const Map< scalar > &regionVolume, const volScalarField &alpha) const
Write volfields with the parts of alpha which are not.
Definition: regionSizeDistribution.C:86
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::divide
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::regionSizeDistribution::execute
virtual void execute()
Execute, currently does nothing.
Definition: regionSizeDistribution.C:395
Foam::fvPatchField::coupled
virtual bool coupled() const
Return true if this patch field is coupled.
Definition: fvPatchField.H:342
Foam::fvMesh::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:199
Foam::label
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
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::labelField
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::minEqOp
Definition: ops.H:78
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
regionSplit.H
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::fvPatchField::patchInternalField
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:208
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:119
Foam::regionSplit::nRegions
label nRegions() const
Return total number of regions.
Definition: regionSplit.H:201
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
fld
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))
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::coordSet
Holds list of sampling positions.
Definition: coordSet.H:49
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::polyPatch::faceCells
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:340
Foam::regionSizeDistribution::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: regionSizeDistribution.C:401
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::fvPatchField::patchNeighbourField
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
Definition: fvPatchField.H:408
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
internalField
conserve internalField()+
Foam::HashTable::iteratorBase::key
const Key & key() const
Return the Key corresponding to the iterator.
Definition: HashTableI.H:260
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::findStrings
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::functionObjectFile
Base class for output file data handling.
Definition: functionObjectFile.H:57
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::regionSizeDistribution::read
virtual void read(const dictionary &)
Read the regionSizeDistribution data.
Definition: regionSizeDistribution.C:363
fvcVolumeIntegrate.H
Volume integrate volField creating a volField.
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
stringListOps.H
Operations on lists of strings.
Foam::OFstream::name
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:153
Foam::regionSizeDistribution::write
virtual void write()
Calculate the regionSizeDistribution and write.
Definition: regionSizeDistribution.C:413
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::mkDir
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:419
Foam::regionSizeDistribution::writeGraph
void writeGraph(const coordSet &coords, const word &valueName, const scalarField &values) const
Definition: regionSizeDistribution.C:62
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
Foam::coordinateSystem
Base class for other coordinate system specifications.
Definition: coordinateSystem.H:85
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::writer::New
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.
Definition: writer.C:35