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::PhiForZeroQ
60 (
61  const scalarField& x
62 ) const
63 {
64  static const scalar sqrt2 = sqrt(scalar(2));
65  return (1 + standardNormal::approxErf((x - mu_)/(sigma_*sqrt2)))/2;
66 }
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
72 (
73  const unitConversion& units,
74  const dictionary& dict,
75  const label sampleQ,
77 )
78 :
80  (
81  typeName,
82  units,
83  dict,
84  sampleQ,
85  std::move(rndGen)
86  ),
87  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"}, units)),
88  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"}, units)),
89  mu_(dict.lookupBackwardsCompatible<scalar>({"mu", "expectation"}, units)),
90  sigma_(dict.lookup<scalar>("sigma", units))
91 {
92  validateBounds(dict);
93  if (q() != 0) validatePositive(dict);
94  mean();
95 }
96 
97 
99 (
100  const label Q,
101  const label sampleQ,
103  const label n,
104  const scalar min,
105  const scalar max,
106  const scalar mu,
107  const scalar sigma
108 )
109 :
111  (
112  Q,
113  sampleQ,
114  std::move(rndGen),
115  n
116  ),
117  min_(min),
118  max_(max),
119  mu_(mu),
120  sigma_(sigma)
121 {}
122 
123 
125 :
127  min_(d.min_),
128  max_(d.max_),
129  mu_(d.mu_),
130  sigma_(d.sigma_)
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  return sampleForZeroQ(rndGen_.sample01<scalar>());
145 }
146 
147 
148 Foam::scalar Foam::distributions::normal::sampleForZeroQ(const scalar s) const
149 {
150  static const scalar sqrt2 = sqrt(scalar(2));
151  const Pair<scalar>& Phi01 = this->Phi01();
152  const scalar PhiS = (1 - s)*Phi01[0] + s*Phi01[1];
153  return standardNormal::approxErfInv(2*PhiS - 1)*sigma_*sqrt2 + mu_;
154 }
155 
156 
158 {
159  return min_;
160 }
161 
162 
164 {
165  return max_;
166 }
167 
168 
170 {
171  if (q() == 0)
172  {
173  // !!! This isn't technically correct. The min/max clipping affects the
174  // mean. To calculate it properly we have to integrate the first moment
175  // of the PDF, and that can't be done analytically; we have to hand
176  // over to the unintegrable handling. This probably isn't worth it for
177  // what should be a small deviation. Also, it looks a bit odd in the
178  // log for the distribution to be reporting a different mean to that
179  // which was specified. So, just use the mean value of the un-clipped
180  // distribution for now.
181  return mu_;
182  }
183  else
184  {
186  }
187 }
188 
189 
191 (
192  Ostream& os,
193  const unitConversion& units
194 ) const
195 {
197 
198  writeEntry(os, "min", units, min_);
199  writeEntry(os, "max", units, max_);
200  writeEntry(os, "mu", units, mu_);
201  writeEntry(os, "sigma", units, sigma_);
202 }
203 
204 
207 {
209  tx.ref()[0] = Foam::max(tx.ref()[0], q() < 0 ? min_/2 : -vGreat);
210  return tx;
211 }
212 
213 
214 // ************************************************************************* //
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 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
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:166
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: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
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
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
static scalar approxErfInv(const scalar y)
Approximate error function inverse.
Base class for distributions that have a closed integral form for the cumulative density function (CD...
virtual scalar mean() const
Return the mean value.
Definition: unintegrable.C:434
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:197
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))
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
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict
randomGenerator rndGen(653213)