Random.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) 2018-2020 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::Random
26 
27 Description
28  Random number generator
29 
30  This is a clone of the drand48 algorithm. This is significantly quicker
31  than drand48, presumably due to the compiler inlining the sampling methods.
32  It is also significantly quicker than the standard library linear
33  congruential engine, as it does not use Schrage's algorithm to prevent
34  overflow.
35 
36  See <http://pubs.opengroup.org/onlinepubs/007908775/xsh/drand48.html> for
37  details of the seeding and iteration sequence.
38 
39 SourceFiles
40  RandomI.H
41 
42 \*---------------------------------------------------------------------------*/
43 
44 #ifndef Random_H
45 #define Random_H
46 
47 #include "scalar.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 
54 /*---------------------------------------------------------------------------*\
55  Class Random Declaration
56 \*---------------------------------------------------------------------------*/
57 
58 class Random
59 {
60  // Private Typedefs
61 
62  //- Working type
63  typedef uint64_t type;
64 
65 
66  // Private static data
67 
68  //- The parameters of the linear congruential iteration
69  static const type A = 0x5DEECE66D, C = 0xB, M = type(1) << 48;
70 
71 
72  // Private Data
73 
74  //- The stored integer
75  type x_;
76 
77  //- Is a normal scalar sample stored?
78  bool scalarNormalStored_;
79 
80  //- A stored normal scalar sample
81  scalar scalarNormalValue_;
82 
83 
84  // Private Member Functions
85 
86  //- Advance the state and return an integer sample
87  inline type sample();
88 
89 
90 public:
91 
92  // Constructors
93 
94  //- Construct from a seed
95  inline Random(const label s);
96 
97 
98  //- Destructor
99  inline ~Random();
100 
101 
102  // Member Functions
103 
104  // Scalars
105 
106  //- Advance the state and return a scalar sample from a uniform
107  // distribution between zero and one
108  inline scalar scalar01();
109 
110  //- Advance the state and return a scalar sample from a uniform
111  // distribution between two limits
112  inline scalar scalarAB(const scalar a, const scalar b);
113 
114  //- Advance the state and return a scalar sample from a normal
115  // distribution with mean zero and standard deviation one
116  scalar scalarNormal();
117 
118 
119  // Other types
120 
121  //- Advance the state and return a sample of a given type from a
122  // uniform distribution between zero and one
123  template<class Type>
124  inline Type sample01();
125 
126  //- Advance the state and return a sample of a given type from a
127  // uniform distribution between two limits
128  template<class Type>
129  inline Type sampleAB(const Type& a, const Type& b);
130 
131  //- Advance the state and return a sample of a given type from a
132  // normal distribution with mean zero and standard deviation one
133  template<class Type>
134  inline Type sampleNormal();
135 
136 
137  // Global scalars
138 
139  //- Advance the state and return a scalar sample from a uniform
140  // distribution between zero and one. Synchronises across all
141  // cores. Use of this is discouraged. It is expensive and
142  // introduces non-randomness in all cores other then the master.
143  scalar globalScalar01();
144 };
145 
146 
147 template<>
148 inline scalar Random::sample01();
149 
150 template<>
151 inline label Random::sample01();
152 
153 template<>
154 inline scalar Random::sampleAB(const scalar& a, const scalar& b);
155 
156 template<>
157 inline label Random::sampleAB(const label& a, const label& b);
158 
159 template<>
160 inline scalar Random::sampleNormal();
161 
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace Foam
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 #include "RandomI.H"
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 #endif
174 
175 // ************************************************************************* //
Graphite solid properties.
Definition: C.H:48
scalar scalarNormal()
Advance the state and return a scalar sample from a normal.
Definition: Random.C:31
Type sampleNormal()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:119
scalar scalarAB(const scalar a, const scalar b)
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:63
Type sampleAB(const Type &a, const Type &b)
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:98
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:70
Random(const label s)
Construct from a seed.
Definition: RandomI.H:41
scalar globalScalar01()
Advance the state and return a scalar sample from a uniform.
Definition: Random.C:94
Random number generator.
Definition: Random.H:57
~Random()
Destructor.
Definition: RandomI.H:51
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:57
Namespace for OpenFOAM.