unintegrable.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2023 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 "unintegrable.H"
27 #include "quadraticEqn.H"
28 #include "SubField.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
34 (
35  const scalarField& x,
36  const scalarField& y
37 )
38 {
39  tmp<scalarField> tResult(new scalarField(x.size()));
40  scalarField& result = tResult.ref();
41 
42  result[0] = 0;
43 
44  for (label i = 1; i < x.size(); ++ i)
45  {
46  result[i] = result[i - 1] + (x[i] - x[i - 1])*(y[i] + y[i - 1])/2;
47  }
48 
49  return tResult;
50 }
51 
52 
54 (
55  const scalarField& x,
56  const scalarField& y
57 )
58 {
59  tmp<scalarField> tResult(new scalarField(x.size()));
60  scalarField& result = tResult.ref();
61 
62  result[0] = 0;
63 
64  for (label i = 1; i < x.size(); ++ i)
65  {
66  const scalar x0 = x[i - 1], x1 = x[i];
67  const scalar y0 = y[i - 1], y1 = y[i];
68 
69  result[i] =
70  result[i - 1] + (x1 - x0)*(2*x0*y0 + x0*y1 + x1*y0 + 2*x1*y1)/6;
71  }
72 
73  return tResult;
74 }
75 
76 
78 (
79  const Pair<scalar>& x,
80  const Pair<scalar>& Phi,
81  const scalar s
82 )
83 {
84  // Do a linear interpolation
85  scalar f = (s - Phi[0])/(Phi[1] - Phi[0]);
86 
87  // Interpolate
88  return (1 - f)*x[0] + f*x[1];
89 }
90 
91 
93 (
94  const Pair<scalar>& x,
95  const Pair<scalar>& phi,
96  const Pair<scalar>& Phi,
97  const scalar s
98 )
99 {
100  // Do a linear interpolation
101  scalar f = (s - Phi[0])/(Phi[1] - Phi[0]);
102 
103  // Attempt to improve by doing a quadratic interpolation
104  const scalar a = (x[1] - x[0])*(phi[1] - phi[0])/2;
105  const scalar b = Phi[1] - Phi[0] - a;
106  const scalar c = Phi[0] - s;
107  const quadraticEqn eqn(a, b, c);
108  const Roots<2> roots = eqn.roots();
109  forAll(roots, rooti)
110  {
111  if
112  (
113  roots.type(rooti) == rootType::real
114  && roots[rooti] > 0
115  && roots[rooti] < 1
116  )
117  {
118  f = roots[rooti];
119  }
120  }
121 
122  // Interpolate
123  return (1 - f)*x[0] + f*x[1];
124 }
125 
126 
128 (
129  const scalarField& x,
130  const scalarField& Phi,
131  const scalar s
132 )
133 {
134  // Find the interval
135  label i = 0;
136  for (; i < x.size() && Phi[i + 1] < s; ++ i);
137 
138  // Sample
139  return
140  sampleInterval
141  (
142  {x[i], x[i + 1]},
143  {Phi[i], Phi[i + 1]},
144  s
145  );
146 }
147 
148 
150 (
151  const scalarField& x,
152  const scalarField& phi,
153  const scalarField& Phi,
154  const scalar s
155 )
156 {
157  // Find the interval
158  label i = 0;
159  for (; i < x.size() && Phi[i + 1] < s; ++ i);
160 
161  // Sample
162  return
163  sampleInterval
164  (
165  {x[i], x[i + 1]},
166  {phi[i], phi[i + 1]},
167  {Phi[i], Phi[i + 1]},
168  s
169  );
170 }
171 
172 
173 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
174 
175 const Foam::scalarField& Foam::distributions::unintegrable::x() const
176 {
177  if (xPtr_.valid()) return xPtr_();
178 
179  xPtr_.set(new scalarField(n_));
180  scalarField& x = xPtr_();
181 
182  forAll(x, i)
183  {
184  const scalar f = scalar(i)/(n_ + 1);
185  x[i] = (1 - f)*this->min() + f*this->max();
186  }
187 
188  static const scalar tol = sqrt(sqrt(pow3(small)));
189 
190  for (label iter = 0; iter < ceil(std::log2(1/tol)); ++ iter)
191  {
192  const scalarField phi(this->phi(this->q(), x));
193  const scalarField Phi(integrate(x, phi));
194 
195  const scalar dPhi = (Phi.last() - Phi.first())/(n_ - 1);
196 
197  if (distribution::debug)
198  {
199  scalar error = -vGreat;
200  for (label i = 0; i < n_ - 1; ++ i)
201  {
202  error = Foam::max(error, (1 - (Phi[i + 1] - Phi[i])/dPhi)/n_);
203  }
204 
205  Info<< indent << "Interval spacing iteration #" << iter
206  << ", error=" << error << endl;
207  }
208 
209  const scalarField xPrev(x);
210 
211  x.first() = this->min();
212 
213  for (label i0 = 0, i = 1; i0 < n_ - 1; ++ i0)
214  {
215  while (Phi[i0 + 1] > i*dPhi)
216  {
217  const scalar xNext =
218  sampleInterval
219  (
220  {xPrev[i0], xPrev[i0 + 1]},
221  {phi[i0], phi[i0 + 1]},
222  {Phi[i0], Phi[i0 + 1]},
223  i*dPhi
224  );
225 
226  x[i] = (x[i] + xNext)/2;
227 
228  i ++;
229  }
230  }
231 
232  x.last() = this->max();
233  }
234 
235  return x;
236 }
237 
238 const Foam::scalarField& Foam::distributions::unintegrable::PDF() const
239 {
240  if (PDFPtr_.valid()) return PDFPtr_();
241 
242  const Pair<scalar> Phi01 = this->Phi01();
243 
244  PDFPtr_.set((phi(this->q(), xPtr_)/(Phi01[1] - Phi01[0])).ptr());
245 
246  return PDFPtr_();
247 }
248 
249 
250 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
251 
253 (
254  const label q,
255  const scalarField& x
256 ) const
257 {
258  return integrate(x, phi(q, x));
259 }
260 
261 
263 (
264  const label q
265 ) const
266 {
267  const scalarField Phi(this->Phi(this->q(), x()));
268 
269  return Pair<scalar>(Phi.first(), Phi.last());
270 }
271 
272 
274 (
275  const label q
276 ) const
277 {
278  if (q == 0)
279  {
280  const scalarField x(scalarList({min(), max()}));
281  const scalarField Phi(this->Phi(q, x));
282  return Pair<scalar>(Phi.first(), Phi.last());
283  }
284  else
285  {
286  return unintegrable::Phi01(q);
287  }
288 }
289 
290 
293 {
294  if (Phi01Ptr_.valid()) return Phi01Ptr_();
295 
296  Phi01Ptr_.set(new Pair<scalar>(Phi01(this->q())));
297 
298  return Phi01Ptr_();
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
303 
305 (
306  const word& name,
307  const dictionary& dict,
308  Random& rndGen,
309  const label sampleQ
310 )
311 :
312  distribution(name, dict, rndGen, sampleQ),
313  n_((1<<dict.lookupOrDefault<label>("level", 16)) + 1)
314 {}
315 
316 
318 (
319  Random& rndGen,
320  const label Q,
321  const label sampleQ,
322  const label n
323 )
324 :
325  distribution(rndGen, Q, sampleQ),
326  n_(n)
327 {}
328 
329 
331 (
332  const unintegrable& d,
333  const label sampleQ
334 )
335 :
336  distribution(d, sampleQ),
337  n_(d.n_)
338 {}
339 
340 
341 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
342 
344 {}
345 
346 
347 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
348 
350 {
351  const scalar s = this->rndGen_.template sample01<scalar>();
352 
353  const scalar dCDF = scalar(1)/(n_ - 1);
354  const label samplei = floor(s/dCDF);
355  return
356  sampleInterval
357  (
358  {x()[samplei], x()[samplei + 1]},
359  {PDF()[samplei], PDF()[samplei + 1]},
360  {samplei*dCDF, (samplei + 1)*dCDF},
361  s
362  );
363 }
364 
365 
367 {
368  const Pair<scalar> Phi01 = this->Phi01();
369  const scalarField Mu(Phi(this->q() + 1, x()));
370  const Pair<scalar> Mu01(Mu.first(), Mu.last());
371  return (Mu01[1] - Mu01[0])/(Phi01[1] - Phi01[0]);
372 }
373 
374 
375 Foam::tmp<Foam::scalarField> Foam::distributions::unintegrable::PDF
376 (
377  const scalarField& x
378 ) const
379 {
380  const scalarField phi(this->phi(this->q(), x));
381  const Pair<scalar> Phi01 = this->Phi01();
382  return this->clipPDF(x, phi/(Phi01[1] - Phi01[0]));
383 }
384 
385 
386 // ************************************************************************* //
scalar y
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
T & last()
Return the last element of the list.
Definition: FixedListI.H:231
const Type & first() const
Return first.
Definition: Pair.H:98
Random number generator.
Definition: Random.H:58
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition: Roots.H:67
void type(const direction i, const rootType t)
Set the type of the i-th root.
Definition: RootsI.H:120
T & first()
Return the first element of the list.
Definition: UListI.H:114
T & last()
Return the last element of the list.
Definition: UListI.H:128
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
virtual Pair< scalar > Phi01(const label q) const
Access cached values of the un-normalised CDF at the minimum and.
Definition: unintegrable.C:263
Base class for distributions that do not have a closed integral form for the cumulative density funct...
Definition: unintegrable.H:58
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
Definition: unintegrable.C:253
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:349
unintegrable(const word &name, const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: unintegrable.C:305
static tmp< scalarField > integrate(const scalarField &x, const scalarField &y)
Integrate the values y with respect to the coordinates x.
Definition: unintegrable.C:34
virtual ~unintegrable()
Destructor.
Definition: unintegrable.C:343
static tmp< scalarField > integrateX(const scalarField &x, const scalarField &y)
Integrate the values x*y with respect to the coordinates x.
Definition: unintegrable.C:54
const Pair< scalar > & Phi01() const
Access cached values of the un-normalised CDF at the minimum and.
Definition: unintegrable.C:292
virtual scalar mean() const
Return the mean value.
Definition: unintegrable.C:366
static scalar sampleInterval(const Pair< scalar > &x, const Pair< scalar > &Phi, const scalar s)
Sample an interval, given the interval's bounding x-coordinates,.
Definition: unintegrable.C:78
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:71
Quadratic equation of the form a*x^2 + b*x + c = 0.
Definition: quadraticEqn.H:52
Roots< 2 > roots() const
Get the roots.
Definition: quadraticEqn.C:31
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
volScalarField & b
Definition: createFields.H:27
const dimensionedScalar c
Speed of light in a vacuum.
label log2(label i)
Return the log base 2 by successive bit-shifting of the given label.
Definition: label.H:79
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar y0(const dimensionedScalar &ds)
messageStream Info
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:221
labelList f(nPoints)
dictionary dict
Random rndGen(label(0))