Loading...
Searching...
No Matches
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-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41namespace Foam
42{
43namespace functionObjects
44{
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* ptr = mesh_.getObjectPtr<volVectorField>(scopedName("forceCoeff"));
67
68 if (!ptr)
69 {
70 ptr = new volVectorField
71 (
73 (
74 scopedName("forceCoeff"),
75 time_.timeName(),
76 mesh_.thisDb(),
80 ),
81 mesh_,
83 );
84
86 }
87
88 return *ptr;
89}
90
91
93{
94 auto* ptr = mesh_.getObjectPtr<volVectorField>(scopedName("momentCoeff"));
95
96 if (!ptr)
97 {
98 ptr = new volVectorField
99 (
101 (
102 scopedName("momentCoeff"),
103 time_.timeName(),
104 mesh_.thisDb(),
108 ),
109 mesh_,
111 );
112
114 }
115
116 return *ptr;
117}
118
119
121{
122 Cf_.reset();
123 Cm_.reset();
125 forceCoeff() == Zero;
126 momentCoeff() == 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_)
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_.csorted())
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_.csorted())
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 }
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 iter = coeffs_.find(key);
361
362 if (!iter.good())
363 {
365 << "Unknown coefficient type " << key
366 << " . Valid entries are : " << coeffs_.sortedToc()
367 << exit(FatalIOError);
368 }
369 auto& coeff = iter.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_.csorted())
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// ************************************************************************* //
Istream and Ostream manipulators taking arguments.
#define Log
Definition PDRblock.C:28
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.
Definition HashTable.H:124
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
@ REGISTER
Request registration (bool: true).
@ 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,...
Definition IOobject.H:191
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
const Cmpt & x() const noexcept
Access to the vector x component.
Definition Vector.H:135
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const Type & value() const noexcept
Return const reference to value.
Abstract base-class for Time/database function objects.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
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.
Definition forceCoeffs.C:85
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.
Definition forceCoeffs.C:46
void calcForceCoeffs()
Calculate force coefficients.
void reset()
Reset containers and fields.
void createIntegratedDataFile()
Create the integrated-coefficient file.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
forceCoeffs(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
virtual bool read(const dictionary &)
Read the function-object dictionary.
void writeIntegratedDataFile()
Write integrated coefficients to the integrated-coefficient file.
volVectorField & forceCoeff()
Return access to the force coefficients field.
Definition forceCoeffs.C:57
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition forces.H:321
virtual void calcForcesMoments()
Calculate forces and moments.
Definition forces.C:665
scalar rhoRef_
Reference density needed for incompressible calculations.
Definition forces.H:387
vector sumPatchForcesV_
Sum of patch viscous forces.
Definition forces.H:336
vector sumPatchForcesP_
Sum of patch pressure forces.
Definition forces.H:331
vector sumInternalMoments_
Sum of internal moments.
Definition forces.H:356
volVectorField & moment()
Return access to the moment field.
Definition forces.C:111
vector sumPatchMomentsV_
Sum of patch viscous moments.
Definition forces.H:346
vector sumInternalForces_
Sum of internal forces.
Definition forces.H:351
forces(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
Definition forces.C:510
autoPtr< coordinateSystem > coordSysPtr_
Coordinate system used when evaluating forces and moments.
Definition forces.H:377
vector sumPatchMomentsP_
Sum of patch pressure moments.
Definition forces.H:341
word rhoName_
Name of density field.
Definition forces.H:407
bool writeFields_
Flag to write force and moment fields.
Definition forces.H:427
volVectorField & force()
Return access to the force field.
Definition forces.C:83
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
virtual bool read(const dictionary &)
Read the function-object dictionary.
Definition forces.C:589
const fvMesh & mesh_
Reference to the fvMesh.
Computes the natural logarithm of an input volScalarField.
Definition log.H:212
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
void setResult(const word &entryName, const Type &value)
Add result.
const Time & time_
Reference to the time database.
virtual autoPtr< OFstream > newFileAtStartTime(const word &name) const
Return autoPtr to a new file using the simulation start time.
Definition writeFile.C:156
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
virtual bool writeToFile() const
Flag to allow writing to file.
Definition writeFile.C:286
bool store()
Register object with its registry and transfer ownership to the registry.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
auto & name
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
constexpr auto key(const Type &t) noexcept
Helper function to return the enum value.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
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.
Definition POSIX.C:801
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
const dimensionSet dimForce
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
static bool initialised_(false)
dictionary dict
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.
Definition stdFoam.H:235