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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "forceCoeffs.H"
27 #include "dictionary.H"
28 #include "Time.H"
29 #include "Pstream.H"
30 #include "IOmanip.H"
31 #include "fvMesh.H"
32 #include "dimensionedTypes.H"
33 #include "volFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(forceCoeffs, 0);
40 }
41 
42 
43 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
44 
46 {
47  // Note: Only possible to create bin files after bins have been initialised
48 
49  if (writeToFile() && !coeffFilePtr_.valid())
50  {
51  coeffFilePtr_ = createFile("coefficient");
52  writeIntegratedHeader("Coefficients", coeffFilePtr_());
53 
54  if (nBin_ > 1)
55  {
56  CmBinFilePtr_ = createFile("CmBin");
57  writeBinHeader("Moment coefficient bins", CmBinFilePtr_());
58  CdBinFilePtr_ = createFile("CdBin");
59  writeBinHeader("Drag coefficient bins", CdBinFilePtr_());
60  ClBinFilePtr_ = createFile("ClBin");
61  writeBinHeader("Lift coefficient bins", ClBinFilePtr_());
62  }
63  }
64 }
65 
66 
68 (
69  const word& header,
70  Ostream& os
71 ) const
72 {
73  writeHeader(os, "Force coefficients");
74  writeHeaderValue(os, "liftDir", liftDir_);
75  writeHeaderValue(os, "dragDir", dragDir_);
76  writeHeaderValue(os, "pitchAxis", pitchAxis_);
77  writeHeaderValue(os, "magUInf", magUInf_);
78  writeHeaderValue(os, "lRef", lRef_);
79  writeHeaderValue(os, "Aref", Aref_);
80  writeHeaderValue(os, "CofR", coordSys_.origin());
81  writeHeader(os, "");
82  writeCommented(os, "Time");
83  writeTabbed(os, "Cm");
84  writeTabbed(os, "Cd");
85  writeTabbed(os, "Cl");
86  writeTabbed(os, "Cl(f)");
87  writeTabbed(os, "Cl(r)");
88  os << endl;
89 }
90 
91 
93 (
94  const word& header,
95  Ostream& os
96 ) const
97 {
98  writeHeader(os, header);
99  writeHeaderValue(os, "bins", nBin_);
100  writeHeaderValue(os, "start", binMin_);
101  writeHeaderValue(os, "delta", binDx_);
102  writeHeaderValue(os, "direction", binDir_);
103 
104  vectorField binPoints(nBin_);
105  writeCommented(os, "x co-ords :");
106  forAll(binPoints, pointI)
107  {
108  binPoints[pointI] = (binMin_ + (pointI + 1)*binDx_)*binDir_;
109  os << tab << binPoints[pointI].x();
110  }
111  os << nl;
112 
113  writeCommented(os, "y co-ords :");
114  forAll(binPoints, pointI)
115  {
116  os << tab << binPoints[pointI].y();
117  }
118  os << nl;
119 
120  writeCommented(os, "z co-ords :");
121  forAll(binPoints, pointI)
122  {
123  os << tab << binPoints[pointI].z();
124  }
125  os << nl;
126 
127  writeHeader(os, "");
128  writeCommented(os, "Time");
129 
130  for (label j = 0; j < nBin_; j++)
131  {
132  word jn(Foam::name(j) + ':');
133  writeTabbed(os, jn + "total");
134  writeTabbed(os, jn + "pressure");
135  writeTabbed(os, jn + "viscous");
136 
137  if (porosity_)
138  {
139  writeTabbed(os, jn + "porous");
140  }
141  }
142 
143  os << endl;
144 }
145 
146 
148 (
149  const word& title,
150  const List<Field<scalar> >& coeff
151 ) const
152 {
153  scalar pressure = sum(coeff[0]);
154  scalar viscous = sum(coeff[1]);
155  scalar porous = sum(coeff[2]);
156  scalar total = pressure + viscous + porous;
157 
158  if (log_)
159  {
160  Info<< " " << title << " : " << total << token::TAB
161  << "("
162  << "pressure: " << pressure << token::TAB
163  << "viscous: " << viscous;
164 
165  if (porosity_)
166  {
167  Info<< token::TAB << "porous: " << porous;
168  }
169 
170  Info<< ")" << endl;
171  }
172 }
173 
174 
176 (
177  const List<Field<scalar> > coeffs,
178  Ostream& os
179 ) const
180 {
181  os << obr_.time().value();
182 
183  for (label binI = 0; binI < nBin_; binI++)
184  {
185  scalar total = coeffs[0][binI] + coeffs[1][binI] + coeffs[2][binI];
186 
187  os << tab << total << tab << coeffs[0][binI] << tab << coeffs[1][binI];
188 
189  if (porosity_)
190  {
191  os << tab << coeffs[2][binI];
192  }
193  }
194 
195  os << endl;
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
200 
202 (
203  const word& name,
204  const objectRegistry& obr,
205  const dictionary& dict,
206  const bool loadFromFiles,
207  const bool readFields
208 )
209 :
210  forces(name, obr, dict, loadFromFiles, false),
211  liftDir_(vector::zero),
212  dragDir_(vector::zero),
213  pitchAxis_(vector::zero),
214  magUInf_(0.0),
215  lRef_(0.0),
216  Aref_(0.0),
217  coeffFilePtr_(),
218  CmBinFilePtr_(),
219  CdBinFilePtr_(),
220  ClBinFilePtr_()
221 {
222  if (readFields)
223  {
224  read(dict);
225  if (log_) Info << endl;
226  }
227 }
228 
229 
230 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
231 
233 {}
234 
235 
236 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
237 
239 {
240  if (!active_)
241  {
242  return;
243  }
244 
246 
247  // Directions for lift and drag forces, and pitch moment
248  dict.lookup("liftDir") >> liftDir_;
249  dict.lookup("dragDir") >> dragDir_;
250  dict.lookup("pitchAxis") >> pitchAxis_;
251 
252  // Free stream velocity magnitude
253  dict.lookup("magUInf") >> magUInf_;
254 
255  // Reference length and area scales
256  dict.lookup("lRef") >> lRef_;
257  dict.lookup("Aref") >> Aref_;
258 
259  if (writeFields_)
260  {
261  const fvMesh& mesh = refCast<const fvMesh>(obr_);
262 
263  tmp<volVectorField> tforceCoeff
264  (
265  new volVectorField
266  (
267  IOobject
268  (
269  fieldName("forceCoeff"),
270  mesh.time().timeName(),
271  mesh,
274  ),
275  mesh,
277  )
278  );
279 
280  obr_.store(tforceCoeff.ptr());
281 
282  tmp<volVectorField> tmomentCoeff
283  (
284  new volVectorField
285  (
286  IOobject
287  (
288  fieldName("momentCoeff"),
289  mesh.time().timeName(),
290  mesh,
293  ),
294  mesh,
296  )
297  );
298 
299  obr_.store(tmomentCoeff.ptr());
300  }
301 }
302 
303 
305 {
306  if (!active_)
307  {
308  return;
309  }
310 
312 
313  createFiles();
314 
315  scalar pDyn = 0.5*rhoRef_*magUInf_*magUInf_;
316 
317  // Storage for pressure, viscous and porous contributions to coeffs
318  List<Field<scalar> > momentCoeffs(3);
320  List<Field<scalar> > liftCoeffs(3);
321  forAll(liftCoeffs, i)
322  {
323  momentCoeffs[i].setSize(nBin_);
324  dragCoeffs[i].setSize(nBin_);
325  liftCoeffs[i].setSize(nBin_);
326  }
327 
328  // Calculate coefficients
329  scalar CmTot = 0;
330  scalar CdTot = 0;
331  scalar ClTot = 0;
332  forAll(liftCoeffs, i)
333  {
334  momentCoeffs[i] = (moment_[i] & pitchAxis_)/(Aref_*pDyn*lRef_);
335  dragCoeffs[i] = (force_[i] & dragDir_)/(Aref_*pDyn);
336  liftCoeffs[i] = (force_[i] & liftDir_)/(Aref_*pDyn);
337 
338  CmTot += sum(momentCoeffs[i]);
339  CdTot += sum(dragCoeffs[i]);
340  ClTot += sum(liftCoeffs[i]);
341  }
342 
343  scalar ClfTot = ClTot/2.0 + CmTot;
344  scalar ClrTot = ClTot/2.0 - CmTot;
345 
346  if (log_) Info
347  << type() << " " << name_ << " output:" << nl
348  << " Coefficients" << nl;
349 
350  writeIntegratedData("Cm", momentCoeffs);
351  writeIntegratedData("Cd", dragCoeffs);
352  writeIntegratedData("Cl", liftCoeffs);
353 
354  if (log_) Info
355  << " Cl(f) : " << ClfTot << nl
356  << " Cl(r) : " << ClrTot << nl
357  << endl;
358 
359  if (writeToFile())
360  {
361  writeTime(coeffFilePtr_());
362  coeffFilePtr_()
363  << tab << CmTot << tab << CdTot
364  << tab << ClTot << tab << ClfTot << tab << ClrTot << endl;
365 
366 
367  if (nBin_ > 1)
368  {
369  if (binCumulative_)
370  {
371  forAll(liftCoeffs, i)
372  {
373  for (label binI = 1; binI < nBin_; binI++)
374  {
375  liftCoeffs[i][binI] += liftCoeffs[i][binI-1];
376  dragCoeffs[i][binI] += dragCoeffs[i][binI-1];
377  momentCoeffs[i][binI] += momentCoeffs[i][binI-1];
378  }
379  }
380  }
381 
382  writeBinData(dragCoeffs, CdBinFilePtr_());
383  writeBinData(liftCoeffs, ClBinFilePtr_());
384  writeBinData(momentCoeffs, CmBinFilePtr_());
385  }
386  }
387 
388  // Write state/results information
389  {
390  setResult("Cm", CmTot);
391  setResult("Cd", CdTot);
392  setResult("Cl", ClTot);
393  setResult("Cl(f)", ClfTot);
394  setResult("Cl(r)", ClrTot);
395  }
396 
397  if (writeFields_)
398  {
399  const volVectorField& force =
400  obr_.lookupObject<volVectorField>(fieldName("force"));
401 
402  const volVectorField& moment =
403  obr_.lookupObject<volVectorField>(fieldName("moment"));
404 
405  volVectorField& forceCoeff =
406  const_cast<volVectorField&>
407  (
408  obr_.lookupObject<volVectorField>(fieldName("forceCoeff"))
409  );
410 
411  volVectorField& momentCoeff =
412  const_cast<volVectorField&>
413  (
414  obr_.lookupObject<volVectorField>(fieldName("momentCoeff"))
415  );
416 
417  dimensionedScalar f0("f0", dimForce, Aref_*pDyn);
418  dimensionedScalar m0("m0", dimForce*dimLength, Aref_*lRef_*pDyn);
419 
420  forceCoeff == force/f0;
421  momentCoeff == moment/m0;
422  }
423 }
424 
425 
427 {
428  // Do nothing
429 }
430 
431 
433 {
434  // Do nothing
435 }
436 
437 
439 {
440  if (!active_)
441  {
442  return;
443  }
444 
445  if (writeFields_)
446  {
447  const volVectorField& forceCoeff =
448  obr_.lookupObject<volVectorField>(fieldName("forceCoeff"));
449 
450  const volVectorField& momentCoeff =
451  obr_.lookupObject<volVectorField>(fieldName("momentCoeff"));
452 
453  forceCoeff.write();
454  momentCoeff.write();
455  }
456 }
457 
458 
459 // ************************************************************************* //
Foam::token::TAB
@ TAB
Definition: token.H:96
volFields.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
forceCoeffs.H
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::jn
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:296
Foam::forces::nBin_
label nBin_
Number of bins.
Definition: forces.H:321
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::forceCoeffs::writeIntegratedHeader
void writeIntegratedHeader(const word &header, Ostream &os) const
Write header for integrated data.
Definition: forceCoeffs.C:68
Foam::forceCoeffs::read
virtual void read(const dictionary &)
Read the forces data.
Definition: forceCoeffs.C:238
dragCoeffs
autoPtr< multiphaseSystem::dragCoeffFields > dragCoeffs(fluid.dragCoeffs())
Foam::dimForce
const dimensionSet dimForce
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::forceCoeffs::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: forceCoeffs.C:432
Foam::forceCoeffs::writeBinHeader
void writeBinHeader(const word &header, Ostream &os) const
Write header for binned data.
Definition: forceCoeffs.C:93
IOmanip.H
Istream and Ostream manipulators taking arguments.
Foam::forceCoeffs::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: forceCoeffs.C:426
pDyn
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimPressure, 0.0))
Foam::forceCoeffs::writeIntegratedData
void writeIntegratedData(const word &title, const List< Field< scalar > > &coeff) const
Write integrated data.
Definition: forceCoeffs.C:148
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::forceCoeffs::createFiles
void createFiles()
Create the output files.
Definition: forceCoeffs.C:45
Foam::forceCoeffs::~forceCoeffs
virtual ~forceCoeffs()
Destructor.
Definition: forceCoeffs.C:232
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Pstream.H
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::forces::calcForcesMoment
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:1063
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::tmp::ptr
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:146
coeff
const scalar coeff[]
Definition: Test-Polynomial.C:38
Foam::forceCoeffs::coeffFilePtr_
autoPtr< OFstream > coeffFilePtr_
Integrated coefficients.
Definition: forceCoeffs.H:233
Foam::functionObjectFile::writeToFile
bool writeToFile() const
Return true if can write to file.
Definition: functionObjectFile.C:206
Foam::forceCoeffs::CdBinFilePtr_
autoPtr< OFstream > CdBinFilePtr_
Drag coefficient.
Definition: forceCoeffs.H:239
Foam::forceCoeffs::write
virtual void write()
Write the forces.
Definition: forceCoeffs.C:438
Foam::forces
This function object calculates the forces and moments by integrating the pressure and skin-friction ...
Definition: forces.H:234
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::tab
static const char tab
Definition: Ostream.H:259
Foam::functionObjectFile::createFile
virtual autoPtr< OFstream > createFile(const word &name) const
Return an autoPtr to a new file.
Definition: functionObjectFile.C:83
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::forceCoeffs::forceCoeffs
forceCoeffs(const forceCoeffs &)
Disallow default bitwise copy construct.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dictionary.H
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
dimensionedTypes.H
Foam::forceCoeffs::execute
virtual void execute()
Execute, currently does nothing.
Definition: forceCoeffs.C:304
Foam::forceCoeffs::ClBinFilePtr_
autoPtr< OFstream > ClBinFilePtr_
Lift coefficient.
Definition: forceCoeffs.H:242
Foam::forces::read
virtual void read(const dictionary &)
Read the forces data.
Definition: forces.C:858
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::forceCoeffs::CmBinFilePtr_
autoPtr< OFstream > CmBinFilePtr_
Moment coefficient.
Definition: forceCoeffs.H:236
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::forceCoeffs::writeBinData
void writeBinData(const List< Field< scalar > > coeffs, Ostream &os) const
Write binned data.
Definition: forceCoeffs.C:176