40 namespace functionObjects
63 if (dict.
found(
"binData"))
66 const label nb = binDict.lookup<
label>(
"nBin");
70 names.
append(forceType +
"_bins");
83 ?
"(pressure viscous porous)" 84 :
"(pressure viscous)" 89 case fileID::mainFile:
94 writeHeaderValue(
file(i),
"CofR", coordSys_.origin());
95 writeCommented(
file(i),
"Time");
98 <<
"forces" << forceTypes <<
tab 99 <<
"moments" << forceTypes;
105 <<
"localForces" << forceTypes <<
tab 106 <<
"localMoments" << forceTypes;
111 case fileID::binsFile:
116 writeHeaderValue(
file(i),
"bins", nBin_);
117 writeHeaderValue(
file(i),
"start", binMin_);
118 writeHeaderValue(
file(i),
"delta", binDx_);
119 writeHeaderValue(
file(i),
"direction", binDir_);
122 writeCommented(
file(i),
"x co-ords :");
125 binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
126 file(i) <<
tab << binPoints[pointi].x();
130 writeCommented(
file(i),
"y co-ords :");
133 file(i) <<
tab << binPoints[pointi].y();
137 writeCommented(
file(i),
"z co-ords :");
140 file(i) <<
tab << binPoints[pointi].z();
144 writeCommented(
file(i),
"Time");
146 for (
label j = 0; j < nBin_; j++)
149 const word f(
"forces" + jn + forceTypes);
150 const word m(
"moments" + jn + forceTypes);
156 for (
label j = 0; j < nBin_; j++)
159 const word f(
"localForces" + jn + forceTypes);
160 const word m(
"localMoments" + jn + forceTypes);
171 <<
"Unhandled file index: " << i
187 if (directForceDensity_)
192 <<
"Could not find " << fDName_ <<
" in database." 206 <<
"Could not find " << UName_ <<
", " << pName_
217 <<
"Could not find " << rhoName_
234 const word& modelName = momentumTransportModel::typeName;
235 const word phaseModelName =
240 if (obr_.foundObject<icoModel>(modelName))
243 obr_.lookupObject<icoModel>(modelName);
247 else if (obr_.foundObject<cmpModel>(modelName))
249 const cmpModel& model =
250 obr_.lookupObject<cmpModel>(modelName);
252 return alpha()*model.devTau();
254 else if (obr_.foundObject<phaseIcoModel>(phaseModelName))
256 const phaseIcoModel& model =
257 obr_.lookupObject<phaseIcoModel>(phaseModelName);
259 return rho()*model.devSigma();
261 else if (obr_.foundObject<phaseCmpModel>(phaseModelName))
263 const phaseCmpModel& model =
264 obr_.lookupObject<phaseCmpModel>(phaseModelName);
266 return model.devTau();
268 else if (obr_.foundObject<
dictionary>(
"physicalProperties"))
273 obr_.lookupObject<
dictionary>(
"physicalProperties");
279 physicalProperties.
lookup(
"nu")
289 <<
"No valid model for viscous stress calculation" 304 const word& modelName = momentumTransportModel::typeName;
305 const word phaseModelName =
310 if (obr_.foundObject<icoModel>(modelName))
313 obr_.lookupObject<icoModel>(modelName);
315 return rho()*model.
nu();
317 else if (obr_.foundObject<cmpModel>(modelName))
319 const cmpModel& model =
320 obr_.lookupObject<cmpModel>(modelName);
322 return model.rho()*model.nu();
324 else if (obr_.foundObject<phaseIcoModel>(phaseModelName))
326 const phaseIcoModel& model =
327 obr_.lookupObject<phaseIcoModel>(phaseModelName);
329 return rho()*model.nu();
331 else if (obr_.foundObject<phaseCmpModel>(phaseModelName))
333 const phaseCmpModel& model =
334 obr_.lookupObject<phaseCmpModel>(phaseModelName);
336 return model.rho()*model.nu();
338 else if (obr_.foundObject<
dictionary>(
"physicalProperties"))
343 obr_.lookupObject<
dictionary>(
"physicalProperties");
349 physicalProperties.
lookup(
"nu")
357 <<
"No valid model for dynamic viscosity calculation" 367 if (rhoName_ ==
"rhoInf")
391 if (rhoName_ !=
"rhoInf")
394 <<
"Dynamic pressure is expected but kinematic is provided." 441 ).boundaryField()[
patchi];
457 force_[0][0] +=
sum(fN);
458 force_[1][0] +=
sum(fT);
459 force_[2][0] +=
sum(fP);
460 moment_[0][0] +=
sum(Md^fN);
461 moment_[1][0] +=
sum(Md^fT);
462 moment_[2][0] +=
sum(Md^fP);
470 label bini =
min(
max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
472 force_[0][bini] += fN[i];
473 force_[1][bini] += fT[i];
474 force_[2][bini] += fP[i];
475 moment_[0][bini] += Md[i]^fN[i];
476 moment_[1][bini] += Md[i]^fT[i];
477 moment_[2][bini] += Md[i]^fP[i];
486 <<
" sum of forces:" <<
nl 487 <<
" pressure : " <<
sum(force_[0]) <<
nl 488 <<
" viscous : " <<
sum(force_[1]) <<
nl 489 <<
" porous : " <<
sum(force_[2]) <<
nl 490 <<
" sum of moments:" <<
nl 491 <<
" pressure : " <<
sum(moment_[0]) <<
nl 492 <<
" viscous : " <<
sum(moment_[1]) <<
nl 493 <<
" porous : " <<
sum(moment_[2])
496 writeTime(
file(fileID::mainFile));
501 <<
sum(force_[0]) <<
setw(1) <<
' ' 502 <<
sum(force_[1]) <<
setw(1) <<
' ' 503 <<
sum(force_[2]) <<
setw(3) <<
") (" 504 <<
sum(moment_[0]) <<
setw(1) <<
' ' 505 <<
sum(moment_[1]) <<
setw(1) <<
' ' 506 <<
sum(moment_[2]) <<
setw(1) <<
')';
511 <<
sum(force_[0]) <<
setw(1) <<
' ' 512 <<
sum(force_[1]) <<
setw(3) <<
") (" 513 <<
sum(moment_[0]) <<
setw(1) <<
' ' 514 <<
sum(moment_[1]) <<
setw(1) <<
')';
522 <<
sum(coordSys_.localVector(force_[0])) <<
setw(1) <<
' ' 523 <<
sum(coordSys_.localVector(force_[1])) <<
setw(1) <<
' ' 524 <<
sum(coordSys_.localVector(force_[2])) <<
setw(3) <<
") (" 525 <<
sum(coordSys_.localVector(moment_[0])) <<
setw(1) <<
' ' 526 <<
sum(coordSys_.localVector(moment_[1])) <<
setw(1) <<
' ' 527 <<
sum(coordSys_.localVector(moment_[2])) <<
setw(1) <<
')';
532 <<
sum(coordSys_.localVector(force_[0])) <<
setw(1) <<
' ' 533 <<
sum(coordSys_.localVector(force_[1])) <<
setw(3) <<
") (" 534 <<
sum(coordSys_.localVector(moment_[0])) <<
setw(1) <<
' ' 535 <<
sum(coordSys_.localVector(moment_[1])) <<
setw(1) <<
')';
557 f[0][i] += f[0][i-1];
558 f[1][i] += f[1][i-1];
559 f[2][i] += f[2][i-1];
561 m[0][i] += m[0][i-1];
562 m[1][i] += m[1][i-1];
563 m[2][i] += m[2][i-1];
567 writeTime(
file(fileID::binsFile));
574 file(fileID::binsFile)
576 << f[0][i] <<
setw(1) <<
' ' 577 << f[1][i] <<
setw(1) <<
' ' 578 << f[2][i] <<
setw(3) <<
") (" 579 << m[0][i] <<
setw(1) <<
' ' 580 << m[1][i] <<
setw(1) <<
' ' 581 << m[2][i] <<
setw(1) <<
')';
585 file(fileID::binsFile)
587 << f[0][i] <<
setw(1) <<
' ' 588 << f[1][i] <<
setw(3) <<
") (" 589 << m[0][i] <<
setw(1) <<
' ' 590 << m[1][i] <<
setw(1) <<
')';
598 lf[0] = coordSys_.localVector(force_[0]);
599 lf[1] = coordSys_.localVector(force_[1]);
600 lf[2] = coordSys_.localVector(force_[2]);
601 lm[0] = coordSys_.localVector(moment_[0]);
602 lm[1] = coordSys_.localVector(moment_[1]);
603 lm[2] = coordSys_.localVector(moment_[2]);
609 lf[0][i] += lf[0][i-1];
610 lf[1][i] += lf[1][i-1];
611 lf[2][i] += lf[2][i-1];
612 lm[0][i] += lm[0][i-1];
613 lm[1][i] += lm[1][i-1];
614 lm[2][i] += lm[2][i-1];
622 file(fileID::binsFile)
624 << lf[0][i] <<
setw(1) <<
' ' 625 << lf[1][i] <<
setw(1) <<
' ' 626 << lf[2][i] <<
setw(3) <<
") (" 627 << lm[0][i] <<
setw(1) <<
' ' 628 << lm[1][i] <<
setw(1) <<
' ' 629 << lm[2][i] <<
setw(1) <<
')';
633 file(fileID::binsFile)
635 << lf[0][i] <<
setw(1) <<
' ' 636 << lf[1][i] <<
setw(3) <<
") (" 637 << lm[0][i] <<
setw(1) <<
' ' 638 << lm[1][i] <<
setw(1) <<
')';
665 directForceDensity_(
false),
677 binCumulative_(
true),
700 directForceDensity_(
false),
712 binCumulative_(
true),
731 initialised_ =
false;
735 directForceDensity_ = dict.
lookupOrDefault(
"directForceDensity",
false);
741 if (directForceDensity_)
772 if (rhoName_ ==
"rhoInf")
774 dict.
lookup(
"rhoInf") >> rhoRef_;
783 if (dict.
found(
"CofR"))
786 localSystem_ =
false;
797 Log <<
" Including porosity effects" <<
endl;
801 Log <<
" Not including porosity effects" <<
endl;
804 if (dict.
found(
"binData"))
807 binDict.
lookup(
"nBin") >> nBin_;
812 <<
"Number of bins (nBin) must be zero or greater" 815 else if ((nBin_ == 0) || (nBin_ == 1))
820 force_[i].setSize(1);
821 moment_[i].setSize(1);
827 binDict.lookup(
"direction") >> binDir_;
828 binDir_ /=
mag(binDir_);
831 scalar binMax = -great;
834 const label patchi = iter.key();
837 binMin_ =
min(
min(d), binMin_);
838 binMax =
max(
max(d), binMax);
845 binMax = 1.0001*(binMax - binMin_) + binMin_;
847 binDx_ = (binMax - binMin_)/scalar(nBin_);
850 binPoints_.setSize(nBin_);
853 binPoints_[i] = (i + 0.5)*binDir_*binDx_;
856 binDict.lookup(
"cumulative") >> binCumulative_;
861 force_[i].setSize(nBin_);
862 moment_[i].setSize(nBin_);
870 force_[0].setSize(1);
871 force_[1].setSize(1);
872 force_[2].setSize(1);
873 moment_[0].setSize(1);
874 moment_[1].setSize(1);
875 moment_[2].setSize(1);
878 resetNames(createFileNames(dict));
896 if (directForceDensity_)
901 mesh_.Sf().boundaryField();
905 const label patchi = iter.key();
909 mesh_.C().boundaryField()[
patchi] - coordSys_.origin()
929 applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[
patchi]);
937 mesh_.Sf().boundaryField();
941 tdevTau().boundaryField();
944 const scalar
pRef = pRef_/
rho(p);
948 const label patchi = iter.key();
952 mesh_.C().boundaryField()[
patchi] - coordSys_.origin()
963 const vectorField fT(Sfb[patchi] & devTaub[patchi]);
967 applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[
patchi]);
983 <<
"Porosity effects requested, but no porosity models found " 998 const label zoneI = cellZoneIDs[i];
999 const cellZone& cZone = mesh_.cellZones()[zoneI];
1007 applyBins(Md, fDummy, fDummy, fP, d);
1021 return sum(force_[0]) +
sum(force_[1]) +
sum(force_[2]);
1027 return sum(moment_[0]) +
sum(moment_[1]) +
sum(moment_[2]);
virtual vector momentEff() const
Return the total moment.
Base class for other coordinate system specifications.
const dimensionSet dimViscosity
virtual void writeFileHeader(const label i)
Output file header information.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Templated abstract base class for multiphase compressible turbulence models.
virtual bool write()
Write function.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
errorManipArg< error, int > exit(error &err, const int errNo=1)
A list of keyword definitions, which are a keyword followed by any number of values (e...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
forces(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An abstract base class for physical properties.
const dimensionSet dimPressure
void size(const label)
Override size to be inconsistent with allocated storage.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
bool empty() const
Return true if the hash table is empty.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static bool master(const label communicator=0)
Am I the master process.
virtual vector forceEff() const
Return the total force.
Vector< scalar > vector
A scalar version of the templated Vector.
const dimensionSet dimless
Templated abstract base class for multiphase incompressible turbulence models.
const labelList & cellZoneIDs() const
Return const access to the cell zone IDs.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
virtual bool read(const dictionary &)
Read the forces data.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Macros for easy insertion into run-time selection tables.
void writeBins()
Helper function to write bin data.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu) const
Return the force over the cell zone(s)
virtual bool execute()
Execute, currently does nothing.
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
Calculate the gradient of the given field.
A class for handling words, derived from string.
virtual tmp< volSymmTensorField > devSigma() const
Return the effective stress tensor.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< volScalarField > mu() const
Dynamic viscosity field.
static word groupName(Name name, const word &group)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
const dimensionSet dimDensity
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
static const word null
An empty word.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
virtual void calcForcesMoment()
Calculate the forces and moments.
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
void applyBins(const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP, const vectorField &d)
Accumulate bin data.
An STL-conforming hash table.
errorManip< error > abort(error &err)
const dimensionedScalar mu
Atomic mass unit.
void initialise()
Initialise the fields.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
defineTypeNameAndDebug(Qdot, 0)
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
void writeForces()
Helper function to write force data.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< volSymmTensorField > devTau() const
Return the effective viscous stress (laminar + turbulent).
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
tmp< volScalarField > alpha() const
Get the volume fraction field.
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Omanip< int > setw(const int i)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
A class for managing temporary objects.
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
virtual ~forces()
Destructor.
Base class for single-phase incompressible turbulence models.
fileID
Enumeration for ensuring the right file is accessed.
static const Vector< scalar > zero
Base class for single-phase compressible turbulence models.
Top level model for porosity models.
virtual bool write()
Write the forces.
wordList createFileNames(const dictionary &dict) const
Create file names for forces and bins.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
const vectorField::subField faceCentres() const
Return face centres.