Go to the documentation of this file.
40 template<
class T,
unsigned Size>
64 const word& valueName,
68 const wordList valNames(1, valueName);
73 OFstream str(outputPath/formatterPtr_().getFileName(coords, valNames));
76 <<
"Writing distribution of " << valueName <<
" to " << str.
name()
81 formatterPtr_().write(coords, valNames, valPtrs, str);
93 const scalar maxDropletVol = 1.0/6.0*
pow(maxDiam_, 3);
108 alphaName_ +
"_liquidCore",
109 obr_.time().timeName(),
121 alphaName_ +
"_background",
122 obr_.time().timeName(),
134 label regionI = regions[cellI];
135 if (patchRegions.found(regionI))
137 backgroundAlpha[cellI] = 0;
141 liquidCore[cellI] = 0;
143 scalar regionVol = regionVolume[regionI];
144 if (regionVol < maxDropletVol)
146 backgroundAlpha[cellI] = 0;
155 Info<<
" Volume of liquid-core = "
158 Info<<
" Volume of background = "
164 <<
"Writing liquid-core field to " << liquidCore.name() <<
endl;
168 <<
"Writing background field to " << backgroundAlpha.name() <<
endl;
169 backgroundAlpha.write();
185 label nPatchFaces = 0;
188 nPatchFaces +=
mesh.boundaryMesh()[iter.key()].size();
204 regions[faceCells[i]],
213 Pstream::mapCombineScatter(patchRegions);
232 result[i] = num[i]/denom[i];
245 const word& fieldName,
252 if (Pstream::master())
258 binSum[indices[i]] += sortedField[i];
267 binSqrSum[indices[i]] +=
Foam::sqr(sortedField[i]);
275 writeGraph(coords, fieldName +
"_sum", binSum);
277 writeGraph(coords, fieldName +
"_avg", binAvg);
279 writeGraph(coords, fieldName +
"_dev", binDev);
286 const word& fieldName,
298 Map<scalar> regionField(regionSum(regions, cellField));
329 const bool loadFromFiles
336 alphaName_(
dict.lookup(
"field")),
337 patchNames_(
dict.lookup(
"patches")),
341 if (isA<fvMesh>(obr_))
349 <<
"No fvMesh available, deactivating " << name_ <<
nl
369 log_.readIfPresent(
"log",
dict);
388 <<
"Transforming all vectorFields with coordinate system "
389 << coordSysPtr_().name() <<
endl;
417 if (log_)
Info <<
type() <<
" " << name_ <<
" output:" <<
nl;
419 const fvMesh&
mesh = refCast<const fvMesh>(obr_);
424 if (log_)
Info <<
" Looking up field " << alphaName_ <<
endl;
428 if (log_)
Info <<
" Reading field " << alphaName_ <<
endl;
455 <<
" Volume of alpha = "
460 const scalar maxDropletVol = 1.0/6.0*
pow(maxDiam_, 3);
461 const scalar
delta = (maxDiam_-minDiam_)/nBins_;
464 <<
" Mesh volume = " << meshVol <<
nl
465 <<
" Maximum droplet diameter = " << maxDiam_ <<
nl
466 <<
" Maximum droplet volume = " << maxDropletVol
482 (ownVal < threshold_ && neiVal > threshold_)
483 || (ownVal > threshold_ && neiVal < threshold_)
486 blockedFace[faceI] =
true;
506 scalar ownVal = ownFld[i];
507 scalar neiVal = nbrFld[i];
511 (ownVal < threshold_ && neiVal > threshold_)
512 || (ownVal > threshold_ && neiVal < threshold_)
515 blockedFace[start+i] =
true;
527 <<
" Determined " << regions.
nRegions()
528 <<
" disconnected regions" <<
endl;
548 <<
" Dumping region as " << volScalarField::typeName
549 <<
" to " << region.name() <<
endl;
553 region[cellI] = regions[cellI];
568 Map<scalar> allRegionAlphaVolume(regionSum(regions, alphaVol));
587 scalar meshSumVol = 0.0;
588 scalar alphaSumVol = 0.0;
597 vIter != allRegionVolume.end()
598 && aIter != allRegionAlphaVolume.end();
599 ++vIter, ++aIter, ++numIter
609 meshSumVol += vIter();
610 alphaSumVol += aIter();
626 Info<<
" Patch connected regions (liquid core):" <<
nl
634 label regionI = iter.key();
647 Info<<
" Background regions:" <<
nl
659 vIter != allRegionVolume.end()
660 && aIter != allRegionAlphaVolume.end();
666 !patchRegions.found(vIter.
key())
667 && vIter() >= maxDropletVol
687 writeAlphaFields(regions, patchRegions, allRegionVolume,
alpha);
696 label regionI = vIter.key();
699 patchRegions.found(regionI)
700 || vIter() >= maxDropletVol
703 allRegionVolume.erase(vIter);
704 allRegionAlphaVolume.erase(regionI);
705 allRegionNumCells.erase(regionI);
709 if (allRegionVolume.size())
721 const coordSet coords(
"diameter",
"x", xBin,
mag(xBin));
727 const labelList sortedRegions = allRegionAlphaVolume.sortedToc();
740 forAll(sortedDiameters, i)
750 labelList indices(sortedDiameters.size());
751 forAll(sortedDiameters, i)
753 indices[i] = (sortedDiameters[i]-minDiam_)/
delta;
758 forAll(sortedDiameters, i)
760 binCount[indices[i]] += 1.0;
766 writeGraph(coords,
"count", binCount);
804 wordList scalarNames(obr_.names(volScalarField::typeName));
809 const word& fldName = scalarNames[selected[i]];
810 if (log_)
Info <<
" Scalar field " << fldName <<
endl;
833 wordList vectorNames(obr_.names(volVectorField::typeName));
838 const word& fldName = vectorNames[selected[i]];
839 if (log_)
Info <<
" Vector field " << fldName <<
endl;
846 if (coordSysPtr_.valid())
849 <<
"Transforming vector field " << fldName
850 <<
" with coordinate system "
851 << coordSysPtr_().name() <<
endl;
853 fld = coordSysPtr_().localVector(
fld);
864 alphaVol*
fld.component(cmp),
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Map< label > findPatchRegions(const polyMesh &, const regionSplit &) const
Mark all regions starting at patches.
regionSizeDistribution(const regionSizeDistribution &)
Disallow default bitwise copy construct.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
word format(conversionProperties.lookup("format"))
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
static tmp< scalarField > divide(const scalarField &, const scalarField &)
Helper: divide if denom != 0.
A class for handling words, derived from string.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
A class for handling file names.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Plus op for FixedList<scalar>
dimensioned< Type > domainIntegrate(const GeometricField< Type, fvPatchField, volMesh > &vf)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
bool read(const char *, int32_t &)
autoPtr< volScalarField > alphaPtr
An STL-conforming const_iterator.
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const word & name() const
Return const reference to name.
@ nComponents
Number of components in this vector space.
void read(const dictionary &dict)
Read.
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
Type gSum(const FieldField< Field, Type > &f)
dimensioned< scalar > mag(const dimensioned< Type > &)
static const char * componentNames[]
virtual ~regionSizeDistribution()
void writeAlphaFields(const regionSplit ®ions, const Map< label > &keepRegions, const Map< scalar > ®ionVolume, const volScalarField &alpha) const
Write volfields with the parts of alpha which are not.
Mesh consisting of general polyhedral cells.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
void divide(FieldField< Field, Type > &f, const FieldField< Field, Type > &f1, const FieldField< Field, scalar > &f2)
Registry of regIOobjects.
virtual void execute()
Execute, currently does nothing.
virtual bool coupled() const
Return true if this patch field is coupled.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Pre-declare SubField and related Field type.
Field< label > labelField
Specialisation of Field<T> for label.
A patch is a list of labels that address the faces in the global face list.
virtual const labelList & faceOwner() const
Return face owner.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
A list of keyword definitions, which are a keyword followed by any number of values (e....
label nInternalFaces() const
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))
Generic dimensioned Type class.
Holds list of sampling positions.
Mesh data needed to do the Finite Volume discretisation.
const labelUList & faceCells() const
Return face-cell addressing.
virtual void end()
Execute at the final time-loop, currently does nothing.
label start() const
Return start label of this patch in the polyMesh face list.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
void correctBoundaryConditions()
Correct boundary field.
static bool master(const label communicator=0)
Am I the master process.
conserve internalField()+
const Key & key() const
Return the Key corresponding to the iterator.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField scalarField(fieldObject, mesh)
Base class for output file data handling.
A 1D vector of objects of type <T> with a fixed size <Size>.
virtual void read(const dictionary &)
Read the regionSizeDistribution data.
Volume integrate volField creating a volField.
const polyPatch & patch() const
Return the polyPatch.
const Time & time() const
Return the top-level database.
Operations on lists of strings.
const fileName & name() const
Return the name of the stream.
const fvPatch & patch() const
Return patch.
vector point
Point is a vector.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
dimensionedScalar cbrt(const dimensionedScalar &ds)
virtual void write()
Calculate the regionSizeDistribution and write.
Generic GeometricField class.
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
void writeGraph(const coordSet &coords, const word &valueName, const scalarField &values) const
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Report a warning using Foam::Warning.
word name(const complex &)
Return a string representation of a complex.
virtual const labelList & faceNeighbour() const
Return face neighbour.
Base class for other coordinate system specifications.
static autoPtr< writer > New(const word &writeFormat)
Return a reference to the selected writer.