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-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::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}} \exp \\
33  \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 
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 
92  // Private Member Functions
93 
94  //- Return values of the un-normalised PDF for the given size exponent
95  // and x-coordinates.
96  virtual tmp<scalarField> phi
97  (
98  const label q,
99  const scalarField& x
100  ) const;
101 
102  //- Return values of the un-normalised CDF for zero effective size
103  // exponent and given x-coordinates
104  virtual tmp<scalarField> PhiForZeroQ(const scalarField& x) const;
105 
106 
107 public:
108 
109  //- Runtime type information
110  TypeName("normal");
111 
112 
113  //- Permit the multiNormal distribution to use private parts of this class
114  friend class multiNormal;
115 
116 
117  // Constructors
118 
119  //- Construct from a dictionary
120  normal
121  (
122  const unitConversion& units,
123  const dictionary& dict,
124  const label sampleQ,
126  );
127 
128  //- Construct from components
129  normal
130  (
131  const label Q,
132  const label sampleQ,
134  const label n,
135  const scalar min,
136  const scalar max,
137  const scalar mu,
138  const scalar sigma
139  );
140 
141  //- Construct copy
142  normal(const normal& d, const label sampleQ);
143 
144  //- Construct and return a clone
145  virtual autoPtr<distribution> clone(const label sampleQ) const
146  {
147  return autoPtr<distribution>(new normal(*this, sampleQ));
148  }
149 
150 
151  //- Destructor
152  virtual ~normal();
153 
154 
155  // Member Functions
156 
157  //- Sample the distribution for zero effective size exponent
158  virtual scalar sampleForZeroQ() const;
159 
160  //- Sample the distribution for zero effective size exponent
161  scalar sampleForZeroQ(const scalar s) const;
162 
163  //- Sample the distribution
165 
166  //- Return the minimum value
167  virtual scalar min() const;
168 
169  //- Return the maximum value
170  virtual scalar max() const;
171 
172  //- Return the mean value
173  virtual scalar mean() const;
174 
175  //- Return the mean value
176  inline scalar mu() const
177  {
178  return mu_;
179  }
180 
181  //- Return the standard deviation
182  inline scalar sigma() const
183  {
184  return sigma_;
185  }
186 
187  //- Write to a stream
188  virtual void write(Ostream& os, const unitConversion& units) const;
189 
190  //- Return coordinates to plot across the range of the distribution
191  virtual tmp<scalarField> plotX(const label n) const;
192 };
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 } // End namespace distributions
198 } // End namespace Foam
199 
200 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
201 
202 #endif
203 
204 // ************************************************************************* //
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
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:143
scalar sigma() const
Return the standard deviation.
Definition: normal.H:180
virtual scalar min() const
Return the minimum value.
Definition: normal.C:157
normal(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
Definition: normal.C:72
virtual ~normal()
Destructor.
Definition: normal.C:136
scalar mu() const
Return the mean value.
Definition: normal.H:174
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
Definition: normal.C:191
virtual scalar sampleForZeroQ() const
Sample the distribution for zero effective size exponent.
Definition: normal.C:142
TypeName("normal")
Runtime type information.
virtual scalar max() const
Return the maximum value.
Definition: normal.C:163
virtual scalar mean() const
Return the mean value.
Definition: normal.C:169
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: normal.C:206
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...
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.
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.
dictionary dict
randomGenerator rndGen(653213)