50void Foam::functionObjects::momentum::purgeFields()
58template<
class GeoField>
60Foam::functionObjects::momentum::newField
86void Foam::functionObjects::momentum::calc()
99 autoPtr<volVectorField> tmomentum, tAngularMom, tAngularVel;
103 const auto&
U = lookupObject<volVectorField>(UName_);
104 const auto*
rhoPtr = findObject<volScalarField>(rhoName_);
122 auto* pmomentum = getObjectPtr<volVectorField>(scopedName(
"momentum"));
127 pmomentum = tmomentum.get();
137 momentum.ref() = (
U * mesh_.V() * rhoRef);
139 momentum.correctBoundaryConditions();
146 getObjectPtr<volVectorField>(scopedName(
"angularMomentum"));
148 if (hasCsys_ && !pAngularMom)
152 pAngularMom = tAngularMom.get();
154 else if (!pAngularMom)
158 pAngularMom = pmomentum;
160 auto& angularMom = *pAngularMom;
167 getObjectPtr<volVectorField>(scopedName(
"angularVelocity"));
174 newField<volVectorField>(
"angularVelocity",
dimVelocity);
175 pAngularVel = tAngularVel.get();
177 auto& angularVel = *pAngularVel;
182 angularVel.primitiveFieldRef() =
183 csys_.invTransform(mesh_.cellCentres(),
U.internalField());
185 angularVel.correctBoundaryConditions();
189 angularMom.ref() = (angularVel * mesh_.V() * (*
rhoPtr));
193 angularMom.ref() = (angularVel * mesh_.V() * rhoRef);
196 angularMom.correctBoundaryConditions();
203 sumAngularMom_ =
Zero;
207 for (label celli=0; celli < mesh_.nCells(); ++celli)
210 sumAngularMom_ += angularMom[celli];
215 for (
const label celli : cellIDs())
218 sumAngularMom_ += angularMom[celli];
231 if (!writeToFile() || writtenHeader_)
239 writeHeaderValue(
os,
"origin", csys_.origin());
240 writeHeaderValue(
os,
"axis", csys_.e3());
252 "Selection " + regionTypeNames_[regionType_]
253 +
" = " + regionName_
258 writeCommented(
os,
"Time");
259 writeTabbed(
os,
"(momentum_x momentum_y momentum_z)");
263 writeTabbed(
os,
"(momentum_r momentum_rtheta momentum_axis)");
280 if (!foundObject<volVectorField>(UName_))
283 <<
"Could not find U: " << UName_ <<
" in database"
288 const auto* pPtr = cfindObject<volScalarField>(pName_);
294 if (!foundObject<volScalarField>(rhoName_))
297 <<
"Could not find rho:" << rhoName_
312 Info<<
" Sum of Momentum";
316 Info<<
' ' << regionTypeNames_[regionType_]
317 <<
' ' << regionName_;
321 <<
" linear : " << sumMomentum_ <<
nl;
325 Info<<
" angular : " << sumAngularMom_ <<
nl;
333 writeCurrentTime(
os);
335 os <<
tab << sumMomentum_;
339 os <<
tab << sumAngularMom_;
358 volRegion(fvMeshFunctionObject::mesh_,
dict),
361 sumAngularMom_(
Zero),
368 writeMomentum_(false),
369 writeVelocity_(false),
370 writePosition_(false),
389 fvMeshFunctionObject(
name, obr,
dict),
390 volRegion(fvMeshFunctionObject::mesh_,
dict),
393 sumAngularMom_(
Zero),
400 writeMomentum_(false),
401 writeVelocity_(false),
402 writePosition_(false),
426 UName_ =
dict.getOrDefault<
word>(
"U",
"U");
427 pName_ =
dict.getOrDefault<
word>(
"p",
"p");
428 rhoName_ =
dict.getOrDefault<
word>(
"rho",
"rho");
429 rhoRef_ =
dict.getOrDefault<scalar>(
"rhoRef", 1.0);
430 hasCsys_ =
dict.getOrDefault(
"cylindrical",
false);
437 writeMomentum_ =
dict.getOrDefault(
"writeMomentum",
false);
438 writeVelocity_ =
dict.getOrDefault(
"writeVelocity",
false);
439 writePosition_ =
dict.getOrDefault(
"writePosition",
false);
441 Info<<
"Integrating for selection: "
442 << regionTypeNames_[regionType_]
443 <<
" (" << regionName_ <<
")" <<
nl;
447 Info<<
" Momentum fields will be written" <<
endl;
467 Info<<
" Angular velocity will be written" <<
endl;
471 newField<volVectorField>(
"angularVelocity",
dimVelocity)
477 Info<<
" Angular position will be written" <<
endl;
491 writeFileHeader(file());
499 setResult(
"momentum_x", sumMomentum_[0]);
500 setResult(
"momentum_y", sumMomentum_[1]);
501 setResult(
"momentum_z", sumMomentum_[2]);
503 setResult(
"momentum_r", sumAngularMom_[0]);
513 if (writeMomentum_ || (hasCsys_ && (writeVelocity_ || writePosition_)))
515 Log <<
"Writing fields" <<
nl;
519 fieldPtr = findObject<volVectorField>(scopedName(
"momentum"));
520 if (fieldPtr) fieldPtr->
write();
522 fieldPtr = findObject<volVectorField>(scopedName(
"angularMomentum"));
523 if (fieldPtr) fieldPtr->
write();
525 fieldPtr = findObject<volVectorField>(scopedName(
"angularVelocity"));
526 if (fieldPtr) fieldPtr->
write();
528 if (hasCsys_ && writePosition_)
533 auto cyl_r = newField<volScalarField>(
"cyl_r",
dimLength);
534 auto cyl_t = newField<volScalarField>(
"cyl_theta",
dimless);
535 auto cyl_z = newField<volScalarField>(
"cyl_z",
dimLength);
539 const auto&
pts = mesh_.cellCentres();
540 const label len =
pts.size();
546 for (label i=0; i < len; ++i)
557 const polyBoundaryMesh&
pbm = mesh_.boundaryMesh();
566 const auto&
pts =
pbm[patchi].faceCentres();
567 const label len =
pts.size();
569 UList<scalar>& r = cyl_r->boundaryFieldRef(
false)[patchi];
570 UList<scalar>& t = cyl_t->boundaryFieldRef(
false)[patchi];
571 UList<scalar>& z = cyl_z->boundaryFieldRef(
false)[patchi];
573 for (label i=0; i < len; ++i)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
T * get() noexcept
Return pointer to managed object without nullptr checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the natural logarithm of an input volScalarField.
Computes linear/angular momentum, reporting integral values and optionally writing the fields.
word pName_
The pressure field name (optional).
scalar rhoRef_
Reference density (for incompressible).
void initialise()
Initialise the fields.
word UName_
The velocity field name (optional).
bool hasCsys_
Are we using the cylindrical coordinate system?
bool writePosition_
Write fields flag.
bool writeMomentum_
Write fields flag.
bool initialised_
Initialised flag.
void writeValues(Ostream &os)
Write momentum data.
momentum(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
vector sumMomentum_
Integral (linear) momentum.
word rhoName_
The density field name (optional).
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
vector sumAngularMom_
Integral angular momentum.
virtual void writeFileHeader(Ostream &os)
Output file header information.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
coordSystem::cylindrical csys_
Coordinate system for evaluating angular momentum.
bool writeVelocity_
Write fields flag.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
const ObjectType * cfindObject(const word &fieldName) const
Return const pointer to the object (eg, a field) in the (sub) objectRegistry.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
const objectRegistry & obr_
Reference to the region objectRegistry.
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
const ObjectType * findObject(const word &fieldName) const
Return const pointer to the object (eg, a field) in the (sub) objectRegistry.
void setResult(const word &entryName, const Type &value)
Add result.
const Time & time_
Reference to the time database.
bool useAllCells() const noexcept
Use all cells, not the cellIDs.
scalar V() const
Return total volume of the selected region.
volRegion(const fvMesh &mesh, const dictionary &dict)
Construct from fvMesh and dictionary.
wordRe regionName_
Region name (cellSet, cellZone, ...).
virtual bool read(const dictionary &dict)
Read from dictionary.
static const Enum< regionTypes > regionTypeNames_
Region type names.
bool update()
Update the cached values as required.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
regionTypes regionType_
Region type.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
virtual bool read(const dictionary &dict)
Read.
bool writtenHeader_
Flag to identify whether the header has been written.
virtual OFstream & file()
Return access to the file (if only 1).
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
virtual bool writeToFile() const
Flag to allow writing to file.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
bool checkOut(regIOobject *io) const
Remove a regIOobject from registry and free memory if the object is ownedByRegistry....
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
bool store()
Register object with its registry and transfer ownership to the registry.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
const dimensionSet dimPressure
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
static void writeHeader(Ostream &os, const word &fieldName)
const dimensionSet dimVelocity
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.
dimensionedScalar log(const dimensionedScalar &ds)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
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.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
constexpr char tab
The tab '\t' character(0x09).
static bool initialised_(false)
Info<< "Reading mechanical properties\n"<< endl;IOdictionary mechanicalProperties(IOobject("mechanicalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));const dictionary &rhoDict(mechanicalProperties.subDict("rho"));word rhoType(rhoDict.get< word >("type"));autoPtr< volScalarField > rhoPtr
#define forAll(list, i)
Loop across all elements in list.