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 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 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 :
43  x_((type(s) << 16) + 0x330E),
44  scalarNormalStored_(false),
45  scalarNormalValue_(0)
46 {}
47 
48 
49 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
50 
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
57 inline Foam::scalar Foam::Random::scalar01()
58 {
59  return scalar(sample())/(M >> 17);
60 }
61 
62 
63 inline Foam::scalar Foam::Random::scalarAB(const scalar a, const scalar b)
64 {
65  return a + scalar01()*(b - a);
66 }
67 
68 
69 template<class Type>
71 {
72  Type value;
73 
74  for (direction i = 0; i < pTraits<Type>::nComponents; ++ i)
75  {
76  value.component(i) = scalar01();
77  }
78 
79  return value;
80 }
81 
82 
83 template<>
84 inline Foam::scalar Foam::Random::sample01()
85 {
86  return scalar01();
87 }
88 
89 
90 template<>
92 {
93  return sample() % 2;
94 }
95 
96 
97 template<class Type>
98 inline Type Foam::Random::sampleAB(const Type& a, const Type& b)
99 {
100  return a + cmptMultiply(sample01<Type>(), b - a);
101 }
102 
103 
104 template<>
105 inline Foam::scalar Foam::Random::sampleAB(const scalar& a, const scalar& b)
106 {
107  return scalarAB(a, b);
108 }
109 
110 
111 template<>
113 {
114  return a + sample() % (b - a);
115 }
116 
117 
118 template<class Type>
120 {
121  Type value;
122 
123  for (direction i = 0; i < pTraits<Type>::nComponents; ++ i)
124  {
125  value.component(i) = scalarNormal();
126  }
127 
128  return value;
129 }
130 
131 
132 template<>
133 inline Foam::scalar Foam::Random::sampleNormal()
134 {
135  return scalarNormal();
136 }
137 
138 
139 // ************************************************************************* //
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 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
uint8_t direction
Definition: direction.H:45
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
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
~Random()
Destructor.
Definition: RandomI.H:51
scalar scalar01()
Advance the state and return a scalar sample from a uniform.
Definition: RandomI.H:57