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 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  //- Sample an interval, given the interval's bounding x-coordinates,
142  // values Phi (the un-normalised CDF), and a random sample. First
143  // order interpolation.
144  static scalar sampleInterval
145  (
146  const Pair<scalar>& x,
147  const Pair<scalar>& Phi,
148  const scalar s
149  );
150 
151  //- Sample an interval, given the interval's bounding x-coordinates,
152  // values phi and Phi (the un-normalised PDF and CDF), and a random
153  // sample. Second order interpolation.
154  static scalar sampleInterval
155  (
156  const Pair<scalar>& x,
157  const Pair<scalar>& phi,
158  const Pair<scalar>& Phi,
159  const scalar s
160  );
161 
162  //- Sample a discretised distribution, given the x-coordinates,
163  // values Phi (the un-normalised CDF), and a random sample. First
164  // order interpolation.
165  static scalar sample
166  (
167  const scalarField& x,
168  const scalarField& Phi,
169  const scalar s
170  );
171 
172  //- Sample a discretised distribution, given the x-coordinates, values
173  // phi and Phi (the un-normalised PDF and CDF), and a random sample.
174  // Second order interpolation.
175  static scalar sample
176  (
177  const scalarField& x,
178  const scalarField& phi,
179  const scalarField& Phi,
180  const scalar s
181  );
182 
183 
184  // Constructors
185 
186  //- Construct from a dictionary
188  (
189  const word& name,
190  const dictionary& dict,
191  Random& rndGen,
192  const label sampleQ
193  );
194 
195  //- Construct from components
197  (
198  Random& rndGen,
199  const label Q,
200  const label sampleQ,
201  const label n
202  );
203 
204  //- Construct copy
205  unintegrable(const unintegrable& d, const label sampleQ);
206 
207 
208  //- Destructor
209  virtual ~unintegrable();
210 
211 
212  // Member Functions
213 
214  //- Sample the distribution
215  scalar sample() const;
216 
217  //- Sample the distribution
218  using distribution::sample;
219 
220  //- Return the mean value
221  virtual scalar mean() const;
222 
223  //- Return coordinates to plot across the range of the distribution
224  using distribution::x;
225 
226  //- Return the distribution probability density function
227  virtual tmp<scalarField> PDF(const scalarField& x) const;
228 };
229 
230 
231 /*---------------------------------------------------------------------------*\
232  Class unintegrable Declaration
233 \*---------------------------------------------------------------------------*/
234 
236 :
237  public unintegrable
238 {
239 protected:
240 
241  // Protected Member Functions
242 
243  //- Return values of the un-normalised CDF at the minimum and maximum
244  // x-coordinates for the given size exponent.
245  virtual Pair<scalar> Phi01(const label q) const;
246 
247  //- Access cached values of the un-normalised CDF at the minimum and
248  // maximum x-coordinates.
249  using unintegrable::Phi01;
250 
251 
252 public:
253 
254  // Constructors
255 
256  //- Inherit constructors
258 };
259 
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 } // End namespace distributions
264 } // End namespace Foam
265 
266 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
267 
268 #endif
269 
270 // ************************************************************************* //
scalar y
label n
Random number generator.
Definition: Random.H:58
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 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 scalar sample() const =0
Sample the distribution.
label q() const
Return the effective distribution size exponent.
Definition: distribution.H:106
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:140
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
virtual tmp< scalarField > phi(const label q, const scalarField &x) const =0
Return values of the un-normalised PDF for the given size exponent.
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
A class for managing temporary objects.
Definition: tmp.H:55
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))
Namespace for OpenFOAM.
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dictionary dict
Random rndGen(label(0))