uniform.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-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 \*---------------------------------------------------------------------------*/
25 
26 #include "uniform.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace distributions
34 {
37 }
38 }
39 
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 Foam::scalar Foam::distributions::uniform::Phi(const scalar x) const
44 {
45  if (q() == -1)
46  {
47  return log(x);
48  }
49  else
50  {
51  return integerPow(x, 1 + q())/(1 + q());
52  }
53 }
54 
55 
56 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
57 
59 (
60  const dictionary& dict,
61  Random& rndGen,
62  const label sampleQ
63 )
64 :
65  FieldDistribution<distribution, uniform>(typeName, dict, rndGen, sampleQ),
66  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"})),
67  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"})),
68  Phi0_(Phi(min_)),
69  Phi1_(Phi(max_))
70 {
71  validateBounds(dict);
72  if (q() != 0) validatePositive(dict);
73  report();
74 }
75 
76 
78 :
80  min_(d.min_),
81  max_(d.max_),
82  Phi0_(Phi(min_)),
83  Phi1_(Phi(max_))
84 {}
85 
86 
87 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
88 
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  const scalar s = rndGen_.sample01<scalar>();
98 
99  if (q() == -1)
100  {
101  return min_*pow(max_/min_, s);
102  }
103  else
104  {
105  const scalar PhiS = (1 - s)*Phi0_ + s*Phi1_;
106 
107  return integerRoot((1 + q())*PhiS, 1 + q());
108  }
109 }
110 
111 
113 {
114  return min_;
115 }
116 
117 
119 {
120  return max_;
121 }
122 
123 
125 {
126  if (q() == -2)
127  {
128  return Foam::log(max_/min_)/(Phi1_ - Phi0_);
129  }
130  else
131  {
132  const scalar Mu0 = integerPow(min_, 2 + q())/(2 + q());
133  const scalar Mu1 = integerPow(max_, 2 + q())/(2 + q());
134 
135  return (Mu1 - Mu0)/(Phi1_ - Phi0_);
136  }
137 }
138 
139 
142 {
143  if (q() == -1)
144  {
145  return clipPDF(x, 1/x/Foam::log(max_/min_));
146  }
147  else
148  {
149  return clipPDF(x, integerPow(x, q())/(Phi1_ - Phi0_));
150  }
151 }
152 
153 
154 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Random number generator.
Definition: Random.H:58
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
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
Distribution in which all values between a specified minimum and maximum have the same probability.
Definition: uniform.H:68
virtual scalar min() const
Return the minimum value.
Definition: uniform.C:112
virtual scalar sample() const
Sample the distribution.
Definition: uniform.C:95
virtual tmp< scalarField > PDF(const scalarField &x) const
Return the distribution probability density function.
Definition: uniform.C:141
uniform(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: uniform.C:59
virtual ~uniform()
Destructor.
Definition: uniform.C:89
virtual scalar max() const
Return the maximum value.
Definition: uniform.C:118
virtual scalar mean() const
Return the mean value.
Definition: uniform.C:124
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))
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
Namespace for OpenFOAM.
scalar integerRoot(const scalar x, const label e)
Compute the power of the number x to the reciprocal integer 1/e.
Definition: scalarI.H:55
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
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dictionary dict
Random rndGen(label(0))