exponential.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 "exponential.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::tmp<Foam::scalarField> Foam::distributions::exponential::phi
44 (
45  const label q,
46  const scalarField& x
47 ) const
48 {
49  return integerPow(x, q)*lambda_*exp(- lambda_*x);
50 }
51 
52 
53 Foam::tmp<Foam::scalarField> Foam::distributions::exponential::Phi
54 (
55  const label q,
56  const scalarField& x
57 ) const
58 {
59  if (q == 0)
60  {
61  return - exp(- lambda_*x);
62  }
63  else
64  {
66  }
67 }
68 
69 
70 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71 
73 (
74  const dictionary& dict,
75  Random& rndGen,
76  const label sampleQ
77 )
78 :
80  (
81  typeName,
82  dict,
83  rndGen,
84  sampleQ
85  ),
86  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"})),
87  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"})),
88  lambda_(dict.lookup<scalar>("lambda"))
89 {
90  validateBounds(dict);
91  validatePositive(dict);
92  mean();
93  report();
94 }
95 
96 
98 (
99  const exponential& d,
100  const label sampleQ
101 )
102 :
104  min_(d.min_),
105  max_(d.max_),
106  lambda_(d.lambda_)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
111 
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  if (q() == 0)
121  {
122  const scalar s = rndGen_.sample01<scalar>();
123  const Pair<scalar>& Phi01 = this->Phi01();
124  const scalar PhiS = (1 - s)*Phi01[0] + s*Phi01[1];
125  return - 1/lambda_*log(- PhiS);
126  }
127  else
128  {
130  }
131 }
132 
133 
135 {
136  return min_;
137 }
138 
139 
141 {
142  return max_;
143 }
144 
145 
148 {
150  tx.ref()[0] = Foam::max(tx.ref()[0], q() < 0 ? min_/2 : rootVSmall);
151  return tx;
152 }
153 
154 
155 // ************************************************************************* //
label n
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
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: distribution.C:140
Exponential distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: exponential.H:73
exponential(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: exponential.C:73
virtual scalar min() const
Return the minimum value.
Definition: exponential.C:134
virtual scalar sample() const
Sample the distribution.
Definition: exponential.C:118
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: exponential.C:147
virtual scalar max() const
Return the maximum value.
Definition: exponential.C:140
virtual ~exponential()
Destructor.
Definition: exponential.C:112
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
Definition: unintegrable.C:253
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:349
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
volScalarField scalarField(fieldObject, mesh)
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.
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.
dimensionedScalar log(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dictionary dict
Random rndGen(label(0))