randomGeneratorI.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-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 "randomGenerator.H"
27 #include "Hash.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 inline uint64_t Foam::randomGenerator::seed::x(const bool global) const
32 {
33  const uint64_t localS =
34  global ? s_ : s_ + (UINT64_MAX/Pstream::nProcs())*Pstream::myProcNo();
35 
36  return (localS << 16) + 0x330E;
37 }
38 
39 
40 inline void Foam::randomGenerator::checkSync() const
41 {
42  if (global_)
43  {
44  uint64_t xMaster = x_;
45  Pstream::scatter(xMaster);
46  if (xMaster != x_)
47  {
49  << "Global random number generator is not synchronised"
50  << exit(FatalError);
51  }
52  }
53 }
54 
55 
56 inline uint64_t Foam::randomGenerator::sample()
57 {
58  x_ = (A*x_ + C) % M;
59  return x_ >> 17;
60 }
61 
62 
63 inline Foam::scalar Foam::randomGenerator::scalar01NoCheckSync()
64 {
65  return scalar(sample())/(M >> 17);
66 }
67 
68 
69 inline Foam::scalar Foam::randomGenerator::scalarABNoCheckSync
70 (
71  const scalar a,
72  const scalar b
73 )
74 {
75  return a + scalar01NoCheckSync()*(b - a);
76 }
77 
78 
79 template<class Type>
80 inline Type Foam::randomGenerator::sample01NoCheckSync()
81 {
82  Type value;
83  for (direction i = 0; i < pTraits<Type>::nComponents; ++ i)
84  {
85  value.component(i) = scalar01NoCheckSync();
86  }
87  return value;
88 }
89 
90 
91 template<>
92 inline Foam::scalar Foam::randomGenerator::sample01NoCheckSync()
93 {
94  return scalar01NoCheckSync();
95 }
96 
97 
98 template<>
99 inline Foam::label Foam::randomGenerator::sample01NoCheckSync()
100 {
101  return sample() % 2;
102 }
103 
104 
105 template<class Type>
106 inline Type Foam::randomGenerator::sampleABNoCheckSync
107 (
108  const Type& a,
109  const Type& b
110 )
111 {
112  return a + cmptMultiply(sample01NoCheckSync<Type>(), b - a);
113 }
114 
115 
116 template<>
117 inline Foam::scalar Foam::randomGenerator::sampleABNoCheckSync
118 (
119  const scalar& a,
120  const scalar& b
121 )
122 {
123  return scalarABNoCheckSync(a, b);
124 }
125 
126 
127 template<>
128 inline Foam::label Foam::randomGenerator::sampleABNoCheckSync
129 (
130  const label& a,
131  const label& b
132 )
133 {
134  return a + sample() % (b - a);
135 }
136 
137 
138 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
139 
141 :
142  s_(s)
143 {}
144 
145 
147 :
148  s_(Hash<word>()(s))
149 {}
150 
151 
152 inline Foam::randomGenerator::randomGenerator(const seed s, const bool global)
153 :
154  global_(global),
155  x_(s.x(global))
156 {
157  checkSync();
158 }
159 
160 
162 :
163  global_(rndGen.global_),
164  x_(rndGen.x_)
165 {
166  checkSync();
167 }
168 
169 
171 :
172  randomGenerator(static_cast<const randomGenerator&>(rndGen))
173 {}
174 
175 
176 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
177 
179 {}
180 
181 
182 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
183 
184 inline Foam::scalar Foam::randomGenerator::scalar01()
185 {
186  #ifdef FULLDEBUG
187  checkSync();
188  #endif
189 
190  return scalar01NoCheckSync();
191 }
192 
193 
195 (
196  const label n
197 )
198 {
199  #ifdef FULLDEBUG
200  checkSync();
201  #endif
202 
203  tmp<scalarField> tValues(new scalarField(n));
204  for (label i = 0; i < n; ++ i)
205  {
206  tValues.ref()[i] = scalar01NoCheckSync();
207  }
208  return tValues;
209 }
210 
211 
213 (
214  const scalar a,
215  const scalar b
216 )
217 {
218  #ifdef FULLDEBUG
219  checkSync();
220  #endif
221 
222  return scalarABNoCheckSync(a, b);
223 }
224 
225 
227 (
228  const label n,
229  const scalar a,
230  const scalar b
231 )
232 {
233  #ifdef FULLDEBUG
234  checkSync();
235  #endif
236 
237  tmp<scalarField> tValues(new scalarField(n));
238  for (label i = 0; i < n; ++ i)
239  {
240  tValues.ref()[i] = scalarABNoCheckSync(a, b);
241  }
242  return tValues;
243 }
244 
245 
246 template<class Type>
248 {
249  #ifdef FULLDEBUG
250  checkSync();
251  #endif
252 
253  return sample01NoCheckSync<Type>();
254 }
255 
256 
257 template<class Type>
259 (
260  const label n
261 )
262 {
263  #ifdef FULLDEBUG
264  checkSync();
265  #endif
266 
267  tmp<Field<Type>> tValues(new Field<Type>(n));
268  for (label i = 0; i < n; ++ i)
269  {
270  tValues.ref()[i] = sample01NoCheckSync<Type>();
271  }
272  return tValues;
273 }
274 
275 
276 template<class Type>
277 inline Type Foam::randomGenerator::sampleAB(const Type& a, const Type& b)
278 {
279  #ifdef FULLDEBUG
280  checkSync();
281  #endif
282 
283  return sampleABNoCheckSync(a, b);
284 }
285 
286 
287 template<class Type>
289 (
290  const label n,
291  const Type& a,
292  const Type& b
293 )
294 {
295  #ifdef FULLDEBUG
296  checkSync();
297  #endif
298 
299  tmp<Field<Type>> tValues(new Field<Type>(n));
300  for (label i = 0; i < n; ++ i)
301  {
302  tValues.ref()[i] = sampleABNoCheckSync(a, b);
303  }
304  return tValues;
305 }
306 
307 
308 template<class Container>
309 inline void Foam::randomGenerator::permute(Container& l)
310 {
311  #ifdef FULLDEBUG
312  checkSync();
313  #endif
314 
315  for (label i = 0; i < l.size(); ++ i)
316  {
317  Swap(l[i], l[sampleABNoCheckSync<label>(i, l.size())]);
318  }
319 }
320 
321 
323 {
324  return randomGenerator(sample(), global_);
325 }
326 
327 
328 // ************************************************************************* //
label n
Pre-declare SubField and related Field type.
Definition: Field.H:83
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
Definition: Hash.H:53
static void scatter(const List< commsStruct > &comms, T &Value, const int tag, const label comm)
Scatter data. Distribute without modification. Reverse of gather.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
seed(const label s)
Construct from a label.
Random number generator.
scalar scalarAB(const scalar a, const scalar b)
Return a scalar uniformly distributed between two limits.
randomGenerator(const seed s, const bool global=false)
Construct from a seed.
scalar scalar01()
Return a scalar uniformly distributed between zero and one.
randomGenerator generator()
Create a randomly seeded generator.
Type sample01()
Return a type with components uniformly distributed between.
void permute(Container &l)
Randomly permute a list.
Type sampleAB(const Type &a, const Type &b)
Return a type with components uniformly distributed between two.
~randomGenerator()
Destructor.
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
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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))
volScalarField & b
Definition: createFields.H:25
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
void Swap(T &a, T &b)
Definition: Swap.H:43
uint8_t direction
Definition: direction.H:45
randomGenerator rndGen(653213)