forceCoeffs.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "forceCoeffs.H"
30 #include "dictionary.H"
31 #include "Time.H"
32 #include "Pstream.H"
33 #include "IOmanip.H"
34 #include "fvMesh.H"
35 #include "dimensionedTypes.H"
36 #include "volFields.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace functionObjects
44 {
45  defineTypeNameAndDebug(forceCoeffs, 0);
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
52 
54 {
55  if (initialised_)
56  {
57  return;
58  }
59 
60  initialised_ = true;
61 }
62 
63 
65 {
66  auto* coeffPtr =
67  mesh_.getObjectPtr<volVectorField>(scopedName("forceCoeff"));
68 
69  if (!coeffPtr)
70  {
71  coeffPtr = new volVectorField
72  (
73  IOobject
74  (
75  scopedName("forceCoeff"),
76  time_.timeName(),
77  mesh_,
80  ),
81  mesh_,
83  );
84 
85  mesh_.objectRegistry::store(coeffPtr);
86  }
87 
88  return *coeffPtr;
89 }
90 
91 
93 {
94  auto* coeffPtr =
95  mesh_.getObjectPtr<volVectorField>(scopedName("momentCoeff"));
96 
97  if (!coeffPtr)
98  {
99  coeffPtr = new volVectorField
100  (
101  IOobject
102  (
103  scopedName("momentCoeff"),
104  time_.timeName(),
105  mesh_,
108  ),
109  mesh_,
111  );
112 
113  mesh_.objectRegistry::store(coeffPtr);
114  }
115 
116  return *coeffPtr;
117 }
118 
119 
121 {
122  Cf_.reset();
123  Cm_.reset();
125  forceCoeff() == dimensionedVector(dimless, Zero);
126  momentCoeff() == dimensionedVector(dimless, Zero);
127 }
128 
129 
132 {
133  HashTable<coeffDesc> coeffs(16);
134 
135  coeffs.insert("Cd", coeffDesc("Drag force", "Cd", 0, 0));
136  coeffs.insert("Cs", coeffDesc("Side force", "Cs", 1, 2));
137  coeffs.insert("Cl", coeffDesc("Lift force", "Cl", 2, 1));
138 
139  // Add front/rear options
140  const auto frontRearCoeffs(coeffs);
141  forAllConstIters(frontRearCoeffs, iter)
142  {
143  const auto& fr = iter.val();
144  coeffs.insert(fr.frontName(), fr.front());
145  coeffs.insert(fr.rearName(), fr.rear());
146  }
147 
148  // Add moments
149  coeffs.insert("CmRoll", coeffDesc("Roll moment", "CmRoll", 0, -1));
150  coeffs.insert("CmPitch", coeffDesc("Pitch moment", "CmPitch", 1, -1));
151  coeffs.insert("CmYaw", coeffDesc("Yaw moment", "CmYaw", 2, -1));
152 
153  return coeffs;
154 }
155 
156 
158 {
159  // Calculate scaling factors
160  const scalar pDyn = 0.5*rhoRef_*sqr(magUInf_);
161  const dimensionedScalar forceScaling
162  (
164  scalar(1)/(Aref_*pDyn + SMALL)
165  );
166 
167  const auto& coordSys = coordSysPtr_();
168 
169  // Calculate force coefficients
170  forceCoeff() = forceScaling*force();
171 
172  Cf_.reset
173  (
174  forceScaling.value()*coordSys.localVector(sumPatchForcesP_),
175  forceScaling.value()*coordSys.localVector(sumPatchForcesV_),
176  forceScaling.value()*coordSys.localVector(sumInternalForces_)
177  );
178 }
179 
180 
182 {
183  // Calculate scaling factors
184  const scalar pDyn = 0.5*rhoRef_*sqr(magUInf_);
185  const dimensionedScalar momentScaling
186  (
188  scalar(1)/(Aref_*pDyn*lRef_ + SMALL)
189  );
190 
191  const auto& coordSys = coordSysPtr_();
192 
193  // Calculate moment coefficients
194  momentCoeff() = momentScaling*moment();
195 
196  Cm_.reset
197  (
198  momentScaling.value()*coordSys.localVector(sumPatchMomentsP_),
199  momentScaling.value()*coordSys.localVector(sumPatchMomentsV_),
200  momentScaling.value()*coordSys.localVector(sumInternalMoments_)
201  );
202 }
203 
204 
206 {
207  if (!coeffFilePtr_.valid())
208  {
209  coeffFilePtr_ = newFileAtStartTime("coefficient");
210  writeIntegratedDataFileHeader("Coefficients", coeffFilePtr_());
211  }
212 }
213 
214 
216 (
217  const word& header,
218  OFstream& os
219 ) const
220 {
221  const auto& coordSys = coordSysPtr_();
222 
223  writeHeader(os, "Force and moment coefficients");
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());
234  writeHeader(os, "");
235  writeCommented(os, "Time");
236 
237  for (const auto& iter : coeffs_.sorted())
238  {
239  const auto& coeff = iter.val();
240 
241  if (!coeff.active_) continue;
242 
243  writeTabbed(os, coeff.name_);
244  }
245 
246  os << endl;
247 }
248 
249 
251 {
252  OFstream& os = coeffFilePtr_();
253 
254  writeCurrentTime(os);
255 
256  for (const auto& iter : coeffs_.sorted())
257  {
258  const auto& coeff = iter.val();
259 
260  if (!coeff.active_) continue;
261 
262  const vector coeffValue = coeff.value(Cf_, Cm_);
263 
264  os << tab << (coeffValue.x() + coeffValue.y() + coeffValue.z());
265  }
266 
267  coeffFilePtr_() << endl;
268 }
269 
270 
271 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
272 
274 (
275  const word& name,
276  const Time& runTime,
277  const dictionary& dict,
278  const bool readFields
279 )
280 :
281  forces(name, runTime, dict, false),
282  Cf_(),
283  Cm_(),
284  coeffFilePtr_(),
285  magUInf_(Zero),
286  lRef_(Zero),
287  Aref_(Zero),
288  initialised_(false)
289 {
290  if (readFields)
291  {
292  read(dict);
293  setCoordinateSystem(dict, "liftDir", "dragDir");
295  }
296 }
297 
298 
299 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
300 
302 {
303  if (!forces::read(dict))
304  {
305  return false;
306  }
307 
308  initialised_ = false;
309 
310  // Reference scales
311  dict.readEntry("magUInf", magUInf_);
312  dict.readEntry("lRef", lRef_);
313  dict.readEntry("Aref", Aref_);
314 
315  // If case is compressible we must read rhoInf
316  // (stored in rhoRef_) to calculate the reference dynamic pressure
317  // Note: for incompressible, rhoRef_ is already initialised to 1
318  if (rhoName_ != "rhoInf")
319  {
320  dict.readEntry("rhoInf", rhoRef_);
321  }
322 
323  Info<< " magUInf: " << magUInf_ << nl
324  << " lRef : " << lRef_ << nl
325  << " Aref : " << Aref_ << nl
326  << " rhoInf : " << rhoRef_ << endl;
327 
328  if (min(magUInf_, rhoRef_) < SMALL || min(lRef_, Aref_) < SMALL)
329  {
331  << "Non-zero values are required for reference scales"
332  << exit(FatalIOError);
333 
334  return false;
335  }
336 
337  if (!dict.found("coefficients"))
338  {
339  Info<< " Selecting all coefficients" << nl;
340 
341  coeffs_ = selectCoeffs();
342 
343  for (auto& iter : coeffs_.sorted())
344  {
345  auto& coeff = iter.val();
346  coeff.active_ = true;
347  Info<< " - " << coeff << nl;
348  }
349  }
350  else
351  {
352  const wordHashSet coeffs(dict.get<wordHashSet>("coefficients"));
353 
354  coeffs_ = selectCoeffs();
355 
356  Info<< " Selecting coefficients:" << nl;
357 
358  for (const word& key : coeffs)
359  {
360  auto coeffIter = coeffs_.find(key);
361 
362  if (!coeffIter.good())
363  {
365  << "Unknown coefficient type " << key
366  << " . Valid entries are : " << coeffs_.sortedToc()
367  << exit(FatalIOError);
368  }
369  auto& coeff = coeffIter.val();
370  coeff.active_ = true;
371  Info<< " - " << coeff << nl;
372  }
373  }
374 
376 
377  return true;
378 }
379 
380 
381 
383 {
385 
386  initialise();
387 
388  reset();
389 
390  Log << type() << " " << name() << " write:" << nl
391  << " " << "Coefficient"
392  << tab << "Total"
393  << tab << "Pressure"
394  << tab << "Viscous"
395  << tab << "Internal"
396  << nl;
397 
398  calcForceCoeffs();
399 
400  calcMomentCoeffs();
401 
402  auto logValues = [](const word& name, const vector& coeff, Ostream& os)
403  {
404  os << " " << name << ":"
405  << tab << coeff.x() + coeff.y() + coeff.z()
406  << tab << coeff.x()
407  << tab << coeff.y()
408  << tab << coeff.z()
409  << nl;
410  };
411 
412  // Always setting all results
413  for (const auto& iter : coeffs_.sorted())
414  {
415  const word& coeffName = iter.key();
416  const auto& coeff = iter.val();
417 
418  // Vectors for x:pressure, y:viscous, z:internal
419  const vector coeffValue = coeff.value(Cf_, Cm_);
420 
421  if (log && coeff.active_)
422  {
423  logValues(coeffName, coeffValue, Info);
424  }
425 
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());
430  }
432  Log << endl;
433 
434  return true;
435 }
436 
437 
439 {
440  if (writeToFile())
441  {
442  Log << " writing force and moment coefficient files." << endl;
443 
444  createIntegratedDataFile();
445 
446  writeIntegratedDataFile();
447  }
448 
449  if (writeFields_)
450  {
451  forceCoeff().write();
452  momentCoeff().write();
453  }
454 
455  Log << endl;
456 
457  return true;
458 }
459 
460 
461 // ************************************************************************* //
const Type & value() const noexcept
Return const reference to value.
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
void calcForceCoeffs()
Calculate force coefficients.
Definition: forceCoeffs.C:150
virtual bool write()
Write to data files/fields and to streams.
Definition: forceCoeffs.C:431
static void writeHeader(Ostream &os, const word &fieldName)
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Output to file stream, using an OSstream.
Definition: OFstream.H:49
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:48
void writeIntegratedDataFile()
Write integrated coefficients to the integrated-coefficient file.
Definition: forceCoeffs.C:243
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
volVectorField & momentCoeff()
Return access to the moment coefficients field.
Definition: forceCoeffs.C:85
Abstract base-class for Time/database function objects.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
HashTable< coeffDesc > selectCoeffs() const
Return the operand coefficients.
Definition: forceCoeffs.C:124
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
forceCoeffs(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from Time and dictionary.
Definition: forceCoeffs.C:267
virtual bool read(const dictionary &dict)
Read the dictionary.
Definition: forces.C:572
Macros for easy insertion into run-time selection tables.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:752
virtual void calcForcesMoments()
Calculate forces and moments.
Definition: forces.C:649
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
A class for handling words, derived from Foam::string.
Definition: word.H:63
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimPressure, Zero))
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
void initialise()
Initialise containers and fields.
Definition: forceCoeffs.C:46
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
void reset()
Reset containers and fields.
Definition: forceCoeffs.C:113
void calcMomentCoeffs()
Calculate moment coefficients.
Definition: forceCoeffs.C:174
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
virtual bool read(const dictionary &dict)
Read the dictionary.
Definition: forceCoeffs.C:294
Istream and Ostream manipulators taking arguments.
const dimensionSet dimForce
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:151
OBJstream os(runTime.globalPath()/outputName)
volVectorField & forceCoeff()
Return access to the force coefficients field.
Definition: forceCoeffs.C:57
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
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.
Definition: forces.C:45
Computes force and moment coefficients over a given list of patches, and optionally over given porous...
Definition: forceCoeffs.H:433
void writeIntegratedDataFileHeader(const word &header, OFstream &os) const
Write header to the integrated-coefficient file.
Definition: forceCoeffs.C:209
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:607
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
Nothing to be read.
#define Log
Definition: PDRblock.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition: forces.H:317
void createIntegratedDataFile()
Create the integrated-coefficient file.
Definition: forceCoeffs.C:198
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
virtual bool execute()
Execute the function object.
Definition: forceCoeffs.C:375
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157