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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2015 OpenFOAM Foundation
9 Copyright (C) 2020 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 "polynomialFunction.H"
31
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
37}
38
39
40// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
41
42Foam::polynomialFunction Foam::polynomialFunction::cloneIntegral
43(
44 const polynomialFunction& poly,
45 const scalar intConstant
46)
47{
48 polynomialFunction newPoly(poly.size()+1);
49
50 newPoly[0] = intConstant;
51 forAll(poly, i)
52 {
53 newPoly[i+1] = poly[i]/(i + 1);
54 }
55
56 return newPoly;
57}
58
59
60Foam::polynomialFunction Foam::polynomialFunction::cloneIntegralMinus1
61(
62 const polynomialFunction& poly,
63 const scalar intConstant
64)
65{
66 polynomialFunction newPoly(poly.size()+1);
67
68 if (poly[0] > VSMALL)
69 {
70 newPoly.logActive_ = true;
71 newPoly.logCoeff_ = poly[0];
72 }
73
74 newPoly[0] = intConstant;
75 for (label i=1; i < poly.size(); ++i)
76 {
77 newPoly[i] = poly[i]/i;
78 }
79
80 return newPoly;
81}
82
83
84// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
85
86void Foam::polynomialFunction::checkSize() const
87{
88 if (this->empty())
89 {
91 << "polynomialFunction coefficients are invalid (empty)"
92 << nl << exit(FatalError);
93 }
94}
95
96
97// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98
100:
101 scalarList(1, Zero),
102 logActive_(false),
103 logCoeff_(0)
104{}
105
106
108:
109 scalarList(order, Zero),
110 logActive_(false),
111 logCoeff_(0)
112{
113 checkSize();
114}
115
116
118(
119 const std::initializer_list<scalar> coeffs
120)
121:
122 scalarList(coeffs),
123 logActive_(false),
124 logCoeff_(0)
125{
126 checkSize();
127}
128
129
131:
132 scalarList(coeffs),
133 logActive_(false),
134 logCoeff_(0)
135{
136 checkSize();
137}
138
139
141:
142 scalarList(is),
143 logActive_(false),
144 logCoeff_(0)
145{
146 checkSize();
147}
148
149
150// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151
153{
154 return logActive_;
155}
156
157
159{
160 return logCoeff_;
161}
162
163
164Foam::scalar Foam::polynomialFunction::value(const scalar x) const
165{
166 const scalarList& coeffs = *this;
167 scalar val = coeffs[0];
168
169 // Avoid costly pow() in calculation
170 scalar powX = x;
171 for (label i=1; i<coeffs.size(); ++i)
172 {
173 val += coeffs[i]*powX;
174 powX *= x;
175 }
176
177 if (logActive_)
178 {
179 val += this->logCoeff_*log(x);
180 }
181
182 return val;
183}
184
185
187(
188 const scalar x1,
189 const scalar x2
190) const
191{
192 const scalarList& coeffs = *this;
193
194 if (logActive_)
195 {
197 << "Cannot integrate polynomial with logarithmic coefficients"
198 << nl << abort(FatalError);
199 }
200
201 // Avoid costly pow() in calculation
202 scalar powX1 = x1;
203 scalar powX2 = x2;
204
205 scalar val = coeffs[0]*(powX2 - powX1);
206 for (label i=1; i<coeffs.size(); ++i)
207 {
208 val += coeffs[i]/(i + 1)*(powX2 - powX1);
209 powX1 *= x1;
210 powX2 *= x2;
211 }
212
213 return val;
214}
215
216
218Foam::polynomialFunction::integral(const scalar intConstant) const
219{
220 return cloneIntegral(*this, intConstant);
221}
222
223
225Foam::polynomialFunction::integralMinus1(const scalar intConstant) const
226{
227 return cloneIntegralMinus1(*this, intConstant);
228}
229
230
231// * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
232
234{
235 return
236 (
237 scalarList::operator==(static_cast<const scalarList&>(rhs))
238 && logActive_ == rhs.logActive_
239 && (!logActive_ || (logCoeff_ == rhs.logCoeff_))
240 );
241}
242
243
246{
247 scalarList& coeffs = *this;
248
249 if (coeffs.size() < poly.size())
250 {
251 coeffs.resize(poly.size(), Zero);
252 }
253
254 forAll(poly, i)
255 {
256 coeffs[i] += poly[i];
257 }
258
259 return *this;
260}
261
262
265{
266 scalarList& coeffs = *this;
267
268 if (coeffs.size() < poly.size())
269 {
270 coeffs.resize(poly.size(), Zero);
271 }
272
273 forAll(poly, i)
274 {
275 coeffs[i] -= poly[i];
276 }
277
278 return *this;
279}
280
281
284{
285 scalarList& coeffs = *this;
286 forAll(coeffs, i)
287 {
288 coeffs[i] *= s;
289 }
290
291 return *this;
292}
293
294
297{
298 scalarList& coeffs = *this;
299 forAll(coeffs, i)
300 {
301 coeffs[i] /= s;
302 }
303
304 return *this;
305}
306
307
308// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
309
311{
312 // Log handling may be unreliable
313 poly.logActive_ = false;
314 poly.logCoeff_ = 0;
315
316 is >> static_cast<scalarList&>(poly);
317
318 poly.checkSize();
319
320 return is;
321}
322
323
325{
326 // output like VectorSpace
328
329 forAll(poly, i)
330 {
331 if (i) os << token::SPACE;
332 os << poly[i];
333 }
335
337
338 return os;
339}
340
341
342// * * * * * * * * * * * * * * * Global Operators * * * * * * * * * * * * * //
343
345Foam::operator+
346(
347 const polynomialFunction& p1,
348 const polynomialFunction& p2
349)
350{
351 polynomialFunction poly(p1);
352 return poly += p2;
353}
354
355
357Foam::operator-
358(
359 const polynomialFunction& p1,
360 const polynomialFunction& p2
361)
362{
363 polynomialFunction poly(p1);
364 return poly -= p2;
365}
366
367
369Foam::operator*
370(
371 const scalar s,
372 const polynomialFunction& p
373)
374{
375 polynomialFunction poly(p);
376 return poly *= s;
377}
378
379
381Foam::operator/
382(
383 const scalar s,
384 const polynomialFunction& p
385)
386{
387 polynomialFunction poly(p);
388 return poly /= s;
389}
390
391
393Foam::operator*
394(
395 const polynomialFunction& p,
396 const scalar s
397)
398{
399 polynomialFunction poly(p);
400 return poly *= s;
401}
402
403
405Foam::operator/
406(
407 const polynomialFunction& p,
408 const scalar s
409)
410{
411 polynomialFunction poly(p);
412 return poly /= s;
413}
414
415
416// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:58
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
friend Ostream & operator(Ostream &, const faMatrix< Type > &)
Polynomial function representation.
scalar logCoeff() const
The log coefficient.
scalar integrate(const scalar x1, const scalar x2) const
Integrate between two values.
polynomialFunction & operator/=(const scalar)
scalar value(const scalar x) const
Return polynomial value.
polynomialFunction & operator+=(const polynomialFunction &)
polynomialFunction()
Default construct as size 1 with coefficient == 0.
polynomialFunction & operator-=(const polynomialFunction &)
polynomialFunction & operator*=(const scalar)
bool logActive() const
True if the log term is active.
polynomialFunction integralMinus1(const scalar intConstant=0) const
Return integral coefficients when lowest order is -1.
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
@ BEGIN_LIST
Begin list [isseparator].
Definition: token.H:155
@ END_LIST
End list [isseparator].
Definition: token.H:156
@ SPACE
Space [isspace].
Definition: token.H:125
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
#define FUNCTION_NAME
Namespace for OpenFOAM.
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
Istream & operator>>(Istream &, directionInfo &)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333