multiNormal.C
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) 2011-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 "multiNormal.H"
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace distributions
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * Private Static Member Functions * * * * * * * * * * //
43 
44 Foam::scalarList Foam::distributions::multiNormal::readCumulativeStrengths
45 (
46  const dictionary& dict
47 )
48 {
49  const scalarList s(dict.lookup<scalarList>("strength"));
50 
51  const scalarList sHat(s/sum(s));
52 
53  scalarList cumSHat(s.size() + 1, scalar(0));
54  forAll(sHat, i)
55  {
56  cumSHat[i + 1] = cumSHat[i] + sHat[i];
57  }
58 
59  return cumSHat;
60 }
61 
62 
63 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
64 
65 Foam::tmp<Foam::scalarField> Foam::distributions::multiNormal::phi
66 (
67  const label q,
68  const scalarField& x
69 ) const
70 {
71  scalarField phiQ0(x.size(), 0);
72 
73  forAll(distributions_, i)
74  {
75  phiQ0 +=
76  (cumulativeStrengths_[i + 1] - cumulativeStrengths_[i])
77  *distributions_[i].phi(0, x);
78  }
79 
80  return integerPow(x, q)*phiQ0;
81 }
82 
83 
84 Foam::tmp<Foam::scalarField> Foam::distributions::multiNormal::Phi
85 (
86  const label q,
87  const scalarField& x
88 ) const
89 {
90  if (q == 0)
91  {
92  tmp<scalarField> tPhiQ0(new scalarField(x.size(), 0));
93  scalarField& PhiQ0 = tPhiQ0.ref();
94 
95  forAll(distributions_, i)
96  {
97  PhiQ0 +=
98  (cumulativeStrengths_[i + 1] - cumulativeStrengths_[i])
99  *distributions_[i].Phi(0, x);
100  }
101 
102  return tPhiQ0;
103  }
104  else
105  {
106  return unintegrableForNonZeroQ::Phi(q, x);
107  }
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
114 (
115  const unitConversion& units,
116  const dictionary& dict,
117  const label sampleQ,
119 )
120 :
122  (
123  typeName,
124  units,
125  dict,
126  sampleQ,
127  std::move(rndGen)
128  ),
129  cumulativeStrengths_(readCumulativeStrengths(dict))
130 {
131  const scalar min
132  (
133  dict.lookupBackwardsCompatible<scalar>({"min", "minValue"}, units)
134  );
135  const scalar max
136  (
137  dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"}, units)
138  );
139  const scalarList mu
140  (
141  dict.lookupBackwardsCompatible<scalarList>({"mu", "expectation"}, units)
142  );
143  const scalarList sigma
144  (
145  dict.lookup<scalarList>("sigma", units)
146  );
147 
148  if
149  (
150  mu.size() != cumulativeStrengths_.size() - 1
151  || sigma.size() != cumulativeStrengths_.size() - 1
152  )
153  {
155  << type() << ": Differing numbers of means, standard deviations "
156  << "and strengths were given:" << nl << " mu = " << mu
157  << ", sigma = " << sigma << ", strength = "
158  << dict.lookup<scalarList>("strength") << abort(FatalIOError);
159  }
160 
161  distributions_.resize(mu.size());
162  forAll(mu, i)
163  {
164  distributions_.set
165  (
166  i,
167  new normal
168  (
169  0,
170  0,
171  rndGen_.generator(),
172  -1,
173  min,
174  max,
175  mu[i],
176  sigma[i]
177  )
178  );
179  }
180 
181  validateBounds(dict);
182  if (q() != 0) validatePositive(dict);
183  mean();
184 }
185 
186 
188 (
189  const multiNormal& d,
190  const label sampleQ
191 )
192 :
194  cumulativeStrengths_(d.cumulativeStrengths_),
195  distributions_(d.distributions_.size())
196 {
197  forAll(d.distributions_, i)
198  {
199  distributions_.set
200  (
201  i,
202  new normal
203  (
204  0,
205  0,
206  randomGenerator(d.distributions_[i].rndGen_),
207  -1,
208  d.distributions_[i].min_,
209  d.distributions_[i].max_,
210  d.distributions_[i].mu_,
211  d.distributions_[i].sigma_
212  )
213  );
214  }
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
219 
221 {}
222 
223 
224 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225 
227 {
228  if (q() == 0)
229  {
230  const scalar S = rndGen_.sample01<scalar>();
231 
232  const label n = cumulativeStrengths_.size() - 1;
233  label i = 0;
234  for (; i < n && cumulativeStrengths_[i + 1] < S; ++ i);
235 
236  const scalar s =
237  (S - cumulativeStrengths_[i])
238  /(cumulativeStrengths_[i + 1] - cumulativeStrengths_[i]);
239 
240  return distributions_[i].sampleForZeroQ(s);
241  }
242  else
243  {
245  }
246 }
247 
248 
250 {
251  return distributions_[0].min();
252 }
253 
254 
256 {
257  return distributions_[0].max();
258 }
259 
260 
262 (
263  Ostream& os,
264  const unitConversion& units
265 ) const
266 {
268 
269  writeEntry(os, "min", units, distributions_[0].min());
270  writeEntry(os, "max", units, distributions_[0].max());
271 
272  scalarList mu(distributions_.size()), sigma(distributions_.size());
273  forAll(distributions_, i)
274  {
275  mu[i] = distributions_[i].mu();
276  sigma[i] = distributions_[i].sigma();
277  }
278 
279  writeEntry(os, "mu", units, mu);
280  writeEntry(os, "sigma", units, sigma);
281 }
282 
283 
286 {
287  return distributions_[0].x(n);
288 }
289 
290 
291 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Base class for statistical distributions.
Definition: distribution.H:76
Multiple superimposed normal distributions.
Definition: multiNormal.H:75
virtual scalar min() const
Return the minimum value.
Definition: multiNormal.C:249
virtual scalar sample() const
Sample the distribution.
Definition: multiNormal.C:226
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: multiNormal.C:285
virtual ~multiNormal()
Destructor.
Definition: multiNormal.C:220
virtual void write(Ostream &os, const unitConversion &units) const
Write stream.
Definition: multiNormal.C:262
multiNormal(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
Definition: multiNormal.C:114
virtual scalar max() const
Return the maximum value.
Definition: multiNormal.C:255
Normal distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: normal.H:74
virtual tmp< scalarField > Phi(const label q, const scalarField &x) const
Return values of the un-normalised CDF for the given size exponent.
Definition: unintegrable.C:254
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:351
Random number generator.
A class for managing temporary objects.
Definition: tmp.H:55
Unit conversion structure. Contains the associated dimensions and the multiplier with which to conver...
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
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))
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar mu
Atomic mass unit.
defineTypeNameAndDebug(exponential, 0)
addToRunTimeSelectionTable(distribution, exponential, dictionary)
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
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 > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
scalar integerPow(const scalar x, const label e)
Compute the power of the number x to the integer e.
Definition: scalarI.H:30
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:50
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const HashTable< unitConversion > & units()
Get the table of unit conversions.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
randomGenerator rndGen(653213)