exponential.H
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 Class
25  Foam::distributions::exponential
26 
27 Description
28  Exponential distribution, scaled so that it spans between a specified
29  minimum and maximum value, rather than from zero to infinity
30 
31  \f[
32  PDF(x) = \lambda \exp(- \lambda x)
33  \f]
34 
35 Usage
36  Example usage:
37  \verbatim
38  {
39  type exponential;
40  Q 3;
41  min 0.01;
42  max 0.5;
43  lambda 3;
44  }
45  \endverbatim
46 
47 SourceFiles
48  exponential.C
49 
50 See also
51  Foam::distribution
52 
53 \*---------------------------------------------------------------------------*/
54 
55 #ifndef exponential_H
56 #define exponential_H
57 
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 namespace Foam
63 {
64 namespace distributions
65 {
66 
67 /*---------------------------------------------------------------------------*\
68  Class exponential Declaration
69 \*---------------------------------------------------------------------------*/
70 
71 class exponential
72 :
73  public FieldDistribution<unintegrableForNonZeroQ, exponential>
74 {
75  // Private Data
76 
77  //- Minimum value
78  const scalar min_;
79 
80  //- Maximum value
81  const scalar max_;
82 
83  //- Parameter
84  const scalar lambda_;
85 
86 
87  // Private Member Functions
88 
89  //- Return values of the un-normalised PDF for the given size exponent
90  // and x-coordinates.
91  virtual tmp<scalarField> phi
92  (
93  const label q,
94  const scalarField& x
95  ) const;
96 
97  //- Return values of the un-normalised CDF for zero effective size
98  // exponent and given x-coordinates
99  virtual tmp<scalarField> PhiForZeroQ(const scalarField& x) const;
100 
101 
102 public:
103 
104  //- Runtime type information
105  TypeName("exponential");
106 
107 
108  // Constructors
109 
110  //- Construct from a dictionary
112  (
113  const unitConversion& units,
114  const dictionary& dict,
115  const label sampleQ,
117  );
118 
119  //- Construct copy
120  exponential(const exponential& d, const label sampleQ);
121 
122  //- Construct and return a clone
123  virtual autoPtr<distribution> clone(const label sampleQ) const
124  {
125  return autoPtr<distribution>(new exponential(*this, sampleQ));
126  }
127 
128 
129  //- Destructor
130  virtual ~exponential();
131 
132 
133  // Member Functions
134 
135  //- Sample the distribution for zero effective size exponent
136  virtual scalar sampleForZeroQ() const;
137 
138  //- Sample the distribution
140 
141  //- Return the minimum value
142  virtual scalar min() const;
143 
144  //- Return the maximum value
145  virtual scalar max() const;
146 
147  //- Write to a stream
148  virtual void write(Ostream& os, const unitConversion& units) const;
149 
150  //- Return coordinates to plot across the range of the distribution
151  virtual tmp<scalarField> plotX(const label n) const;
152 };
153 
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
157 } // End namespace distributions
158 } // End namespace Foam
159 
160 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 
162 #endif
163 
164 // ************************************************************************* //
label n
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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:65
virtual autoPtr< distribution > clone(const label sampleQ) const
Construct and return a clone.
Definition: exponential.H:122
virtual scalar min() const
Return the minimum value.
Definition: exponential.C:120
TypeName("exponential")
Runtime type information.
virtual void write(Ostream &os, const unitConversion &units) const
Write to a stream.
Definition: exponential.C:133
virtual scalar sampleForZeroQ() const
Sample the distribution for zero effective size exponent.
Definition: exponential.C:111
virtual scalar max() const
Return the maximum value.
Definition: exponential.C:126
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: exponential.C:147
virtual ~exponential()
Destructor.
Definition: exponential.C:105
Random number generator.
A class for managing temporary objects.
Definition: tmp.H:55
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
Namespace for OpenFOAM.
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
const HashTable< unitConversion > & units()
Get the table of unit conversions.
dictionary dict
randomGenerator rndGen(653213)