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-2016 OpenFOAM Foundation
6  \\/ M anipulation |
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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace functionObjects
34 {
35  defineTypeNameAndDebug(forceCoeffs, 0);
36  addToRunTimeSelectionTable(functionObject, forceCoeffs, dictionary);
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
42 
44 {
45  if (i == 0)
46  {
47  // force coeff data
48 
49  writeHeader(file(i), "Force coefficients");
50  writeHeaderValue(file(i), "liftDir", liftDir_);
51  writeHeaderValue(file(i), "dragDir", dragDir_);
52  writeHeaderValue(file(i), "pitchAxis", pitchAxis_);
53  writeHeaderValue(file(i), "magUInf", magUInf_);
54  writeHeaderValue(file(i), "lRef", lRef_);
55  writeHeaderValue(file(i), "Aref", Aref_);
56  writeHeaderValue(file(i), "CofR", coordSys_.origin());
57  writeCommented(file(i), "Time");
58  writeTabbed(file(i), "Cm");
59  writeTabbed(file(i), "Cd");
60  writeTabbed(file(i), "Cl");
61  writeTabbed(file(i), "Cl(f)");
62  writeTabbed(file(i), "Cl(r)");
63  }
64  else if (i == 1)
65  {
66  // bin coeff data
67 
68  writeHeader(file(i), "Force coefficient bins");
69  writeHeaderValue(file(i), "bins", nBin_);
70  writeHeaderValue(file(i), "start", binMin_);
71  writeHeaderValue(file(i), "delta", binDx_);
72  writeHeaderValue(file(i), "direction", binDir_);
73 
74  vectorField binPoints(nBin_);
75  writeCommented(file(i), "x co-ords :");
76  forAll(binPoints, pointi)
77  {
78  binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
79  file(i) << tab << binPoints[pointi].x();
80  }
81  file(i) << nl;
82 
83  writeCommented(file(i), "y co-ords :");
84  forAll(binPoints, pointi)
85  {
86  file(i) << tab << binPoints[pointi].y();
87  }
88  file(i) << nl;
89 
90  writeCommented(file(i), "z co-ords :");
91  forAll(binPoints, pointi)
92  {
93  file(i) << tab << binPoints[pointi].z();
94  }
95  file(i) << nl;
96 
97  writeCommented(file(i), "Time");
98 
99  for (label j = 0; j < nBin_; j++)
100  {
101  const word jn('(' + Foam::name(j) + ')');
102  writeTabbed(file(i), "Cm" + jn);
103  writeTabbed(file(i), "Cd" + jn);
104  writeTabbed(file(i), "Cl" + jn);
105  }
106  }
107  else
108  {
110  << "Unhandled file index: " << i
111  << abort(FatalError);
112  }
113 
114  file(i)<< endl;
115 }
116 
117 
118 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119 
120 Foam::functionObjects::forceCoeffs::forceCoeffs
121 (
122  const word& name,
123  const Time& runTime,
124  const dictionary& dict
125 )
126 :
127  forces(name, runTime, dict),
128  liftDir_(Zero),
129  dragDir_(Zero),
130  pitchAxis_(Zero),
131  magUInf_(0.0),
132  lRef_(0.0),
133  Aref_(0.0)
134 {
135  read(dict);
136 }
137 
138 
139 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
140 
142 {}
143 
144 
145 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
146 
148 {
149  forces::read(dict);
150 
151  // Directions for lift and drag forces, and pitch moment
152  dict.lookup("liftDir") >> liftDir_;
153  dict.lookup("dragDir") >> dragDir_;
154  dict.lookup("pitchAxis") >> pitchAxis_;
155 
156  // Free stream velocity magnitude
157  dict.lookup("magUInf") >> magUInf_;
158 
159  // Reference length and area scales
160  dict.lookup("lRef") >> lRef_;
161  dict.lookup("Aref") >> Aref_;
162 
163  return true;
164 }
165 
166 
168 {
169  return true;
170 }
171 
172 
174 {
176 
177  if (Pstream::master())
178  {
180 
181  scalar pDyn = 0.5*rhoRef_*magUInf_*magUInf_;
182 
183  Field<vector> totForce(force_[0] + force_[1] + force_[2]);
184  Field<vector> totMoment(moment_[0] + moment_[1] + moment_[2]);
185 
186  List<Field<scalar>> coeffs(3);
187  coeffs[0].setSize(nBin_);
188  coeffs[1].setSize(nBin_);
189  coeffs[2].setSize(nBin_);
190 
191  // lift, drag and moment
192  coeffs[0] = (totForce & liftDir_)/(Aref_*pDyn);
193  coeffs[1] = (totForce & dragDir_)/(Aref_*pDyn);
194  coeffs[2] = (totMoment & pitchAxis_)/(Aref_*lRef_*pDyn);
195 
196  scalar Cl = sum(coeffs[0]);
197  scalar Cd = sum(coeffs[1]);
198  scalar Cm = sum(coeffs[2]);
199 
200  scalar Clf = Cl/2.0 + Cm;
201  scalar Clr = Cl/2.0 - Cm;
202 
203  writeTime(file(0));
204  file(0)
205  << tab << Cm << tab << Cd
206  << tab << Cl << tab << Clf << tab << Clr << endl;
207 
208  Log << type() << " " << name() << " write:" << nl
209  << " Cm = " << Cm << nl
210  << " Cd = " << Cd << nl
211  << " Cl = " << Cl << nl
212  << " Cl(f) = " << Clf << nl
213  << " Cl(r) = " << Clr << endl;
214 
215  if (nBin_ > 1)
216  {
217  if (binCumulative_)
218  {
219  for (label i = 1; i < coeffs[0].size(); i++)
220  {
221  coeffs[0][i] += coeffs[0][i-1];
222  coeffs[1][i] += coeffs[1][i-1];
223  coeffs[2][i] += coeffs[2][i-1];
224  }
225  }
226 
227  writeTime(file(1));
228 
229  forAll(coeffs[0], i)
230  {
231  file(1)
232  << tab << coeffs[2][i]
233  << tab << coeffs[1][i]
234  << tab << coeffs[0][i];
235  }
236 
237  file(1) << endl;
238  }
239 
240  Log << endl;
241  }
242 
243  return true;
244 }
245 
246 
247 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
virtual bool write()
Write the forces.
Definition: forceCoeffs.C:173
static const char tab
Definition: Ostream.H:261
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
virtual bool read(const dictionary &)
Read the forces data.
Definition: forces.C:612
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual bool read(const dictionary &)
Read the forces data.
Definition: forceCoeffs.C:147
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class for handling words, derived from string.
Definition: word.H:59
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:747
static const zero Zero
Definition: zero.H:91
errorManip< error > abort(error &err)
Definition: errorManip.H:131
static const char nl
Definition: Ostream.H:262
virtual void writeFileHeader(const label i)
Output file header information.
Definition: forceCoeffs.C:43
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setSize(const label)
Reset size of List.
Definition: List.C:295
virtual ~forceCoeffs()
Destructor.
Definition: forceCoeffs.C:141
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
addToRunTimeSelectionTable(functionObject, blendingFactor, dictionary)
virtual bool write()
Write function.
Definition: writeFiles.C:197
#define Log
Report write to Foam::Info if the local log switch is true.
volScalarField pDyn(IOobject("pDyn", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, dimensionedScalar("zero", dimPressure, 0.0))
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
This function object calculates the forces and moments by integrating the pressure and skin-friction ...
Definition: forces.H:196
virtual bool execute()
Execute, currently does nothing.
Definition: forceCoeffs.C:167
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451