Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
regionSizeDistribution Class Reference

This function object creates a size distribution via interrogating a continuous phase fraction field. More...

Inheritance diagram for regionSizeDistribution:
Inheritance graph
[legend]
Collaboration diagram for regionSizeDistribution:
Collaboration graph
[legend]

Public Member Functions

 TypeName ("regionSizeDistribution")
 Runtime type information. More...
 
 regionSizeDistribution (const word &name, const objectRegistry &, const dictionary &, const bool loadFromFiles=false)
 Construct for given objectRegistry and dictionary. More...
 
virtual ~regionSizeDistribution ()
 
virtual const wordname () const
 Return name of the set of regionSizeDistribution. More...
 
virtual void read (const dictionary &)
 Read the regionSizeDistribution data. More...
 
virtual void execute ()
 Execute, currently does nothing. More...
 
virtual void end ()
 Execute at the final time-loop, currently does nothing. More...
 
virtual void timeSet ()
 Called when time was set at the end of the Time::operator++. More...
 
virtual void write ()
 Calculate the regionSizeDistribution and write. More...
 
virtual void updateMesh (const mapPolyMesh &)
 Update for changes of mesh. More...
 
virtual void movePoints (const polyMesh &)
 Update for changes of mesh. More...
 
template<class Type >
Foam::Map< Type > regionSum (const regionSplit &regions, const Field< Type > &fld) const
 
template<class Type >
Foam::List< Type > extractData (const UList< label > &keys, const Map< Type > &regionData) const
 
- Public Member Functions inherited from functionObjectFile
 functionObjectFile (const objectRegistry &obr, const word &prefix)
 Construct null. More...
 
 functionObjectFile (const objectRegistry &obr, const word &prefix, const word &fileName, const dictionary &dict)
 Construct from components and read options from dictionary. More...
 
virtual ~functionObjectFile ()
 Destructor. More...
 
void read (const dictionary &dict)
 Read. More...
 
OFstreamfile ()
 Return access to the file (if only 1) More...
 
bool writeToFile () const
 Return true if can write to file. More...
 
void writeCommented (Ostream &os, const string &str) const
 Write a commented string to stream. More...
 
void writeTabbed (Ostream &os, const string &str) const
 Write a tabbed string to stream. More...
 
void writeHeader (Ostream &os, const string &str) const
 Write a commented header to stream. More...
 
void writeTime (Ostream &os) const
 Write the current time to stream. More...
 
template<class Type >
void writeHeaderValue (Ostream &os, const string &property, const Type &value) const
 Write a (commented) header property and value pair. More...
 
label charWidth () const
 Return width of character stream output. More...
 

Private Member Functions

template<class Type >
Map< Type > regionSum (const regionSplit &, const Field< Type > &) const
 
template<class Type >
List< Type > extractData (const UList< label > &keys, const Map< Type > &) const
 Get data in order. More...
 
void writeGraph (const coordSet &coords, const word &valueName, const scalarField &values) const
 
void writeAlphaFields (const regionSplit &regions, const Map< label > &keepRegions, const Map< scalar > &regionVolume, const volScalarField &alpha) const
 Write volfields with the parts of alpha which are not. More...
 
Map< labelfindPatchRegions (const polyMesh &, const regionSplit &) const
 Mark all regions starting at patches. More...
 
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. More...
 
void writeGraphs (const word &fieldName, const scalarField &cellField, const regionSplit &regions, const labelList &sortedRegions, const scalarField &sortedNormalisation, const labelList &indices, const scalarField &binCount, const coordSet &coords) const
 Given per-cell data calculate per-bin average/deviation and graph. More...
 
 regionSizeDistribution (const regionSizeDistribution &)
 Disallow default bitwise copy construct. More...
 
void operator= (const regionSizeDistribution &)
 Disallow default bitwise assignment. More...
 

Static Private Member Functions

static tmp< scalarFielddivide (const scalarField &, const scalarField &)
 Helper: divide if denom != 0. More...
 

Private Attributes

word name_
 Name of this set of regionSizeDistribution objects. More...
 
const objectRegistryobr_
 
bool active_
 on/off switch More...
 
word alphaName_
 Name of field. More...
 
wordReList patchNames_
 Patches to walk from. More...
 
Switch log_
 Switch to send output to Info as well as to file. More...
 
scalar threshold_
 Clip value. More...
 
scalar maxDiam_
 Maximum droplet diameter. More...
 
scalar minDiam_
 Minimum droplet diameter. More...
 
label nBins_
 Mumber of bins. More...
 
wordReList fields_
 Names of fields to sample on regions. More...
 
autoPtr< writer< scalar > > formatterPtr_
 Output formatter to write. More...
 
autoPtr< coordinateSystemcoordSysPtr_
 Optional coordinate system. More...
 

Additional Inherited Members

- Static Public Attributes inherited from functionObjectFile
static const word outputPrefix = "postProcessing"
 Directory prefix. More...
 
static label addChars = 7
 Additional characters for writing. More...
 
- Protected Member Functions inherited from functionObjectFile
virtual void initStream (Ostream &os) const
 Initialise the output stream for writing. More...
 
virtual fileName baseFileDir () const
 Return the base directory for output. More...
 
virtual fileName baseTimeDir () const
 Return the base directory for the current time value. More...
 
virtual autoPtr< OFstreamcreateFile (const word &name) const
 Return an autoPtr to a new file. More...
 
virtual void resetFile (const word &name)
 Reset internal file pointer to new file with new name. More...
 
virtual Omanip< int > valueWidth (const label offset=0) const
 Return the value width when writing to stream with optional offset. More...
 
 functionObjectFile (const functionObjectFile &)
 Disallow default bitwise copy construct. More...
 
void operator= (const functionObjectFile &)
 Disallow default bitwise assignment. More...
 
- Protected Attributes inherited from functionObjectFile
bool writeToFile_
 Flag to enable/disable writing to file. More...
 

Detailed Description

This function object creates a size distribution via interrogating a continuous phase fraction field.

Looks up a phase-fraction (alpha) field and splits the mesh into regions based on where the field is below the threshold value. These regions ("droplets") can now be analysed.

Regions:

Output (volume scalar) fields include:

Histogram:

Example of function object specification:

regionSizeDistribution1
{
    type            regionSizeDistribution;
    functionObjectLibs ("libfieldFunctionObjects.so");
    ...
    field           alpha;
    patches         (inlet);
    threshold       0.4;
    fields          (p U);
    nBins           100;
    maxDiameter     0.5e-4;
    minDiameter     0;
    setFormat       gnuplot;
    coordinateSystem
    {
        type            cartesian;
        origin          (0 0 0);
        e3              (0 1 1);
        e1              (1 0 0);
    }
}


Function object usage

Property Description Required Default value
type type name: regionSizeDistribution yes
field phase field to interrogate yes
patches patches from which the liquid core is identified yes
threshold phase fraction applied to delimit regions yes
fields fields to sample yes
nBins number of bins for histogram yes
maxDiameter maximum region equivalent diameter yes
minDiameter minimum region equivalent diameter no 0
setFormat writing format yes
coordinateSystem transformation for vector fields no
log Log to standard output no yes
See also
Foam::functionObject Foam::OutputFilterFunctionObject
Source files

Definition at line 194 of file regionSizeDistribution.H.

Constructor & Destructor Documentation

◆ regionSizeDistribution() [1/2]

Disallow default bitwise copy construct.

◆ regionSizeDistribution() [2/2]

regionSizeDistribution ( const word name,
const objectRegistry obr,
const dictionary dict,
const bool  loadFromFiles = false 
)

Construct for given objectRegistry and dictionary.

Allow the possibility to load fields from files

Definition at line 325 of file regionSizeDistribution.C.

References dict, Foam::endl(), Foam::nl, Foam::read(), and WarningInFunction.

Here is the call graph for this function:

◆ ~regionSizeDistribution()

~regionSizeDistribution ( )
virtual

Definition at line 357 of file regionSizeDistribution.C.

Member Function Documentation

◆ regionSum() [1/2]

Map<Type> regionSum ( const regionSplit ,
const Field< Type > &   
) const
private

◆ extractData() [1/2]

List<Type> extractData ( const UList< label > &  keys,
const Map< Type > &   
) const
private

Get data in order.

◆ writeGraph()

void writeGraph ( const coordSet coords,
const word valueName,
const scalarField values 
) const
private

Definition at line 62 of file regionSizeDistribution.C.

References Foam::endl(), Foam::Info, Foam::mkDir(), and OFstream::name().

Here is the call graph for this function:

◆ writeAlphaFields()

void writeAlphaFields ( const regionSplit regions,
const Map< label > &  keepRegions,
const Map< scalar > &  regionVolume,
const volScalarField alpha 
) const
private

Write volfields with the parts of alpha which are not.

droplets (liquidCore, backGround)

Definition at line 86 of file regionSizeDistribution.C.

References Foam::constant::atomic::alpha, GeometricField::correctBoundaryConditions(), Foam::fvc::domainIntegrate(), Foam::endl(), forAll, Foam::Info, Foam::pow(), and dimensioned::value().

Here is the call graph for this function:

◆ findPatchRegions()

Foam::Map< Foam::label > findPatchRegions ( const polyMesh mesh,
const regionSplit regions 
) const
private

Mark all regions starting at patches.

Definition at line 174 of file regionSizeDistribution.C.

References polyPatch::faceCells(), forAll, forAllConstIter(), and mesh.

Here is the call graph for this function:

◆ divide()

Foam::tmp< Foam::scalarField > divide ( const scalarField num,
const scalarField denom 
)
staticprivate

Helper: divide if denom != 0.

Definition at line 220 of file regionSizeDistribution.C.

References forAll, and scalarField().

Here is the call graph for this function:

◆ writeGraphs() [1/2]

void writeGraphs ( const word fieldName,
const labelList indices,
const scalarField sortedField,
const scalarField binCount,
const coordSet coords 
) const
private

Given per-region data calculate per-bin average/deviation and graph.

Definition at line 244 of file regionSizeDistribution.C.

References Foam::divide(), forAll, Foam::sqr(), and Foam::sqrt().

Here is the call graph for this function:

◆ writeGraphs() [2/2]

void writeGraphs ( const word fieldName,
const scalarField cellField,
const regionSplit regions,
const labelList sortedRegions,
const scalarField sortedNormalisation,
const labelList indices,
const scalarField binCount,
const coordSet coords 
) const
private

Given per-cell data calculate per-bin average/deviation and graph.

Definition at line 285 of file regionSizeDistribution.C.

◆ operator=()

void operator= ( const regionSizeDistribution )
private

Disallow default bitwise assignment.

◆ TypeName()

TypeName ( "regionSizeDistribution"  )

Runtime type information.

◆ name()

virtual const word& name ( ) const
inlinevirtual

Return name of the set of regionSizeDistribution.

Definition at line 331 of file regionSizeDistribution.H.

References regionSizeDistribution::name_.

◆ read()

void read ( const dictionary dict)
virtual

◆ execute()

void execute ( )
virtual

Execute, currently does nothing.

Definition at line 395 of file regionSizeDistribution.C.

◆ end()

void end ( )
virtual

Execute at the final time-loop, currently does nothing.

Definition at line 401 of file regionSizeDistribution.C.

◆ timeSet()

void timeSet ( )
virtual

Called when time was set at the end of the Time::operator++.

Definition at line 407 of file regionSizeDistribution.C.

◆ write()

void write ( )
virtual

◆ updateMesh()

virtual void updateMesh ( const mapPolyMesh )
inlinevirtual

Update for changes of mesh.

Definition at line 352 of file regionSizeDistribution.H.

◆ movePoints()

virtual void movePoints ( const polyMesh )
inlinevirtual

Update for changes of mesh.

Definition at line 356 of file regionSizeDistribution.H.

◆ regionSum() [2/2]

Foam::Map<Type> regionSum ( const regionSplit regions,
const Field< Type > &  fld 
) const

Definition at line 34 of file regionSizeDistributionTemplates.C.

References fld(), forAll, and regionSplit::nRegions().

Here is the call graph for this function:

◆ extractData() [2/2]

Foam::List<Type> extractData ( const UList< label > &  keys,
const Map< Type > &  regionData 
) const

Definition at line 66 of file regionSizeDistributionTemplates.C.

References forAll, and UList::size().

Here is the call graph for this function:

Field Documentation

◆ name_

word name_
private

Name of this set of regionSizeDistribution objects.

Definition at line 201 of file regionSizeDistribution.H.

Referenced by regionSizeDistribution::name().

◆ obr_

const objectRegistry& obr_
private

Definition at line 203 of file regionSizeDistribution.H.

◆ active_

bool active_
private

on/off switch

Definition at line 206 of file regionSizeDistribution.H.

◆ alphaName_

word alphaName_
private

Name of field.

Definition at line 209 of file regionSizeDistribution.H.

◆ patchNames_

wordReList patchNames_
private

Patches to walk from.

Definition at line 212 of file regionSizeDistribution.H.

◆ log_

Switch log_
private

Switch to send output to Info as well as to file.

Definition at line 215 of file regionSizeDistribution.H.

◆ threshold_

scalar threshold_
private

Clip value.

Definition at line 218 of file regionSizeDistribution.H.

◆ maxDiam_

scalar maxDiam_
private

Maximum droplet diameter.

Definition at line 221 of file regionSizeDistribution.H.

◆ minDiam_

scalar minDiam_
private

Minimum droplet diameter.

Definition at line 224 of file regionSizeDistribution.H.

◆ nBins_

label nBins_
private

Mumber of bins.

Definition at line 227 of file regionSizeDistribution.H.

◆ fields_

wordReList fields_
private

Names of fields to sample on regions.

Definition at line 230 of file regionSizeDistribution.H.

◆ formatterPtr_

autoPtr<writer<scalar> > formatterPtr_
private

Output formatter to write.

Definition at line 233 of file regionSizeDistribution.H.

◆ coordSysPtr_

autoPtr<coordinateSystem> coordSysPtr_
private

Optional coordinate system.

Definition at line 236 of file regionSizeDistribution.H.


The documentation for this class was generated from the following files: