39 namespace functionObjects
65 const label regioni = regions[celli];
66 regionToSum(regioni, Type(
Zero)) +=
fld[celli];
82 sortedData[i] = regionData[keys[i]];
92void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
94 const regionSplit& regions,
95 const Map<label>& patchRegions,
96 const Map<scalar>& regionVolume,
97 const volScalarField&
alpha
100 const scalar maxDropletVol = 1.0/6.0*
pow3(maxDiam_);
121 fvPatchField<scalar>::calculatedType()
134 fvPatchField<scalar>::calculatedType()
141 const label regioni = regions[celli];
142 if (patchRegions.found(regioni))
144 backgroundAlpha[celli] = 0;
148 liquidCore[celli] = 0;
150 const scalar regionVol = regionVolume[regioni];
151 if (regionVol < maxDropletVol)
153 backgroundAlpha[celli] = 0;
157 liquidCore.correctBoundaryConditions();
158 backgroundAlpha.correctBoundaryConditions();
162 Info<<
" Volume of liquid-core = "
165 <<
" Volume of background = "
170 Log <<
" Writing liquid-core field to " << liquidCore.name() <<
endl;
173 Log<<
" Writing background field to " << backgroundAlpha.name() <<
endl;
174 backgroundAlpha.
write();
179Foam::functionObjects::regionSizeDistribution::findPatchRegions
181 const regionSplit& regions
188 const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
190 label nPatchFaces = 0;
191 for (
const label patchi : patchIDs)
193 nPatchFaces += mesh_.boundaryMesh()[patchi].size();
197 Map<label> patchRegions(nPatchFaces);
198 for (
const label patchi : patchIDs)
200 const polyPatch& pp = mesh_.boundaryMesh()[patchi];
203 const labelList& faceCells = pp.faceCells();
205 for (
const label celli : faceCells)
224Foam::functionObjects::regionSizeDistribution::divide
231 auto& result = tresult.ref();
237 result[i] = num[i]/denom[i];
248void Foam::functionObjects::regionSizeDistribution::writeGraphs
250 const word& fieldName,
255 const coordSet& coords
264 binSum[indices[i]] += sortedField[i];
273 binSqrSum[indices[i]] +=
Foam::sqr(sortedField[i]);
281 auto&
writer = formatterPtr_();
311 Log <<
" Writing distribution of "
312 << fieldName <<
" to " <<
writer.path() <<
endl;
314 writer.write(fieldName +
"_sum", binSum);
315 writer.write(fieldName +
"_avg", binAvg);
316 writer.write(fieldName +
"_dev", binDev);
322void Foam::functionObjects::regionSizeDistribution::writeGraphs
324 const word& fieldName,
326 const regionSplit& regions,
332 const coordSet& coords
336 Map<scalar> regionField(
regionSum(regions, cellField));
368 alphaName_(
dict.get<
word>(
"field")),
370 isoPlanes_(
dict.getOrDefault(
"isoPlanes", false))
383 dict.readEntry(
"nBins", nBins_);
384 dict.readEntry(
"field", alphaName_);
385 dict.readEntry(
"threshold", threshold_);
386 dict.readEntry(
"maxDiameter", maxDiam_);
388 dict.readIfPresent(
"minDiameter", minDiam_);
389 dict.readEntry(
"patches", patchNames_);
390 dict.readEntry(
"fields", fields_);
396 dict.subOrEmptyDict(
"formatOptions").optionalSubDict(
setFormat)
399 if (
dict.found(coordinateSystem::typeName_()))
406 Info<<
"Transforming all vectorFields with coordinate system "
407 << csysPtr_->name() <<
endl;
416 dict.readEntry(
"origin", origin_);
417 dict.readEntry(
"direction", direction_);
418 dict.readEntry(
"maxD", maxDiameter_);
419 dict.readEntry(
"nDownstreamBins", nDownstreamBins_);
420 dict.readEntry(
"maxDownstream", maxDownstream_);
421 direction_.normalise();
442 Log <<
" Looking up field " << alphaName_ <<
endl;
446 Info<<
" Reading field " << alphaName_ <<
endl;
454 mesh_.time().timeName(),
463 const auto&
alpha = talpha();
465 Log <<
" Volume of alpha = "
469 const scalar meshVol =
gSum(mesh_.V());
470 const scalar maxDropletVol = 1.0/6.0*
pow3(maxDiam_);
471 const scalar
delta = (maxDiam_-minDiam_)/nBins_;
473 Log <<
" Mesh volume = " << meshVol <<
nl
474 <<
" Maximum droplet diameter = " << maxDiam_ <<
nl
475 <<
" Maximum droplet volume = " << maxDropletVol
480 boolList blockedFace(mesh_.nFaces(),
false);
484 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
486 scalar ownVal =
alpha[mesh_.faceOwner()[facei]];
487 scalar neiVal =
alpha[mesh_.faceNeighbour()[facei]];
491 (ownVal < threshold_ && neiVal > threshold_)
492 || (ownVal > threshold_ && neiVal < threshold_)
495 blockedFace[facei] =
true;
508 const auto& ownFld = townFld();
509 const auto& nbrFld = tnbrFld();
515 scalar ownVal = ownFld[i];
516 scalar neiVal = nbrFld[i];
520 (ownVal < threshold_ && neiVal > threshold_)
521 || (ownVal > threshold_ && neiVal < threshold_)
524 blockedFace[start+i] =
true;
536 <<
" disconnected regions" <<
endl;
546 mesh_.time().timeName(),
554 Info<<
" Dumping region as volScalarField to " << region.
name()
559 region[celli] = regions[celli];
567 Map<label> patchRegions(findPatchRegions(regions));
590 scalar meshSumVol = 0.0;
591 scalar alphaSumVol = 0.0;
594 auto vIter = allRegionVolume.
cbegin();
595 auto aIter = allRegionAlphaVolume.
cbegin();
596 auto numIter = allRegionNumCells.
cbegin();
600 vIter.good() && aIter.good();
601 ++vIter, ++aIter, ++numIter
610 meshSumVol += vIter();
611 alphaSumVol += aIter();
624 Info<<
" Patch connected regions (liquid core):" <<
nl;
632 const label regioni = iter.key();
643 Info<<
" Background regions:" <<
nl;
649 auto vIter = allRegionVolume.
cbegin();
650 auto aIter = allRegionAlphaVolume.
cbegin();
655 vIter.good() && aIter.good();
661 !patchRegions.
found(vIter.key())
662 && vIter() >= maxDropletVol
681 writeAlphaFields(regions, patchRegions, allRegionVolume,
alpha);
692 const label regioni = vIter.key();
695 patchRegions.
found(regioni)
696 || vIter() >= maxDropletVol
697 || (allRegionAlphaVolume[regioni]/vIter() < threshold_)
700 allRegionVolume.
erase(vIter);
701 allRegionAlphaVolume.
erase(regioni);
702 allRegionNumCells.
erase(regioni);
706 if (allRegionVolume.
size())
720 const coordSet coords(
"diameter",
"x", xBin,
mag(xBin));
740 (
alpha.primitiveField()*mesh_.V())
741 *(mesh_.C().primitiveField() - origin_)()
759 centroids = sortedMoment/sortedVols + origin_;
762 scalarField distToPlane((centroids - origin_) & direction_);
766 (centroids - origin_) - (distToPlane*direction_)
769 const scalar deltaX = maxDownstream_/nDownstreamBins_;
775 (
mag(radialDistToOrigin[i]) < maxDiameter_)
776 && (distToPlane[i] < maxDownstream_)
779 downstreamIndices[i] = distToPlane[i]/deltaX;
786 if (downstreamIndices[i] != -1)
788 binDownCount[downstreamIndices[i]] += 1.0;
799 scalar
x = 0.5*deltaX;
807 const coordSet coords(
"distance",
"x", xBin,
mag(xBin));
809 auto&
writer = formatterPtr_();
818 writer.write(
"isoPlanes", binDownCount);
825 Info<<
" Iso-planes Bins:" <<
nl
832 forAll(binDownCount, bini)
846 forAll(sortedDiameters, i)
857 forAll(sortedDiameters, i)
859 indices[i] = (sortedDiameters[i]-minDiam_)/
delta;
864 forAll(sortedDiameters, i)
866 binCount[indices[i]] += 1.0;
872 auto&
writer = formatterPtr_();
881 writer.write(
"count", binCount);
927 Log <<
" Scalar field " << fldName <<
endl;
933 const auto&
fld = tfld();
958 Log <<
" Vector field " << fldName <<
endl;
967 Log <<
"Transforming vector field " << fldName
968 <<
" with coordinate system "
972 tfld = csysPtr_->localVector(tfld());
974 const auto&
fld = tfld();
983 alphaVol*
fld.component(cmpt),
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
vtk::internalMeshWriter writer(topoMesh, topoCells, vtk::formatType::INLINE_ASCII, runTime.path()/"blockTopology")
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Generic templated field type.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correctBoundaryConditions()
Correct boundary field.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
const_iterator cbegin() const
const_iterator set to the beginning of the HashTable
bool found(const Key &key) const
Return true if hashed entry is found in table.
label size() const noexcept
The number of elements in table.
bool erase(const iterator &iter)
Erase an entry specified by given iterator.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
const word & name() const noexcept
Return the object name.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A HashTable to objects of type <T> with a label key.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nProcs() const noexcept
Number of ranks associated with PstreamBuffers.
static void mapCombineAllGather(const List< commsStruct > &comms, Container &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
virtual bool read()
Re-read model coefficients if they have changed.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
static word timeName(const scalar t, const int precision=precision_)
void size(const label n)
Older name for setAddressableSize.
static word suffix(const word &fldName, const word &fileExt=word::null)
Name suffix based on fieldName (underscore separator)
Holds list of sampling positions.
const word & name() const noexcept
The coord-set name.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the natural logarithm of an input volScalarField.
Computes the magnitude of an input field.
const objectRegistry & obr_
Reference to the region objectRegistry.
Creates a droplet size distribution via interrogating a continuous phase fraction field.
virtual bool execute()
Do nothing.
virtual bool write()
Calculate the regionSizeDistribution and write.
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
Base class for writing single files from the function objects.
fileName baseTimeDir() const
Return the base directory for the current time value.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
virtual bool coupled() const
Return true if this patch field is coupled.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
const fvPatch & patch() const
Return patch.
const polyPatch & patch() const
Return the polyPatch.
const Time & time() const noexcept
Return time registry.
static const char *const componentNames[]
static constexpr direction nComponents
Number of components in bool is 1.
label start() const
Return start label of this patch in the polyMesh face list.
int myProcNo() const noexcept
Return processor number.
virtual bool write(const bool valid=true) const
Write using setting from DB.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
label nRegions() const
Return total number of regions.
splitCell * master() const
A class for managing temporary objects.
void reset(tmp< T > &&other) noexcept
Clear existing and transfer ownership.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
word outputName("finiteArea-edges.obj")
Volume integrate volField creating a volField.
constexpr scalar pi(M_PI)
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
List< word > wordList
A List of words.
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
static Map< Type > regionSum(const regionSplit ®ions, const Field< Type > &fld)
Field< label > labelField
Specialisation of Field<T> for label.
static constexpr const zero Zero
Global zero (0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
static List< Type > extractData(const labelUList &keys, const Map< Type > ®ionData)
dimensionedScalar cbrt(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.
#define forAllIters(container, iter)
Iterate across all elements in the container object.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
word setFormat(propsDict.getOrDefault< word >("setFormat", "vtk"))