normal.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) 2011-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::normal
26 
27 Description
28  Normal distribution, scaled so that it spans between a specified
29  minimum and maximum value, rather than from zero to infinity
30 
31  \f[
32  PDF(x) = \frac{1}{\sigma \sqrt{2 \pi}} \\
33  \exp \left( \frac{1}{2} \left( \frac{x - \mu}{\sigma} \right)^2 \right)
34  \f]
35 
36 Usage
37  Example usage:
38  \verbatim
39  {
40  type normal;
41  Q 0;
42  min 0.001;
43  max 0.019;
44  mu 0.011;
45  sigma 0.003;
46  }
47  \endverbatim
48 
49 SourceFiles
50  normal.C
51 
52 See also
53  Foam::distribution
54 
55 \*---------------------------------------------------------------------------*/
56 
57 #ifndef normal_H
58 #define normal_H
59 
60 #include "unintegrable.H"
61 
62 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64 namespace Foam
65 {
66 namespace distributions
67 {
68 
69 /*---------------------------------------------------------------------------*\
70  Class normal Declaration
71 \*---------------------------------------------------------------------------*/
72 
73 class normal
74 :
75  public FieldDistribution<unintegrableForNonZeroQ, normal>
76 {
77  // Private Data
78 
79  //- Minimum value
80  const scalar min_;
81 
82  //- Maximum value
83  const scalar max_;
84 
85  //- Mean
86  const scalar mu_;
87 
88  //- Standard deviation
89  const scalar sigma_;
90 
91  //- Constant for approximate error function and inverse error function
92  static const scalar a_;
93 
94 
95  // Private Member Functions
96 
97  //- Approximate error function
98  static tmp<scalarField> approxErf(const scalarField& x);
99 
100  //- Approximate error function inverse
101  static scalar approxErfInv(const scalar y);
102 
103  //- Return values of the un-normalised PDF for the given size exponent
104  // and x-coordinates.
105  virtual tmp<scalarField> phi
106  (
107  const label q,
108  const scalarField& x
109  ) const;
110 
111  //- Return values of the un-normalised CDF for the given size exponent
112  // and x-coordinates.
113  virtual tmp<scalarField> Phi
114  (
115  const label q,
116  const scalarField& x
117  ) const;
118 
119  //- Sample the distribution for zero effective size exponent
120  scalar sampleForZeroQ(const scalar s) const;
121 
122 
123 public:
124 
125  //- Runtime type information
126  TypeName("normal");
127 
128 
129  //- Permit the multiNormal distribution to use private parts of this class
130  friend class multiNormal;
131 
132 
133  // Constructors
134 
135  //- Construct from a dictionary
136  normal
137  (
138  const dictionary& dict,
139  Random& rndGen,
140  const label sampleQ
141  );
142 
143  //- Construct from components
144  normal
145  (
146  Random& rndGen,
147  const label Q,
148  const label sampleQ,
149  const label n,
150  const scalar min,
151  const scalar max,
152  const scalar mu,
153  const scalar sigma
154  );
155 
156  //- Construct copy
157  normal(const normal& d, const label sampleQ);
158 
159  //- Construct and return a clone
160  virtual autoPtr<distribution> clone(const label sampleQ) const
161  {
162  return autoPtr<distribution>(new normal(*this, sampleQ));
163  }
164 
165 
166  //- Destructor
167  virtual ~normal();
168 
169 
170  // Member Functions
171 
172  //- Sample the distribution
173  virtual scalar sample() const;
174 
175  //- Sample the distribution
177 
178  //- Return the minimum value
179  virtual scalar min() const;
180 
181  //- Return the maximum value
182  virtual scalar max() const;
183 
184  //- Return the mean value
185  virtual scalar mean() const;
186 
187  //- Return coordinates to plot across the range of the distribution
188  virtual tmp<scalarField> x(const label n) const;
189 };
190 
191 
192 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 
194 } // End namespace distributions
195 } // End namespace Foam
196 
197 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
198 
199 #endif
200 
201 // ************************************************************************* //
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
Multiple superimposed normal distributions.
Definition: multiNormal.H:75
Normal distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: normal.H:74
virtual autoPtr< distribution > clone(const label sampleQ) const
Construct and return a clone.
Definition: normal.H:158
virtual scalar min() const
Return the minimum value.
Definition: normal.C:185
virtual scalar sample() const
Sample the distribution.
Definition: normal.C:172
normal(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: normal.C:109
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: normal.C:219
virtual ~normal()
Destructor.
Definition: normal.C:166
TypeName("normal")
Runtime type information.
virtual scalar max() const
Return the maximum value.
Definition: normal.C:191
virtual scalar mean() const
Return the mean value.
Definition: normal.C:197
A class for managing temporary objects.
Definition: tmp.H:55
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 sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
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
dictionary dict
Random rndGen(label(0))