44 for (
label i = 1; i <
x.size(); ++ i)
46 result[i] = result[i - 1] + (
x[i] -
x[i - 1])*(
y[i] +
y[i - 1])/2;
64 for (
label i = 1; i <
x.size(); ++ i)
66 const scalar x0 =
x[i - 1], x1 =
x[i];
67 const scalar
y0 =
y[i - 1],
y1 =
y[i];
70 result[i - 1] + (x1 - x0)*(2*x0*
y0 + x0*
y1 + x1*
y0 + 2*x1*
y1)/6;
85 scalar
f = (
s - Phi[0])/(Phi[1] - Phi[0]);
88 return (1 -
f)*
x[0] +
f*
x[1];
101 scalar
f = (
s - Phi[0])/(Phi[1] - Phi[0]);
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;
123 return (1 -
f)*
x[0] +
f*
x[1];
136 for (; i <
x.size() && Phi[i + 1] <
s; ++ i);
143 {Phi[i], Phi[i + 1]},
159 for (; i <
x.size() && Phi[i + 1] <
s; ++ i);
166 {phi[i], phi[i + 1]},
167 {Phi[i], Phi[i + 1]},
177 if (xPtr_.valid())
return xPtr_();
184 const scalar
f = scalar(i)/(n_ + 1);
185 x[i] = (1 -
f)*this->
min() + f*this->
max();
195 const scalar dPhi = (Phi.
last() - Phi.
first())/(n_ - 1);
197 if (distribution::debug)
199 scalar
error = -vGreat;
200 for (
label i = 0; i < n_ - 1; ++ i)
211 x.first() = this->
min();
213 for (
label i0 = 0, i = 1; i0 < n_ - 1; ++ i0)
215 while (Phi[i0 + 1] > i*dPhi)
220 {xPrev[i0], xPrev[i0 + 1]},
221 {phi[i0], phi[i0 + 1]},
222 {Phi[i0], Phi[i0 + 1]},
226 x[i] = (
x[i] + xNext)/2;
232 x.last() = this->
max();
241 if (PDFPtr_.valid())
return PDFPtr_();
243 const Pair<scalar> Phi01 = this->Phi01();
245 PDFPtr_.set((phi(this->q(), xPtr_)/(Phi01[1] - Phi01[0])).ptr());
259 return integrate(
x, phi(q,
x));
295 if (Phi01Ptr_.valid())
return Phi01Ptr_();
315 n_((1<<
dict.lookupOrDefault<
label>(
"level", 16)) + 1)
353 const scalar
s = this->rndGen_.template sample01<scalar>();
355 const scalar dCDF = scalar(1)/(n_ - 1);
356 const label samplei = floor(
s/dCDF);
360 {
x()[samplei],
x()[samplei + 1]},
361 {PDF()[samplei], PDF()[samplei + 1]},
362 {samplei*dCDF, (samplei + 1)*dCDF},
373 return (Mu01[1] - Mu01[0])/(Phi01[1] - Phi01[0]);
386 label n = n_ - 1, level = 0;
387 while (
n >>= 1) ++ level;
400 return this->clipPDF(
x, phi/(Phi01[1] - Phi01[0]));
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
T & last()
Return the last element of the list.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
const Type & first() const
Return first.
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
void type(const direction i, const rootType t)
Set the type of the i-th root.
T & first()
Return the first element of the list.
T & last()
Return the last element of the list.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Base class for statistical distributions.
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
virtual Pair< scalar > Phi01(const label q) const
Access cached values of the un-normalised CDF at the minimum and.
Base class for distributions that do not have a closed integral form for the cumulative density funct...
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
scalar sample() const
Sample the distribution.
static tmp< scalarField > integrate(const scalarField &x, const scalarField &y)
Integrate the values y with respect to the coordinates x.
unintegrable(const word &name, const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
virtual ~unintegrable()
Destructor.
static tmp< scalarField > integrateX(const scalarField &x, const scalarField &y)
Integrate the values x*y with respect to the coordinates x.
const Pair< scalar > & Phi01() const
Access cached values of the un-normalised CDF at the minimum and.
virtual scalar mean() const
Return the mean value.
static scalar sampleInterval(const Pair< scalar > &x, const Pair< scalar > &Phi, const scalar s)
Sample an interval, given the interval's bounding x-coordinates,.
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Quadratic equation of the form a*x^2 + b*x + c = 0.
Roots< 2 > roots() const
Get the roots.
A class for managing temporary objects.
T & ref() const
Return non-const reference or generate a fatal error.
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
A class for handling words, derived from string.
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))
const dimensionedScalar c
Speed of light in a vacuum.
void writeEntryIfDifferent(Ostream &os, const word &entryName, const EntryType &value1, const EntryType &value2)
Helper function to write the keyword and entry only if the.
label log2(label i)
Return the log base 2 by successive bit-shifting of the given 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.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
dimensionedScalar y0(const dimensionedScalar &ds)
List< scalar > scalarList
A List of scalars.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
dimensionedScalar y1(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Ostream & indent(Ostream &os)
Indent stream.
randomGenerator rndGen(653213)