43namespace functionObjects
75 scopedName(
"forceCoeff"),
85 mesh_.objectRegistry::store(coeffPtr);
103 scopedName(
"momentCoeff"),
113 mesh_.objectRegistry::store(coeffPtr);
140 const auto frontRearCoeffs(coeffs);
143 const auto& fr = iter.val();
144 coeffs.
insert(fr.frontName(), fr.front());
145 coeffs.
insert(fr.rearName(), fr.rear());
160 const scalar
pDyn = 0.5*rhoRef_*
sqr(magUInf_);
164 scalar(1)/(Aref_*
pDyn + SMALL)
167 const auto& coordSys = coordSysPtr_();
170 forceCoeff() = forceScaling*force();
174 forceScaling.
value()*coordSys.localVector(sumPatchForcesV_),
175 forceScaling.
value()*coordSys.localVector(sumPatchForcesP_),
176 forceScaling.
value()*coordSys.localVector(sumInternalForces_)
184 const scalar
pDyn = 0.5*rhoRef_*
sqr(magUInf_);
188 scalar(1)/(Aref_*
pDyn*lRef_ + SMALL)
191 const auto& coordSys = coordSysPtr_();
194 momentCoeff() = momentScaling*moment();
198 momentScaling.
value()*coordSys.localVector(sumPatchMomentsP_),
199 momentScaling.
value()*coordSys.localVector(sumPatchMomentsV_),
200 momentScaling.
value()*coordSys.localVector(sumInternalMoments_)
207 if (!coeffFilePtr_.valid())
209 coeffFilePtr_ = createFile(
"coefficient");
210 writeIntegratedDataFileHeader(
"Coefficients", coeffFilePtr_());
221 const auto& coordSys = coordSysPtr_();
224 writeHeaderValue(
os,
"dragDir", coordSys.e1());
225 writeHeaderValue(
os,
"sideDir", coordSys.e2());
226 writeHeaderValue(
os,
"liftDir", coordSys.e3());
227 writeHeaderValue(
os,
"rollAxis", coordSys.e1());
228 writeHeaderValue(
os,
"pitchAxis", coordSys.e2());
229 writeHeaderValue(
os,
"yawAxis", coordSys.e3());
230 writeHeaderValue(
os,
"magUInf", magUInf_);
231 writeHeaderValue(
os,
"lRef", lRef_);
232 writeHeaderValue(
os,
"Aref", Aref_);
233 writeHeaderValue(
os,
"CofR", coordSys.origin());
235 writeCommented(
os,
"Time");
237 for (
const auto& iter : coeffs_.sorted())
239 const auto& coeff = iter.val();
241 if (!coeff.active_)
continue;
243 writeTabbed(
os, coeff.name_);
254 writeCurrentTime(
os);
256 for (
const auto& iter : coeffs_.sorted())
258 const auto& coeff = iter.val();
260 if (!coeff.active_)
continue;
262 const vector coeffValue = coeff.value(Cf_, Cm_);
264 os <<
tab << (coeffValue.
x() + coeffValue.
y() + coeffValue.
z());
267 coeffFilePtr_() <<
endl;
308 initialised_ =
false;
311 dict.readEntry(
"magUInf", magUInf_);
312 dict.readEntry(
"lRef", lRef_);
313 dict.readEntry(
"Aref", Aref_);
318 if (rhoName_ !=
"rhoInf")
320 dict.readEntry(
"rhoInf", rhoRef_);
323 Info<<
" magUInf: " << magUInf_ <<
nl
324 <<
" lRef : " << lRef_ <<
nl
325 <<
" Aref : " << Aref_ <<
nl
326 <<
" rhoInf : " << rhoRef_ <<
endl;
328 if (
min(magUInf_, rhoRef_) < SMALL ||
min(lRef_, Aref_) < SMALL)
331 <<
"Non-zero values are required for reference scales"
337 if (!
dict.found(
"coefficients"))
339 Info<<
" Selecting all coefficients" <<
nl;
341 coeffs_ = selectCoeffs();
343 for (
auto& iter : coeffs_.sorted())
345 auto& coeff = iter.val();
346 coeff.active_ =
true;
347 Info<<
" - " << coeff <<
nl;
354 coeffs_ = selectCoeffs();
356 Info<<
" Selecting coefficients:" <<
nl;
358 for (
const word& key : coeffs)
360 auto coeffIter = coeffs_.find(key);
362 if (!coeffIter.good())
365 <<
"Unknown coefficient type " << key
366 <<
" . Valid entries are : " << coeffs_.sortedToc()
369 auto& coeff = coeffIter.val();
370 coeff.active_ =
true;
371 Info<<
" - " << coeff <<
nl;
391 <<
" " <<
"Coefficient"
405 <<
tab << coeff.x() + coeff.y() + coeff.z()
413 for (
const auto& iter : coeffs_.sorted())
415 const word& coeffName = iter.key();
416 const auto& coeff = iter.val();
419 const vector coeffValue = coeff.value(Cf_, Cm_);
421 if (
log && coeff.active_)
423 logValues(coeffName, coeffValue,
Info);
426 setResult(coeffName, coeffValue.
x() + coeffValue.
y() + coeffValue.
z());
427 setResult(coeffName &
"pressure", coeffValue.
x());
428 setResult(coeffName &
"viscous", coeffValue.
y());
429 setResult(coeffName &
"internal", coeffValue.
z());
442 Log <<
" writing force and moment coefficient files." <<
endl;
444 createIntegratedDataFile();
446 writeIntegratedDataFile();
451 forceCoeff().write();
452 momentCoeff().write();
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.
A HashTable similar to std::unordered_map.
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Output to file stream, using an OSstream.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
const Cmpt & z() const
Access to the vector z component.
const Cmpt & y() const
Access to the vector y component.
const Cmpt & x() const
Access to the vector x component.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const Type & value() const
Return const reference to value.
Abstract base-class for Time/database function objects.
Computes force and moment coefficients over a given list of patches, and optionally over given porous...
volVectorField & momentCoeff()
Return access to the moment coefficients field.
void calcMomentCoeffs()
Calculate moment coefficients.
HashTable< coeffDesc > selectCoeffs() const
Return the operand coefficients.
void writeIntegratedDataFileHeader(const word &header, OFstream &os) const
Write header to the integrated-coefficient file.
void initialise()
Initialise containers and fields.
void calcForceCoeffs()
Calculate force coefficients.
virtual bool read(const dictionary &dict)
Read the dictionary.
void reset()
Reset containers and fields.
void createIntegratedDataFile()
Create the integrated-coefficient file.
virtual bool execute()
Execute the function object.
virtual bool write()
Write to data files/fields and to streams.
void writeIntegratedDataFile()
Write integrated coefficients to the integrated-coefficient file.
volVectorField & forceCoeff()
Return access to the force coefficients field.
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 setCoordinateSystem(const dictionary &dict, const word &e3Name=word::null, const word &e1Name=word::null)
Set the co-ordinate system from dictionary and axes names.
Computes the natural logarithm of an input volScalarField.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
OBJstream os(runTime.globalPath()/outputName)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
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.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a)
constexpr char tab
The tab '\t' character(0x09)
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimPressure, Zero))
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.