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