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-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 "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 unitConversion& units,
75  const dictionary& dict,
76  const label sampleQ,
78 )
79 :
81  (
82  typeName,
83  units,
84  dict,
85  sampleQ,
86  std::move(rndGen)
87  ),
88  min_(dict.lookupBackwardsCompatible<scalar>({"min", "minValue"}, units)),
89  max_(dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"}, units)),
90  lambda_(dict.lookup<scalar>("lambda", unitless))
91 {
92  validateBounds(dict);
93  validatePositive(dict);
94  mean();
95 }
96 
97 
99 (
100  const exponential& d,
101  const label sampleQ
102 )
103 :
105  min_(d.min_),
106  max_(d.max_),
107  lambda_(d.lambda_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
112 
114 {}
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 {
121  if (q() == 0)
122  {
123  const scalar s = rndGen_.sample01<scalar>();
124  const Pair<scalar>& Phi01 = this->Phi01();
125  const scalar PhiS = (1 - s)*Phi01[0] + s*Phi01[1];
126  return - 1/lambda_*log(- PhiS);
127  }
128  else
129  {
131  }
132 }
133 
134 
136 {
137  return min_;
138 }
139 
140 
142 {
143  return max_;
144 }
145 
146 
148 (
149  Ostream& os,
150  const unitConversion& units
151 ) const
152 {
154 
155  writeEntry(os, "min", units, min_);
156  writeEntry(os, "max", units, max_);
157  writeEntry(os, "lambda", unitless, lambda_);
158 }
159 
160 
163 {
165  tx.ref()[0] = Foam::max(tx.ref()[0], q() < 0 ? min_/2 : rootVSmall);
166  return tx;
167 }
168 
169 
170 // ************************************************************************* //
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
Exponential distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: exponential.H:73
exponential(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
Definition: exponential.C:73
virtual scalar min() const
Return the minimum value.
Definition: exponential.C:135
virtual scalar sample() const
Sample the distribution.
Definition: exponential.C:119
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: exponential.C:162
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
Definition: exponential.C:148
virtual scalar max() const
Return the maximum value.
Definition: exponential.C:141
virtual ~exponential()
Destructor.
Definition: exponential.C:113
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:254
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:351
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...
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)
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.
dimensionedScalar log(const dimensionedScalar &ds)
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)
const unitConversion unitless
dictionary dict
randomGenerator rndGen(653213)