54void Foam::functionObjects::interfaceHeight::writePositions()
60 if (
mag(direction_) > 0.0)
62 gHat = direction_/
mag(direction_);
66 gHat =
g.value()/
mag(
g.value());
72 autoPtr<interpolation<scalar>>
81 files(fileID::positionFile) <<
mesh_.time().timeName() <<
tab;
87 const midPointAndFaceSet
set
93 locations_[li] + gHat*
mesh_.bounds().mag(),
94 locations_[li] - gHat*
mesh_.bounds().mag()
98 scalar hLB =
set.size() ? - gHat & (locations_[li] -
set[0]) : - GREAT;
99 reduce(hLB, maxOp<scalar>());
104 scalar sumLength = 0, sumLengthAlpha = 0;
105 for(label si = 0; si <
set.size() - 1; ++ si)
107 if (
set.segments()[si] !=
set.segments()[si+1])
113 const label c0 =
set.cells()[si],
c1 =
set.cells()[si+1];
114 const label f0 =
set.faces()[si], f1 =
set.faces()[si+1];
115 const scalar
a0 = interpolator->interpolate(
p0, c0, f0);
116 const scalar a1 = interpolator->interpolate(p1, c1, f1);
118 const scalar l = - gHat & (p1 -
p0);
120 sumLengthAlpha += l*(
a0 + a1)/2;
123 reduce(sumLength, sumOp<scalar>());
124 reduce(sumLengthAlpha, sumOp<scalar>());
131 liquid_ ? sumLengthAlpha : sumLength - sumLengthAlpha;
132 const scalar hIL = hIB - hLB;
135 const point p = locations_[li] - gHat*hIL;
139 files(fileID::heightFile) << w << hIB << w << hIL;
140 files(fileID::positionFile) <<
'(' << w <<
p.x() << w <<
p.y()
147 files(fileID::heightFile).endl();
169 case fileID::heightFile:
175 "Interface height above the boundary"
181 "Interface height above the location"
185 case fileID::positionFile:
187 writeHeaderValue(files(fid),
"p",
"Interface position");
192 const Foam::Omanip<int> w = valueWidth(1);
194 writeCommented(files(fid),
"Location");
199 case fileID::heightFile:
200 files(fid) << w << li << w <<
' ';
202 case fileID::positionFile:
203 files(fid) << w << li << w <<
' ' << w <<
' ' <<
" ";
209 writeCommented(files(fid),
"Time");
214 case fileID::heightFile:
215 files(fid) << w <<
"hB" << w <<
"hL";
217 case fileID::positionFile:
218 files(fid) << w <<
"p" << w <<
' ' << w <<
' ' <<
" ";
236 logFiles(obr_,
name),
239 interpolationScheme_(
"cellPoint"),
258 dict.readIfPresent(
"alpha", alphaName_);
259 dict.readIfPresent(
"liquid", liquid_);
260 dict.readEntry(
"locations", locations_);
261 dict.readIfPresent(
"interpolationScheme", interpolationScheme_);
262 dict.readIfPresent(
"direction", direction_);
Istream and Ostream manipulators taking arguments.
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const uniformDimensionedVectorField & g
An Ostream manipulator taking arguments.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
This function object reports the height of the interface above a set of locations.
interfaceHeight(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
PtrList< OFstream > & files()
Inherit logFiles methods.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
virtual void writeFileHeader(const fileID fid)
Output file header information.
virtual bool execute()
Execute the function-object operations.
OFstream & files(const fileID fid)
Return file corresponding to enumeration.
virtual bool write()
Write the function-object results.
virtual bool end()
Execute at the final time-loop.
virtual void resetNames(const wordList &names)
Reset the list of names from a wordList.
virtual bool write()
Write function.
const objectRegistry & obr_
Reference to the region objectRegistry.
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
Omanip< int > valueWidth(const label offset=0) const
Return the value width when writing to stream with optional offset.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
const Time & time() const
Return the top-level database.
static autoPtr< interpolation< Type > > New(const word &interpolationType, const GeometricField< Type, fvPatchField, volMesh > &psi)
Return a reference to the specified interpolation scheme.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & p0
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
const dimensionedScalar a0
Bohr radius: default SI units: [m].
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
vector point
Point is a vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
constexpr char tab
The tab '\t' character(0x09).
#define forAll(list, i)
Loop across all elements in list.