normal.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "normal.H"
27 #include "standardNormal.H"
28 #include "mathematicalConstants.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace distributions
36 {
39 }
40 }
41 
42 
43 using namespace Foam::constant::mathematical;
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 Foam::tmp<Foam::scalarField> Foam::distributions::normal::phi
49 (
50  const label q,
51  const scalarField& x
52 ) const
53 {
54  static const scalar sqrt2Pi = sqrt(2*pi);
55  return integerPow(x, q)/(sigma_*sqrt2Pi)*exp(- sqr((x - mu_)/sigma_)/2);
56 }
57 
58 
59 Foam::tmp<Foam::scalarField> Foam::distributions::normal::Phi
60 (
61  const label q,
62  const scalarField& x
63 ) const
64 {
65  if (q == 0)
66  {
67  static const scalar sqrt2 = sqrt(scalar(2));
68  return (1 + standardNormal::approxErf((x - mu_)/(sigma_*sqrt2)))/2;
69  }
70  else
71  {
72  return unintegrableForNonZeroQ::Phi(q, x);
73  }
74 }
75 
76 
77 Foam::scalar Foam::distributions::normal::sampleForZeroQ(const scalar s) const
78 {
79  static const scalar sqrt2 = sqrt(scalar(2));
80  const Pair<scalar>& Phi01 = this->Phi01();
81  const scalar PhiS = (1 - s)*Phi01[0] + s*Phi01[1];
82  return standardNormal::approxErfInv(2*PhiS - 1)*sigma_*sqrt2 + mu_;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const unitConversion& units,
91  const dictionary& dict,
92  const label sampleQ,
94 )
95 :
97  (
98  typeName,
99  units,
100  dict,
101  sampleQ,
102  std::move(rndGen)
103  ),
104  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"}, units)),
105  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"}, units)),
106  mu_(dict.lookupBackwardsCompatible<scalar>({"mu", "expectation"}, units)),
107  sigma_(dict.lookup<scalar>("sigma", units))
108 {
109  validateBounds(dict);
110  if (q() != 0) validatePositive(dict);
111  mean();
112 }
113 
114 
116 (
117  const label Q,
118  const label sampleQ,
120  const label n,
121  const scalar min,
122  const scalar max,
123  const scalar mu,
124  const scalar sigma
125 )
126 :
128  (
129  Q,
130  sampleQ,
131  std::move(rndGen),
132  n
133  ),
134  min_(min),
135  max_(max),
136  mu_(mu),
137  sigma_(sigma)
138 {}
139 
140 
142 :
144  min_(d.min_),
145  max_(d.max_),
146  mu_(d.mu_),
147  sigma_(d.sigma_)
148 {}
149 
150 
151 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
152 
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 {
161  if (q() == 0)
162  {
163  return sampleForZeroQ(rndGen_.sample01<scalar>());
164  }
165  else
166  {
168  }
169 }
170 
171 
173 {
174  return min_;
175 }
176 
177 
179 {
180  return max_;
181 }
182 
183 
185 {
186  if (q() == 0)
187  {
188  // !!! This isn't technically correct. The min/max clipping affects the
189  // mean. To calculate it properly we have to integrate the first moment
190  // of the PDF, and that can't be done analytically; we have to hand
191  // over to the unintegrable handling. This probably isn't worth it for
192  // what should be a small deviation. Also, it looks a bit odd in the
193  // log for the distribution to be reporting a different mean to that
194  // which was specified. So, just use the mean value of the un-clipped
195  // distribution for now.
196  return mu_;
197  }
198  else
199  {
201  }
202 }
203 
204 
206 (
207  Ostream& os,
208  const unitConversion& units
209 ) const
210 {
212 
213  writeEntry(os, "min", units, min_);
214  writeEntry(os, "max", units, max_);
215  writeEntry(os, "mu", units, mu_);
216  writeEntry(os, "sigma", units, sigma_);
217 }
218 
219 
222 {
224  tx.ref()[0] = Foam::max(tx.ref()[0], q() < 0 ? min_/2 : -vGreat);
225  return tx;
226 }
227 
228 
229 // ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:147
Normal distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: normal.H:74
virtual scalar min() const
Return the minimum value.
Definition: normal.C:172
virtual scalar sample() const
Sample the distribution.
Definition: normal.C:159
normal(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
Definition: normal.C:89
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: normal.C:221
virtual ~normal()
Destructor.
Definition: normal.C:153
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
Definition: normal.C:206
virtual scalar max() const
Return the maximum value.
Definition: normal.C:178
virtual scalar mean() const
Return the mean value.
Definition: normal.C:184
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:351
virtual scalar mean() const
Return the mean value.
Definition: unintegrable.C:368
Random number generator.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
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))
mathematical constants.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
scalar integerPow(const scalar x, const label e)
Compute the power of the number x to the integer e.
Definition: scalarI.H:30
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.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
randomGenerator rndGen(653213)