31 inline uint64_t Foam::randomGenerator::seed::x(
const bool global)
const
33 const uint64_t localS =
36 return (localS << 16) + 0x330E;
40 inline void Foam::randomGenerator::checkSync()
const
44 uint64_t xMaster = x_;
49 <<
"Global random number generator is not synchronised"
56 inline uint64_t Foam::randomGenerator::sample()
63 inline Foam::scalar Foam::randomGenerator::scalar01NoCheckSync()
65 return scalar(sample())/(M >> 17);
69 inline Foam::scalar Foam::randomGenerator::scalarABNoCheckSync
75 return a + scalar01NoCheckSync()*(
b - a);
80 inline Type Foam::randomGenerator::sample01NoCheckSync()
83 for (
direction i = 0; i < pTraits<Type>::nComponents; ++ i)
85 value.component(i) = scalar01NoCheckSync();
92 inline Foam::scalar Foam::randomGenerator::sample01NoCheckSync()
94 return scalar01NoCheckSync();
99 inline Foam::label Foam::randomGenerator::sample01NoCheckSync()
106 inline Type Foam::randomGenerator::sampleABNoCheckSync
117 inline Foam::scalar Foam::randomGenerator::sampleABNoCheckSync
123 return scalarABNoCheckSync(a,
b);
134 return a + sample() % (
b - a);
190 return scalar01NoCheckSync();
204 for (
label i = 0; i <
n; ++ i)
206 tValues.
ref()[i] = scalar01NoCheckSync();
222 return scalarABNoCheckSync(a,
b);
238 for (
label i = 0; i <
n; ++ i)
240 tValues.
ref()[i] = scalarABNoCheckSync(a,
b);
253 return sample01NoCheckSync<Type>();
268 for (
label i = 0; i <
n; ++ i)
270 tValues.
ref()[i] = sample01NoCheckSync<Type>();
283 return sampleABNoCheckSync(a,
b);
300 for (
label i = 0; i <
n; ++ i)
302 tValues.
ref()[i] = sampleABNoCheckSync(a,
b);
308 template<
class Container>
315 for (
label i = 0; i < l.size(); ++ i)
317 Swap(l[i], l[sampleABNoCheckSync<label>(i, l.size())]);
Pre-declare SubField and related Field type.
Hash function class for primitives. All non-primitives used to hash entries on hash tables likely nee...
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.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
seed(const label s)
Construct from a label.
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.
T & ref() const
Return non-const reference or generate a fatal error.
A class for handling words, derived from string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
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))
errorManipArg< error, int > exit(error &err, const int errNo=1)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
randomGenerator rndGen(653213)