unintegrable.H
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-2024 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 Class
25  Foam::distributions::unintegrable
26 
27 Description
28  Base class for distributions that do not have a closed integral form for
29  the cumulative density function (CDF) for some or all effective size
30  exponents.
31 
32 SourceFiles
33  unintegrable.C
34 
35 See also
36  Foam::distribution
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef unintegrable_H
41 #define unintegrable_H
42 
43 #include "distribution.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace distributions
50 {
51 
52 /*---------------------------------------------------------------------------*\
53  Class unintegrable Declaration
54 \*---------------------------------------------------------------------------*/
55 
56 class unintegrable
57 :
58  public distribution
59 {
60 private:
61 
62  // Private Data
63 
64  //- Number of intervals to use over the range of the distribution
65  const label n_;
66 
67  //- X-coordinates. These are spaced to contain equal increments of the
68  // CDF. Sampling is a inversion of the CDF; i.e., solve CDF(x) = s for
69  // x, where s is a random sample between 0 and 1. Uniform increments
70  // of CDF make this inversion trivial.
71  mutable autoPtr<scalarField> xPtr_;
72 
73  //- Values of Phi (i.e., the un-normalised CDF) at the minimum and
74  // maximum x-coordinates
75  mutable autoPtr<Pair<scalar>> Phi01Ptr_;
76 
77  //- Values of the PDF at the X-coordinates
78  mutable autoPtr<scalarField> PDFPtr_;
79 
80 
81  // Private Member Functions
82 
83  //- X-coordinates
84  const scalarField& x() const;
85 
86  //- Values of the PDF at the X-coordinates
87  const scalarField& PDF() const;
88 
89 
90 protected:
91 
92  // Protected Member Functions
93 
94  //- Return values of the un-normalised PDF for the given size exponent
95  // and x-coordinates. Must be provided by derivations.
96  virtual tmp<scalarField> phi
97  (
98  const label q,
99  const scalarField& x
100  ) const = 0;
101 
102  //- Return values of the un-normalised CDF for the given size exponent
103  // and x-coordinates. Can be overloaded by derivations, for example, if
104  // there is a simpler analytic solution for certain effective size
105  // exponents (probably zero).
106  virtual tmp<scalarField> Phi
107  (
108  const label q,
109  const scalarField& x
110  ) const;
111 
112  //- Return values of the un-normalised CDF at the minimum and maximum
113  // x-coordinates for the given size exponent. Can be overloaded by
114  // derivations, for example, if there is a simpler analytic solution
115  // for certain effective size exponents (probably zero).
116  virtual Pair<scalar> Phi01(const label q) const;
117 
118  //- Access cached values of the un-normalised CDF at the minimum and
119  // maximum x-coordinates.
120  const Pair<scalar>& Phi01() const;
121 
122 
123 public:
124 
125  // Static Member Functions
126 
127  //- Integrate the values y with respect to the coordinates x
129  (
130  const scalarField& x,
131  const scalarField& y
132  );
133 
134  //- Integrate the values x*y with respect to the coordinates x
136  (
137  const scalarField& x,
138  const scalarField& y
139  );
140 
141  //- Integrate the values x^e*y with respect to the coordinates x,
142  // and interpolating onto a separate set of x coordinates
144  (
145  const scalarField& xStar,
146  const label e,
147  const scalarField& yStar,
148  const scalarField& x
149  );
150 
151  //- Sample an interval, given the interval's bounding x-coordinates,
152  // values Phi (the un-normalised CDF), and a random sample. First
153  // order interpolation.
154  static scalar sampleInterval
155  (
156  const Pair<scalar>& x,
157  const Pair<scalar>& Phi,
158  const scalar s
159  );
160 
161  //- Sample an interval, given the interval's bounding x-coordinates,
162  // values phi and Phi (the un-normalised PDF and CDF), and a random
163  // sample. Second order interpolation.
164  static scalar sampleInterval
165  (
166  const Pair<scalar>& x,
167  const Pair<scalar>& phi,
168  const Pair<scalar>& Phi,
169  const scalar s
170  );
171 
172  //- Sample a discretised distribution, given the x-coordinates,
173  // values Phi (the un-normalised CDF), and a random sample. First
174  // order interpolation.
175  static scalar sample
176  (
177  const scalarField& x,
178  const scalarField& Phi,
179  const scalar s
180  );
181 
182  //- Sample a discretised distribution, given the x-coordinates, values
183  // phi and Phi (the un-normalised PDF and CDF), and a random sample.
184  // Second order interpolation.
185  static scalar sample
186  (
187  const scalarField& x,
188  const scalarField& phi,
189  const scalarField& Phi,
190  const scalar s
191  );
192 
193 
194  // Constructors
195 
196  //- Construct from a dictionary
198  (
199  const word& name,
200  const unitConversion& units,
201  const dictionary& dict,
202  const label sampleQ,
204  );
205 
206  //- Construct from components
208  (
209  const label Q,
210  const label sampleQ,
212  const label n
213  );
214 
215  //- Construct copy
216  unintegrable(const unintegrable& d, const label sampleQ);
217 
218 
219  //- Destructor
220  virtual ~unintegrable();
221 
222 
223  // Member Functions
224 
225  //- Sample the distribution
226  virtual scalar sample() const;
227 
228  //- Sample the distribution
229  using distribution::sample;
230 
231  //- Return the mean value
232  virtual scalar mean() const;
233 
234  //- Return the integral of the PDF multiplied by an integer power of x
236  (
237  const scalarField& x,
238  const label e,
239  const bool consistent = false
240  ) const;
241 
242  //- Write to a stream
243  virtual void write(Ostream& os, const unitConversion& units) const;
244 
245  //- Return coordinates to plot across the range of the distribution
246  using distribution::plotX;
247 
248  //- Return values to plot the probability density function
249  virtual tmp<scalarField> plotPDF(const scalarField& x) const;
250 };
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace distributions
256 } // End namespace Foam
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 #endif
261 
262 // ************************************************************************* //
scalar y
label n
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
label sampleQ() const
Access the sample size exponent.
Definition: distribution.C:134
virtual scalar sample() const =0
Sample the distribution.
label q() const
Return the effective distribution size exponent.
Definition: distribution.H:107
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:166
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:320
virtual tmp< scalarField > integralPDFxPow(const scalarField &x, const label e, const bool consistent=false) const
Return the integral of the PDF multiplied by an integer power of x.
Definition: unintegrable.C:445
virtual tmp< scalarField > phi(const label q, const scalarField &x) const =0
Return values of the un-normalised PDF for the given size exponent.
virtual scalar sample() const
Sample the distribution.
Definition: unintegrable.C:416
static tmp< scalarField > integrate(const scalarField &x, const scalarField &y)
Integrate the values y with respect to the coordinates x.
Definition: unintegrable.C:34
unintegrable(const word &name, const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
Definition: unintegrable.C:353
virtual tmp< scalarField > plotPDF(const scalarField &x) const
Return values to plot the probability density function.
Definition: unintegrable.C:472
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
Definition: unintegrable.C:456
virtual ~unintegrable()
Destructor.
Definition: unintegrable.C:410
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:340
virtual scalar mean() const
Return the mean value.
Definition: unintegrable.C:434
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:144
static tmp< scalarField > interpolateIntegrateXPow(const scalarField &xStar, const label e, const scalarField &yStar, const scalarField &x)
Integrate the values x^e*y with respect to the coordinates x,.
Definition: unintegrable.C:79
Random number generator.
A class for managing temporary objects.
Definition: tmp.H:55
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
A class for handling words, derived from string.
Definition: word.H:62
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
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
const HashTable< unitConversion > & units()
Get the table of unit conversions.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict
randomGenerator rndGen(653213)