34 namespace functionObjects
41 regionSizeDistribution,
47 template<
class T,
unsigned Size>
68 void Foam::functionObjects::regionSizeDistribution::writeGraph
71 const word& valueName,
75 const wordList valNames(1, valueName);
77 fileName outputPath = file_.baseTimeDir();
80 OFstream str(outputPath/formatterPtr_().getFileName(coords, valNames));
82 Info<<
" Writing distribution of " << valueName <<
" to " << str.
name()
87 formatterPtr_().write(coords, valNames, valPtrs, str);
91 void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
99 const scalar maxDropletVol = 1.0/6.0*
pow(maxDiam_, 3);
114 alphaName_ +
"_liquidCore",
115 obr_.time().timeName(),
127 alphaName_ +
"_background",
128 obr_.time().timeName(),
140 label regionI = regions[celli];
141 if (patchRegions.
found(regionI))
143 backgroundAlpha[celli] = 0;
147 liquidCore[celli] = 0;
149 scalar regionVol = regionVolume[regionI];
150 if (regionVol < maxDropletVol)
152 backgroundAlpha[celli] = 0;
156 liquidCore.correctBoundaryConditions();
157 backgroundAlpha.correctBoundaryConditions();
159 Info<<
" Volume of liquid-core = " 162 Info<<
" Volume of background = " 166 Info<<
" Writing liquid-core field to " << liquidCore.name() <<
endl;
168 Info<<
" Writing background field to " << backgroundAlpha.name() <<
endl;
169 backgroundAlpha.
write();
174 Foam::functionObjects::regionSizeDistribution::findPatchRegions
183 const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
185 label nPatchFaces = 0;
188 nPatchFaces += mesh_.boundaryMesh()[iter.key()].size();
204 regions[faceCells[i]],
220 Foam::functionObjects::regionSizeDistribution::divide
233 result[i] = num[i]/denom[i];
244 void Foam::functionObjects::regionSizeDistribution::writeGraphs
246 const word& fieldName,
259 binSum[indices[i]] += sortedField[i];
268 binSqrSum[indices[i]] +=
Foam::sqr(sortedField[i]);
276 writeGraph(coords, fieldName +
"_sum", binSum);
278 writeGraph(coords, fieldName +
"_avg", binAvg);
280 writeGraph(coords, fieldName +
"_dev", binDev);
285 void Foam::functionObjects::regionSizeDistribution::writeGraphs
287 const word& fieldName,
299 Map<scalar> regionField(regionSum(regions, cellField));
334 alphaName_(dict.
lookup(
"field")),
335 patchNames_(dict.
lookup(
"patches"))
351 dict.
lookup(
"field") >> alphaName_;
352 dict.
lookup(
"patches") >> patchNames_;
353 dict.
lookup(
"threshold") >> threshold_;
354 dict.
lookup(
"maxDiameter") >> maxDiam_;
357 dict.
lookup(
"nBins") >> nBins_;
358 dict.
lookup(
"fields") >> fields_;
363 if (dict.
found(
"coordinateSystem"))
367 Info<<
"Transforming all vectorFields with coordinate system " 368 << coordSysPtr_().name() <<
endl;
388 Info<<
" Looking up field " << alphaName_ <<
endl;
392 Info<<
" Reading field " << alphaName_ <<
endl;
400 mesh_.time().timeName(),
418 Info<<
" Volume of alpha = " 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_;
426 Info<<
" Mesh volume = " << meshVol <<
endl;
427 Info<<
" Maximum droplet diameter = " << maxDiam_ <<
endl;
428 Info<<
" Maximum droplet volume = " << maxDropletVol <<
endl;
432 boolList blockedFace(mesh_.nFaces(),
false);
436 for (
label facei = 0; facei < mesh_.nInternalFaces(); facei++)
438 scalar ownVal = alpha[mesh_.faceOwner()[facei]];
439 scalar neiVal = alpha[mesh_.faceNeighbour()[facei]];
443 (ownVal < threshold_ && neiVal > threshold_)
444 || (ownVal > threshold_ && neiVal < threshold_)
447 blockedFace[facei] =
true;
467 scalar ownVal = ownFld[i];
468 scalar neiVal = nbrFld[i];
472 (ownVal < threshold_ && neiVal > threshold_)
473 || (ownVal > threshold_ && neiVal < threshold_)
476 blockedFace[start+i] =
true;
488 <<
" disconnected regions" <<
endl;
498 mesh_.time().timeName(),
506 Info<<
" Dumping region as volScalarField to " << region.name()
511 region[celli] = regions[celli];
513 region.correctBoundaryConditions();
519 Map<label> patchRegions(findPatchRegions(regions));
524 Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
525 Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
538 <<
tab <<
"Volume(mesh)" 539 <<
tab <<
"Volume(" << alpha.
name() <<
"):" 542 scalar meshSumVol = 0.0;
543 scalar alphaSumVol = 0.0;
552 vIter != allRegionVolume.end()
553 && aIter != allRegionAlphaVolume.
end();
554 ++vIter, ++aIter, ++numIter
563 meshSumVol += vIter();
564 alphaSumVol += aIter();
569 <<
tab << alphaSumVol
579 Info<<
" Patch connected regions (liquid core):" <<
endl;
581 <<
tab <<
"Volume(mesh)" 582 <<
tab <<
"Volume(" << alpha.
name() <<
"):" 586 label regionI = iter.key();
587 Info<<
" " <<
tab << iter.key()
588 <<
tab << allRegionVolume[regionI]
589 <<
tab << allRegionAlphaVolume[regionI] <<
endl;
596 Info<<
" Background regions:" <<
endl;
598 <<
tab <<
"Volume(mesh)" 599 <<
tab <<
"Volume(" << alpha.
name() <<
"):" 607 vIter != allRegionVolume.end()
608 && aIter != allRegionAlphaVolume.
end();
615 && vIter() >= maxDropletVol
634 writeAlphaFields(regions, patchRegions, allRegionVolume, alpha);
643 label regionI = vIter.key();
646 patchRegions.
found(regionI)
647 || vIter() >= maxDropletVol
650 allRegionVolume.erase(vIter);
651 allRegionAlphaVolume.
erase(regionI);
652 allRegionNumCells.
erase(regionI);
656 if (allRegionVolume.size())
664 xBin[i] =
point(x, 0, 0);
668 const coordSet coords(
"diameter",
"x", xBin,
mag(xBin));
687 forAll(sortedDiameters, i)
697 labelList indices(sortedDiameters.size());
698 forAll(sortedDiameters, i)
700 indices[i] = (sortedDiameters[i]-minDiam_)/delta;
705 forAll(sortedDiameters, i)
707 binCount[indices[i]] += 1.0;
713 writeGraph(coords,
"count", binCount);
720 <<
tab <<
"Min diameter" 729 <<
tab << binCount[binI] <<
endl;
753 const word& fldName = scalarNames[selected[i]];
754 Info<<
" Scalar field " << fldName <<
endl;
759 >(fldName).primitiveField();
782 const word& fldName = vectorNames[selected[i]];
783 Info<<
" Vector field " << fldName <<
endl;
788 >(fldName).primitiveField();
790 if (coordSysPtr_.valid())
792 Info<<
"Transforming vector field " << fldName
793 <<
" with coordinate system " 794 << coordSysPtr_().name()
797 fld = coordSysPtr_().localVector(fld);
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Field< label > labelField
Specialisation of Field<T> for label.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
#define forAll(list, i)
Loop across all elements in list.
virtual bool write()
Calculate the regionSizeDistribution and write.
virtual Ostream & write(const char)=0
Write character.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
A class for handling file names.
An STL-conforming const_iterator.
virtual bool execute()
Do nothing.
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
A 1D vector of objects of type <T> with a fixed size <Size>.
A list of keyword definitions, which are a keyword followed by any number of values (e...
static const char *const typeName
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
void size(const label)
Override size to be inconsistent with allocated storage.
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
static const char *const componentNames[]
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Plus op for FixedList<scalar>
static bool master(const label communicator=0)
Am I the master process.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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 dimensionSet dimless
const fileName & name() const
Return the name of the stream.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
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.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
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.
bool read(const char *, int32_t &)
static const direction nComponents
Number of components in this vector space.
const Key & key() const
Return the Key corresponding to the iterator.
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.
Type gSum(const FieldField< Field, Type > &f)
const labelUList & faceCells() const
Return face-cell addressing.
bool found(const Key &) const
Return true if hashedEntry is found in table.
A class for handling words, derived from string.
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.
const polyPatch & patch() const
Return the polyPatch.
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.
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.
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
const fvPatch & patch() const
Return patch.
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
virtual ~regionSizeDistribution()
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
word name(const complex &)
Return a string representation of a complex.
defineTypeNameAndDebug(Qdot, 0)
static autoPtr< setWriter > New(const word &writeFormat)
Return a reference to the selected setWriter.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
static autoPtr< coordinateSystem > New(const objectRegistry &obr, const dictionary &dict)
Select constructed from dictionary and objectRegistry.
vector point
Point is a vector.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< label > sortedToc() const
Return the table of contents as a sorted list.
label start() const
Return start label of this patch in the polyMesh face list.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
label nRegions() const
Return total number of regions.
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...
Specialisation 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.
A patch is a list of labels that address the faces in the global face list.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.