34 namespace functionObjects
47 template<
class T,
unsigned Size>
68 void Foam::functionObjects::regionSizeDistribution::writeAlphaFields
70 const regionSplit& regions,
71 const Map<label>& patchRegions,
72 const Map<scalar>& regionVolume,
76 const scalar maxDropletVol = 1.0/6.0*
pow(maxDiam_, 3);
91 alphaName_ +
"_liquidCore",
104 alphaName_ +
"_background",
117 label regionI = regions[celli];
118 if (patchRegions.found(regionI))
120 backgroundAlpha[celli] = 0;
124 liquidCore[celli] = 0;
126 scalar regionVol = regionVolume[regionI];
127 if (regionVol < maxDropletVol)
129 backgroundAlpha[celli] = 0;
133 liquidCore.correctBoundaryConditions();
134 backgroundAlpha.correctBoundaryConditions();
136 Info<<
" Volume of liquid-core = "
139 Info<<
" Volume of background = "
143 Info<<
" Writing liquid-core field to " << liquidCore.name() <<
endl;
145 Info<<
" Writing background field to " << backgroundAlpha.name() <<
endl;
146 backgroundAlpha.
write();
151 Foam::functionObjects::regionSizeDistribution::findPatchRegions
153 const regionSplit& regions
160 const labelHashSet patchIDs(mesh_.boundaryMesh().patchSet(patchNames_));
162 label nPatchFaces = 0;
165 nPatchFaces += mesh_.boundaryMesh()[iter.key()].size();
169 Map<label> patchRegions(nPatchFaces);
172 const polyPatch& pp = mesh_.boundaryMesh()[iter.key()];
175 const labelList& faceCells = pp.faceCells();
181 regions[faceCells[i]],
198 Foam::functionObjects::regionSizeDistribution::divide
211 result[i] = num[i]/denom[i];
224 void Foam::functionObjects::regionSizeDistribution::generateFields
226 const word& fieldName,
240 binSum[indices[i]] += sortedField[i];
256 void Foam::functionObjects::regionSizeDistribution::generateFields
258 const word& fieldName,
263 PtrList<scalarField>&
fields
272 binSum[indices[i]] += sortedField[i];
282 binSqrSum[indices[i]] +=
Foam::sqr(sortedField[i]);
302 void Foam::functionObjects::regionSizeDistribution::generateFields
304 const word& fieldName,
305 const Field<Type>& cellField,
306 const regionSplit& regions,
312 PtrList<Field<Type>>&
fields
316 Map<Type> regionField(regionSum(regions, cellField));
319 Field<Type> sortedField
321 sortedNormalisation*extractData(sortedRegions, regionField)
348 alphaName_(
dict.lookup(
"alpha")),
349 patchNames_(
dict.lookup(
"patches"))
365 dict.lookup(
"alpha") >> alphaName_;
366 dict.lookup(
"patches") >> patchNames_;
367 dict.lookup(
"threshold") >> threshold_;
368 dict.lookup(
"maxDiameter") >> maxDiam_;
370 dict.readIfPresent(
"minDiameter", minDiam_);
371 dict.lookup(
"nBins") >> nBins_;
372 dict.lookup(
"fields") >> fields_;
383 fields.append(alphaName_);
401 Info<<
" Looking up field " << alphaName_ <<
endl;
405 Info<<
" Reading field " << alphaName_ <<
endl;
431 Info<<
" Volume of alpha = "
435 const scalar meshVol =
gSum(mesh_.V());
436 const scalar maxDropletVol = 1.0/6.0*
pow(maxDiam_, 3);
437 const scalar
delta = (maxDiam_-minDiam_)/nBins_;
439 Info<<
" Mesh volume = " << meshVol <<
endl;
440 Info<<
" Maximum droplet diameter = " << maxDiam_ <<
endl;
441 Info<<
" Maximum droplet volume = " << maxDropletVol <<
endl;
445 boolList blockedFace(mesh_.nFaces(),
false);
448 for (
label facei = 0; facei < mesh_.nInternalFaces(); facei++)
450 scalar ownVal =
alpha[mesh_.faceOwner()[facei]];
451 scalar neiVal =
alpha[mesh_.faceNeighbour()[facei]];
455 (ownVal < threshold_ && neiVal > threshold_)
456 || (ownVal > threshold_ && neiVal < threshold_)
459 blockedFace[facei] =
true;
478 scalar ownVal = ownFld[i];
479 scalar neiVal = nbrFld[i];
483 (ownVal < threshold_ && neiVal > threshold_)
484 || (ownVal > threshold_ && neiVal < threshold_)
487 blockedFace[start+i] =
true;
498 <<
" disconnected regions" <<
endl;
516 Info<<
" Dumping region as volScalarField to " << region.
name()
521 region[celli] = regions[celli];
529 Map<label> patchRegions(findPatchRegions(regions));
534 Map<scalar> allRegionVolume(regionSum(regions, mesh_.V()));
535 Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
548 <<
tab <<
"Volume(mesh)"
549 <<
tab <<
"Volume(" <<
alpha.name() <<
"):"
552 scalar meshSumVol = 0.0;
553 scalar alphaSumVol = 0.0;
562 vIter != allRegionVolume.
end()
563 && aIter != allRegionAlphaVolume.
end();
564 ++vIter, ++aIter, ++numIter
567 Info<<
" " <<
tab << vIter.key()
573 meshSumVol += vIter();
574 alphaSumVol += aIter();
579 <<
tab << alphaSumVol
589 Info<<
" Patch connected regions (liquid core):" <<
endl;
591 <<
tab <<
"Volume(mesh)"
592 <<
tab <<
"Volume(" <<
alpha.name() <<
"):"
596 label regionI = iter.key();
597 Info<<
" " <<
tab << iter.key()
598 <<
tab << allRegionVolume[regionI]
599 <<
tab << allRegionAlphaVolume[regionI] <<
endl;
606 Info<<
" Background regions:" <<
endl;
608 <<
tab <<
"Volume(mesh)"
609 <<
tab <<
"Volume(" <<
alpha.name() <<
"):"
617 vIter != allRegionVolume.
end()
618 && aIter != allRegionAlphaVolume.
end();
624 !patchRegions.
found(vIter.key())
625 && vIter() >= maxDropletVol
628 Info<<
" " <<
tab << vIter.key()
644 writeAlphaFields(regions, patchRegions, allRegionVolume,
alpha);
653 label regionI = vIter.key();
656 patchRegions.
found(regionI)
657 || vIter() >= maxDropletVol
660 allRegionVolume.
erase(vIter);
661 allRegionAlphaVolume.
erase(regionI);
662 allRegionNumCells.
erase(regionI);
666 if (allRegionVolume.
size())
693 forAll(sortedDiameters, i)
704 forAll(sortedDiameters, i)
706 indices[i] = (sortedDiameters[i]-minDiam_)/
delta;
711 forAll(sortedDiameters, i)
713 binCount[indices[i]] += 1.0;
720 <<
tab <<
"Min diameter"
729 <<
tab << binCount[binI] <<
endl;
737 #define DeclareTypeFields(Type, nullArg) \
738 PtrList<Field<Type>> Type##Fields;
740 #undef DeclareTypeFields
744 #define TypeFieldsAppend(Type, nullArg) \
745 appendFields(binCount, Type##Fields);
746 #undef TypeFieldsAppend
764 #define GenerateTypeFields(Type, nullArg) \
766 if (obr_.foundObject<VolField<Type>>(fields_[fieldi])) \
770 const VolField<Type>& field = \
771 obr_.lookupObject<VolField<Type>>(fields_[fieldi]); \
776 (alphaVol*field)(), \
787 #undef GenerateTypeFields
789 if (!
found) cannotFindObject(fields_[fieldi]);
793 #define TypeFieldsExpand(Type, nullArg) \
794 Type##Fields.setSize(fieldNames.size());
796 #undef TypeFieldsAppend
799 formatterPtr_().write
#define forAll(list, i)
Loop across all elements in list.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Macros for easy insertion into run-time selection tables.
Pre-declare SubField and related Field type.
A 1D vector of objects of type <T> with a fixed size <Size>.
Generic GeometricField class.
void correctBoundaryConditions()
Correct boundary field.
List< Key > sortedToc() const
Return the table of contents as a sorted list.
bool erase(const iterator &)
Erase a hashedEntry specified by given iterator.
iterator begin()
Iterator set to the beginning of the HashTable.
label size() const
Return number of elements in table.
bool found(const Key &) const
Return true if hashedEntry is found in table.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
const word & name() const
Return name.
Plus op for FixedList<scalar>
void append(const T &)
Append an element at the end of the list.
void size(const label)
Override size to be inconsistent with allocated storage.
A HashTable to objects of type <T> with a label key.
virtual Ostream & write(const char)=0
Write character.
static void mapCombineScatter(const List< commsStruct > &comms, Container &Values, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void mapCombineGather(const List< commsStruct > &comms, Container &Values, const CombineOp &cop, const int tag, const label comm)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static bool master(const label communicator=0)
Am I the master process.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Holds list of sampling positions.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const word & name() const
Return const reference to name.
Abstract base-class for Time/database functionObjects.
const Time & time_
Reference to time.
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const objectRegistry & obr_
Reference to the region objectRegistry.
Creates a size distribution via interrogating a continuous phase fraction field.
regionSizeDistribution(const word &name, const Time &runTime, const dictionary &)
Construct for given objectRegistry and dictionary.
virtual wordList fields() const
Return the list of fields required.
virtual ~regionSizeDistribution()
virtual bool execute()
Do nothing.
virtual bool write()
Calculate the regionSizeDistribution and write.
virtual bool read(const dictionary &)
Read the regionSizeDistribution data.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
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.
virtual tmp< Field< Type > > patchNeighbourField(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) const
Return patchField on the opposite patch of a coupled patch.
const fvPatch & patch() const
Return patch.
static const word & calculatedType()
Return the type of the calculated for of fvPatchField.
const polyPatch & patch() const
Return the polyPatch.
label start() const
Return start label of this patch in the polyMesh face list.
virtual bool write(const bool write=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.
static autoPtr< setWriter > New(const word &writeType, const IOstream::streamFormat writeFormat=IOstream::ASCII, const IOstream::compressionType writeCompression=IOstream::UNCOMPRESSED)
Select given write options.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
Volume integrate volField creating a volField.
static List< word > fieldNames
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
dimensioned< Type > domainIntegrate(const VolField< Type > &vf)
List< word > wordList
A List of words.
Type gSum(const FieldField< Field, Type > &f)
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
const dimensionSet dimless
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Field< label > labelField
Specialisation of Field<T> for label.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
dimensionedScalar cbrt(const dimensionedScalar &ds)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
#define GenerateTypeFields(Type, nullArg)
#define TypeFieldsParameter(Type, nullArg)
#define DeclareTypeFields(Type, nullArg)
#define TypeFieldsExpand(Type, nullArg)
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable