polynomialFunction.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 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 "polynomialFunction.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineTypeNameAndDebug(polynomialFunction, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38 
39 
40 Foam::polynomialFunction Foam::polynomialFunction::cloneIntegral
41 (
42  const polynomialFunction& poly,
43  const scalar intConstant
44 )
45 {
46  polynomialFunction newPoly(poly.size()+1);
47 
48  newPoly[0] = intConstant;
49  forAll(poly, i)
50  {
51  newPoly[i+1] = poly[i]/(i + 1);
52  }
53 
54  return newPoly;
55 }
56 
57 
58 Foam::polynomialFunction Foam::polynomialFunction::cloneIntegralMinus1
59 (
60  const polynomialFunction& poly,
61  const scalar intConstant
62 )
63 {
64  polynomialFunction newPoly(poly.size()+1);
65 
66  if (poly[0] > VSMALL)
67  {
68  newPoly.logActive_ = true;
69  newPoly.logCoeff_ = poly[0];
70  }
71 
72  newPoly[0] = intConstant;
73  for (label i=1; i < poly.size(); ++i)
74  {
75  newPoly[i] = poly[i]/i;
76  }
77 
78  return newPoly;
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
85 :
86  scalarList(order, 0.0),
87  logActive_(false),
88  logCoeff_(0.0)
89 {
90  if (this->empty())
91  {
93  (
94  "polynomialFunction::polynomialFunction(const label order)"
95  ) << "polynomialFunction coefficients are invalid (empty)"
96  << nl << exit(FatalError);
97  }
98 }
99 
100 
102 :
103  scalarList(poly),
104  logActive_(poly.logActive_),
105  logCoeff_(poly.logCoeff_)
106 {}
107 
108 
110 :
111  scalarList(coeffs),
112  logActive_(false),
113  logCoeff_(0.0)
114 {
115  if (this->empty())
116  {
118  (
119  "polynomialFunction::polynomialFunction(const UList<scalar>&)"
120  ) << "polynomialFunction coefficients are invalid (empty)"
121  << nl << exit(FatalError);
122  }
123 }
124 
125 
127 :
128  scalarList(is),
129  logActive_(false),
130  logCoeff_(0.0)
131 {
132  if (this->empty())
133  {
135  (
136  "polynomialFunction::polynomialFunction(Istream&)"
137  ) << "polynomialFunction coefficients are invalid (empty)"
138  << nl << exit(FatalError);
139  }
140 }
141 
142 
143 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
144 
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
152 {
153  return logActive_;
154 }
155 
156 
158 {
159  return logCoeff_;
160 }
161 
162 
163 Foam::scalar Foam::polynomialFunction::value(const scalar x) const
164 {
165  const scalarList& coeffs = *this;
166  scalar val = coeffs[0];
167 
168  // avoid costly pow() in calculation
169  scalar powX = x;
170  for (label i=1; i<coeffs.size(); ++i)
171  {
172  val += coeffs[i]*powX;
173  powX *= x;
174  }
175 
176  if (logActive_)
177  {
178  val += this->logCoeff_*log(x);
179  }
180 
181  return val;
182 }
183 
184 
186 (
187  const scalar x1,
188  const scalar x2
189 ) const
190 {
191  const scalarList& coeffs = *this;
192 
193  if (logActive_)
194  {
196  (
197  "scalar polynomialFunction::integrate"
198  "("
199  "const scalar, "
200  "const scalar"
201  ") const"
202  ) << "Cannot integrate polynomial with logarithmic coefficients"
203  << nl << abort(FatalError);
204  }
205 
206  // avoid costly pow() in calculation
207  scalar powX1 = x1;
208  scalar powX2 = x2;
209 
210  scalar val = coeffs[0]*(powX2 - powX1);
211  for (label i=1; i<coeffs.size(); ++i)
212  {
213  val += coeffs[i]/(i + 1)*(powX2 - powX1);
214  powX1 *= x1;
215  powX2 *= x2;
216  }
217 
218  return val;
219 }
220 
221 
223 Foam::polynomialFunction::integral(const scalar intConstant) const
224 {
225  return cloneIntegral(*this, intConstant);
226 }
227 
228 
230 Foam::polynomialFunction::integralMinus1(const scalar intConstant) const
231 {
232  return cloneIntegralMinus1(*this, intConstant);
233 }
234 
235 
236 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
237 
240 {
241  scalarList& coeffs = *this;
242 
243  if (coeffs.size() > poly.size())
244  {
245  forAll(poly, i)
246  {
247  coeffs[i] += poly[i];
248  }
249  }
250  else
251  {
252  coeffs.setSize(poly.size(), 0.0);
253 
254  forAll(coeffs, i)
255  {
256  coeffs[i] += poly[i];
257  }
258  }
259 
260  return *this;
261 }
262 
263 
266 {
267  scalarList& coeffs = *this;
268 
269  if (coeffs.size() > poly.size())
270  {
271  forAll(poly, i)
272  {
273  coeffs[i] -= poly[i];
274  }
275  }
276  else
277  {
278  coeffs.setSize(poly.size(), 0.0);
279 
280  forAll(coeffs, i)
281  {
282  coeffs[i] -= poly[i];
283  }
284  }
285 
286  return *this;
287 }
288 
289 
292 {
293  scalarList& coeffs = *this;
294  forAll(coeffs, i)
295  {
296  coeffs[i] *= s;
297  }
298 
299  return *this;
300 }
301 
302 
305 {
306  scalarList& coeffs = *this;
307  forAll(coeffs, i)
308  {
309  coeffs[i] /= s;
310  }
311 
312  return *this;
313 }
314 
315 
316 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
317 
319 {
320  // output like VectorSpace
321  os << token::BEGIN_LIST;
322 
323  if (!poly.empty())
324  {
325  for (int i=0; i<poly.size()-1; i++)
326  {
327  os << poly[i] << token::SPACE;
328  }
329  os << poly.last();
330  }
331  os << token::END_LIST;
332 
333 
334  // Check state of Ostream
335  os.check("operator<<(Ostream&, const polynomialFunction&)");
336 
337  return os;
338 }
339 
340 
341 // * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
342 
344 Foam::operator+
345 (
346  const polynomialFunction& p1,
347  const polynomialFunction& p2
348 )
349 {
350  polynomialFunction poly(p1);
351  return poly += p2;
352 }
353 
354 
356 Foam::operator-
357 (
358  const polynomialFunction& p1,
359  const polynomialFunction& p2
360 )
361 {
362  polynomialFunction poly(p1);
363  return poly -= p2;
364 }
365 
366 
368 Foam::operator*
369 (
370  const scalar s,
371  const polynomialFunction& p
372 )
373 {
374  polynomialFunction poly(p);
375  return poly *= s;
376 }
377 
378 
380 Foam::operator/
381 (
382  const scalar s,
383  const polynomialFunction& p
384 )
385 {
386  polynomialFunction poly(p);
387  return poly /= s;
388 }
389 
390 
392 Foam::operator*
393 (
394  const polynomialFunction& p,
395  const scalar s
396 )
397 {
398  polynomialFunction poly(p);
399  return poly *= s;
400 }
401 
402 
404 Foam::operator/
405 (
406  const polynomialFunction& p,
407  const scalar s
408 )
409 {
410  polynomialFunction poly(p);
411  return poly /= s;
412 }
413 
414 
415 // ************************************************************************* //
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject( name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE ))
bool empty() const
Return true if the UList is empty (ie, size() is zero).
Definition: UListI.H:313
T & last()
Return the last element of the list.
Definition: UListI.H:131
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
bool logActive() const
Return true if the log term is active.
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:133
Namespace for OpenFOAM.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
scalar value(const scalar x) const
Return polynomial value.
dimensionedScalar log(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:260
scalar logCoeff() const
Return the log coefficient.
void setSize(const label)
Reset size of List.
Definition: List.C:318
Polynomial function representation.
volScalarField & p
Definition: createFields.H:51
#define forAll(list, i)
Definition: UList.H:421
polynomialFunction & operator/=(const scalar)
Macros for easy insertion into run-time selection tables.
polynomialFunction(const label)
Construct a particular size, with all coefficients = 0.0.
polynomialFunction & operator+=(const polynomialFunction &)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
error FatalError
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
polynomialFunction & operator*=(const scalar)
virtual ~polynomialFunction()
Destructor.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
polynomialFunction integralMinus1(const scalar intConstant=0.0) const
Return integral coefficients when lowest order is -1.
polynomialFunction integral(const scalar intConstant=0.0) const
Return integral coefficients.
polynomialFunction & operator-=(const polynomialFunction &)
scalar integrate(const scalar x1, const scalar x2) const
Integrate between two values.
defineTypeNameAndDebug(combustionModel, 0)