RandomI.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-2022 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 "Random.H"
27 
28 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
29 
30 inline Foam::Random::type Foam::Random::sample()
31 {
32  x_ = (A*x_ + C) % M;
33 
34  return x_ >> 17;
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 :
42  x_((type(s) << 16) + 0x330E),
43  scalarNormalStored_(false),
44  scalarNormalValue_(0)
45 {}
46 
47 
48 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
49 
51 {}
52 
53 
54 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
55 
56 inline Foam::scalar Foam::Random::scalar01()
57 {
58  return scalar(sample())/(M >> 17);
59 }
60 
61 
62 inline Foam::scalar Foam::Random::scalarAB(const scalar a, const scalar b)
63 {
64  return a + scalar01()*(b - a);
65 }
66 
67 
68 template<class Type>
70 {
71  Type value;
72 
73  for (direction i = 0; i < pTraits<Type>::nComponents; ++ i)
74  {
75  value.component(i) = scalar01();
76  }
77 
78  return value;
79 }
80 
81 
82 template<>
83 inline Foam::scalar Foam::Random::sample01()
84 {
85  return scalar01();
86 }
87 
88 
89 template<>
91 {
92  return sample() % 2;
93 }
94 
95 
96 template<class Type>
97 inline Type Foam::Random::sampleAB(const Type& a, const Type& b)
98 {
99  return a + cmptMultiply(sample01<Type>(), b - a);
100 }
101 
102 
103 template<>
104 inline Foam::scalar Foam::Random::sampleAB(const scalar& a, const scalar& b)
105 {
106  return scalarAB(a, b);
107 }
108 
109 
110 template<>
112 {
113  return a + sample() % (b - a);
114 }
115 
116 
117 template<class Type>
119 {
120  Type value;
121 
122  for (direction i = 0; i < pTraits<Type>::nComponents; ++ i)
123  {
124  value.component(i) = scalarNormal();
125  }
126 
127  return value;
128 }
129 
130 
131 template<>
132 inline Foam::scalar Foam::Random::sampleNormal()
133 {
134  return scalarNormal();
135 }
136 
137 
138 template<class Container>
139 inline void Foam::Random::permute(Container& l)
140 {
141  for (label i = 0; i < l.size(); ++ i)
142  {
143  Swap(l[i], l[sampleAB<label>(i, l.size())]);
144  }
145 }
146 
147 
148 // ************************************************************************* //
#define M(I)
scalar scalarAB(const scalar a, const scalar b)
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:62
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:56
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:69
~Random()
Destructor.
Definition: RandomI.H:50
void permute(Container &l)
Randomly permute a list.
Definition: RandomI.H:139
Type sampleAB(const Type &a, const Type &b)
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:97
Type sampleNormal()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:118
Random(const label s)
Construct from a seed.
Definition: RandomI.H:40
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:27
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
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
void Swap(T &a, T &b)
Definition: Swap.H:43
uint8_t direction
Definition: direction.H:45