41namespace functionObjects
75 if (
dict.
found(coordinateSystem::typeName_()))
83 coordinateSystem::typeName_()
96 auto* forcePtr = mesh_.getObjectPtr<
volVectorField>(scopedName(
"force"));
114 mesh_.objectRegistry::store(forcePtr);
123 auto* momentPtr = mesh_.getObjectPtr<
volVectorField>(scopedName(
"moment"));
131 scopedName(
"moment"),
141 mesh_.objectRegistry::store(momentPtr);
155 if (directForceDensity_)
157 if (!foundObject<volVectorField>(fDName_))
160 <<
"Could not find " << fDName_ <<
" in database"
168 !foundObject<volVectorField>(UName_)
169 || !foundObject<volScalarField>(pName_)
173 <<
"Could not find U: " << UName_
174 <<
" or p:" << pName_ <<
" in database"
178 if (rhoName_ !=
"rhoInf" && !foundObject<volScalarField>(rhoName_))
181 <<
"Could not find rho:" << rhoName_ <<
" in database"
192 sumPatchForcesP_ =
Zero;
193 sumPatchForcesV_ =
Zero;
194 sumPatchMomentsP_ =
Zero;
195 sumPatchMomentsV_ =
Zero;
197 sumInternalForces_ =
Zero;
198 sumInternalMoments_ =
Zero;
200 auto& force = this->force();
201 auto& moment = this->moment();
213 if (foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
216 lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
218 return turb.devRhoReff();
220 else if (foundObject<icoTurbModel>(icoTurbModel::propertiesName))
223 lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
231 const auto&
U = lookupObject<volVectorField>(UName_);
235 else if (foundObject<transportModel>(
"transportProperties"))
237 const auto& laminarT =
238 lookupObject<transportModel>(
"transportProperties");
240 const auto&
U = lookupObject<volVectorField>(UName_);
244 else if (foundObject<dictionary>(
"transportProperties"))
247 lookupObject<dictionary>(
"transportProperties");
251 const auto&
U = lookupObject<volVectorField>(UName_);
258 <<
"No valid model for viscous stress calculation"
274 else if (foundObject<transportModel>(
"transportProperties"))
276 const auto& laminarT =
277 lookupObject<transportModel>(
"transportProperties");
279 return rho()*laminarT.nu();
281 else if (foundObject<dictionary>(
"transportProperties"))
284 lookupObject<dictionary>(
"transportProperties");
293 <<
"No valid model for dynamic viscosity calculation"
303 if (rhoName_ ==
"rhoInf")
310 mesh_.time().timeName(),
318 return (lookupObject<volScalarField>(rhoName_));
329 if (rhoName_ !=
"rhoInf")
332 <<
"Dynamic pressure is expected but kinematic is provided."
348 sumPatchForcesP_ +=
sum(fP);
349 sumPatchForcesV_ +=
sum(fV);
350 force().boundaryFieldRef()[patchi] += fP + fV;
355 sumPatchMomentsP_ +=
sum(mP);
356 sumPatchMomentsV_ +=
sum(mV);
357 moment().boundaryFieldRef()[patchi] += mP + mV;
368 auto& force = this->force();
369 auto& moment = this->moment();
373 const label celli =
cellIDs[i];
375 sumInternalForces_ +=
f[i];
376 force[celli] +=
f[i];
379 sumInternalMoments_ += m;
387 if (!forceFilePtr_.valid())
389 forceFilePtr_ = createFile(
"force");
390 writeIntegratedDataFileHeader(
"Force", forceFilePtr_());
393 if (!momentFilePtr_.valid())
395 momentFilePtr_ = createFile(
"moment");
396 writeIntegratedDataFileHeader(
"Moment", momentFilePtr_());
407 const auto& coordSys = coordSysPtr_();
408 const auto vecDesc = [](
const word& root)->
string
410 return root +
"_x " + root +
"_y " + root +
"_z";
413 writeHeaderValue(
os,
"CofR", coordSys.origin());
415 writeCommented(
os,
"Time");
416 writeTabbed(
os, vecDesc(
"total"));
417 writeTabbed(
os, vecDesc(
"pressure"));
418 writeTabbed(
os, vecDesc(
"viscous"));
422 writeTabbed(
os, vecDesc(
"porous"));
431 const auto& coordSys = coordSysPtr_();
433 writeIntegratedDataFile
435 coordSys.localVector(sumPatchForcesP_),
436 coordSys.localVector(sumPatchForcesV_),
437 coordSys.localVector(sumInternalForces_),
441 writeIntegratedDataFile
443 coordSys.localVector(sumPatchMomentsP_),
444 coordSys.localVector(sumPatchMomentsV_),
445 coordSys.localVector(sumInternalMoments_),
459 writeCurrentTime(
os);
461 writeValue(
os, pres + vis + internal);
462 writeValue(
os, pres);
467 writeValue(
os, internal);
476 const string& descriptor,
487 Log <<
" Sum of " << descriptor.c_str() <<
nl
488 <<
" Total : " << (pres + vis + internal) <<
nl
489 <<
" Pressure : " << pres <<
nl
490 <<
" Viscous : " << vis <<
nl;
494 Log <<
" Porous : " << internal <<
nl;
511 sumPatchForcesP_(
Zero),
512 sumPatchForcesV_(
Zero),
513 sumPatchMomentsP_(
Zero),
514 sumPatchMomentsV_(
Zero),
515 sumInternalForces_(
Zero),
516 sumInternalMoments_(
Zero),
519 coordSysPtr_(nullptr),
527 directForceDensity_(false),
551 sumPatchForcesP_(
Zero),
552 sumPatchForcesV_(
Zero),
553 sumPatchMomentsP_(
Zero),
554 sumPatchMomentsV_(
Zero),
555 sumInternalForces_(
Zero),
556 sumInternalMoments_(
Zero),
559 coordSysPtr_(nullptr),
567 directForceDensity_(false),
590 initialised_ =
false;
595 mesh_.boundaryMesh().patchSet
600 dict.readIfPresent(
"directForceDensity", directForceDensity_);
601 if (directForceDensity_)
604 if (
dict.readIfPresent<
word>(
"fD", fDName_))
612 if (
dict.readIfPresent<
word>(
"p", pName_))
616 if (
dict.readIfPresent<
word>(
"U", UName_))
620 if (
dict.readIfPresent<
word>(
"rho", rhoName_))
622 Info<<
" rho: " << rhoName_ <<
endl;
626 if (rhoName_ ==
"rhoInf")
629 Info<<
" Freestream density (rhoInf) set to " << rhoRef_ <<
endl;
633 if (
dict.readIfPresent<scalar>(
"pRef", pRef_))
635 Info<<
" Reference pressure (pRef) set to " << pRef_ <<
endl;
639 dict.readIfPresent(
"porosity", porosity_);
642 Info<<
" Including porosity effects" <<
endl;
646 Info<<
" Not including porosity effects" <<
endl;
649 writeFields_ =
dict.getOrDefault(
"writeFields",
false);
652 Info<<
" Fields will be written" <<
endl;
666 const point& origin = coordSysPtr_->origin();
668 if (directForceDensity_)
670 const auto& fD = lookupObject<volVectorField>(fDName_);
672 const auto& Sfb = mesh_.Sf().boundaryField();
674 for (
const label patchi : patchSet_)
676 const vectorField& d = mesh_.C().boundaryField()[patchi];
687 Sfb[patchi] & fD.boundaryField()[patchi]
692 const vectorField fV(sA*fD.boundaryField()[patchi] - fP);
694 addToPatchFields(patchi, Md, fP, fV);
699 const auto&
p = lookupObject<volScalarField>(pName_);
701 const auto& Sfb = mesh_.Sf().boundaryField();
704 const auto& devRhoReffb = tdevRhoReff().boundaryField();
707 const scalar rhoRef =
rho(
p);
708 const scalar pRef = pRef_/rhoRef;
710 for (
const label patchi : patchSet_)
712 const vectorField& d = mesh_.C().boundaryField()[patchi];
718 rhoRef*Sfb[patchi]*(
p.boundaryField()[patchi] - pRef)
721 const vectorField fV(Sfb[patchi] & devRhoReffb[patchi]);
723 addToPatchFields(patchi, Md, fP, fV);
729 const auto&
U = lookupObject<volVectorField>(UName_);
738 <<
"Porosity effects requested, "
739 <<
"but no porosity models found in the database"
750 const labelList& cellZoneIDs = pm.cellZoneIDs();
752 for (
const label zonei : cellZoneIDs)
754 const cellZone& cZone = mesh_.cellZones()[zonei];
760 addToInternalField(cZone, Md, fP);
776 return sumPatchForcesP_ + sumPatchForcesV_ + sumInternalForces_;
782 return sumPatchMomentsP_ + sumPatchMomentsV_ + sumInternalMoments_;
792 const auto& coordSys = coordSysPtr_();
794 const auto localFp(coordSys.localVector(sumPatchForcesP_));
795 const auto localFv(coordSys.localVector(sumPatchForcesV_));
796 const auto localFi(coordSys.localVector(sumInternalForces_));
798 logIntegratedData(
"forces", localFp, localFv, localFi);
800 const auto localMp(coordSys.localVector(sumPatchMomentsP_));
801 const auto localMv(coordSys.localVector(sumPatchMomentsV_));
802 const auto localMi(coordSys.localVector(sumInternalMoments_));
804 logIntegratedData(
"moments", localMp, localMv, localMi);
806 setResult(
"pressureForce", localFp);
807 setResult(
"viscousForce", localFv);
808 setResult(
"internalForce", localFi);
809 setResult(
"pressureMoment", localMp);
810 setResult(
"viscousMoment", localMv);
811 setResult(
"internalMoment", localMi);
821 Log <<
" writing force and moment files." <<
endl;
823 createIntegratedDataFiles();
824 writeIntegratedDataFiles();
829 Log <<
" writing force and moment fields." <<
endl;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
compressible::turbulenceModel & turb
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Templated abstract base class for single-phase incompressible turbulence models.
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Output to file stream, using an OSstream.
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read()
Re-read model coefficients if they have changed.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
static const word dictName
A Cartesian coordinate system.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Abstract base-class for Time/database function objects.
void writeIntegratedDataFile()
Write integrated coefficients to the integrated-coefficient file.
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
virtual void calcForcesMoments()
Calculate forces and moments.
void writeIntegratedDataFileHeader(const word &header, OFstream &os) const
Write header for an integrated-data file.
void initialise()
Initialise containers and fields.
void addToPatchFields(const label patchi, const vectorField &Md, const vectorField &fP, const vectorField &fV)
Add patch contributions to force and moment fields.
void writeIntegratedDataFiles()
Write integrated data to files.
volVectorField & moment()
Return access to the moment field.
void createIntegratedDataFiles()
Create the integrated-data files.
tmp< volScalarField > mu() const
Return dynamic viscosity field.
virtual bool read(const dictionary &dict)
Read the dictionary.
virtual vector forceEff() const
Return the total force.
void logIntegratedData(const string &descriptor, const vector &pres, const vector &vis, const vector &internal) const
Write integrated data to stream.
autoPtr< coordinateSystem > coordSysPtr_
Coordinate system used when evaluating forces and moments.
virtual vector momentEff() const
Return the total moment.
void addToInternalField(const labelList &cellIDs, const vectorField &Md, const vectorField &f)
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
void reset()
Reset containers and fields.
virtual bool execute()
Execute the function object.
volVectorField & force()
Return access to the force field.
virtual bool write()
Write to data files/fields and to streams.
tmp< volSymmTensorField > devRhoReff() const
Return the effective stress (viscous + turbulent)
void setCoordinateSystem(const dictionary &dict, const word &e3Name=word::null, const word &e1Name=word::null)
Set the co-ordinate system from dictionary and axes names.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the natural logarithm of an input volScalarField.
Computes the magnitude of an input field.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
const objectRegistry & obr_
Reference to the region objectRegistry.
Base class for writing single files from the function objects.
Registry of regIOobjects.
Top level model for porosity models.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
A List of wordRe with additional matching capabilities.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
const volScalarField & mu
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
Calculate the gradient of the given field.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
const dimensionSet dimPressure
const dimensionSet dimViscosity
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static void writeHeader(Ostream &os, const word &fieldName)
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Ostream & endl(Ostream &os)
Add newline and flush stream.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimForce
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a)
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.