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-2023 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 dictionary& dict,
116  Random& rndGen,
117  const label sampleQ
118 )
119 :
121  (
122  typeName,
123  dict,
124  rndGen,
125  sampleQ
126  ),
127  cumulativeStrengths_(readCumulativeStrengths(dict))
128 {
129  const scalar min
130  (
131  dict.lookupBackwardsCompatible<scalar>({"min", "minValue"})
132  );
133  const scalar max
134  (
135  dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"})
136  );
137  const scalarList mu
138  (
139  dict.lookupBackwardsCompatible<scalarList>({"mu", "expectation"})
140  );
141  const scalarList sigma
142  (
143  dict.lookup<scalarList>("sigma")
144  );
145 
146  if
147  (
148  mu.size() != cumulativeStrengths_.size() - 1
149  || sigma.size() != cumulativeStrengths_.size() - 1
150  )
151  {
153  << type() << ": Differing numbers of means, standard deviations "
154  << "and strengths were given:" << nl << " mu = " << mu
155  << ", sigma = " << sigma << ", strength = "
156  << dict.lookup<scalarList>("strength") << abort(FatalIOError);
157  }
158 
159  distributions_.resize(mu.size());
160  forAll(mu, i)
161  {
162  distributions_.set
163  (
164  i,
165  new normal
166  (
167  rndGen_,
168  0,
169  0,
170  -1,
171  min,
172  max,
173  mu[i],
174  sigma[i]
175  )
176  );
177  }
178 
179  validateBounds(dict);
180  if (q() != 0) validatePositive(dict);
181  mean();
182  report();
183 }
184 
185 
187 (
188  const multiNormal& d,
189  const label sampleQ
190 )
191 :
193  cumulativeStrengths_(d.cumulativeStrengths_),
194  distributions_(d.distributions_.size())
195 {
196  forAll(d.distributions_, i)
197  {
198  distributions_.set
199  (
200  i,
201  new normal
202  (
203  rndGen_,
204  0,
205  0,
206  -1,
207  d.distributions_[i].min_,
208  d.distributions_[i].max_,
209  d.distributions_[i].mu_,
210  d.distributions_[i].sigma_
211  )
212  );
213  }
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
218 
220 {}
221 
222 
223 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
224 
226 {
227  if (q() == 0)
228  {
229  const scalar S = rndGen_.sample01<scalar>();
230 
231  const label n = cumulativeStrengths_.size() - 1;
232  label i = 0;
233  for (; i < n && cumulativeStrengths_[i + 1] < S; ++ i);
234 
235  const scalar s =
236  (S - cumulativeStrengths_[i])
237  /(cumulativeStrengths_[i + 1] - cumulativeStrengths_[i]);
238 
239  return distributions_[i].sampleForZeroQ(s);
240  }
241  else
242  {
244  }
245 }
246 
247 
249 {
250  return distributions_[0].min();
251 }
252 
253 
255 {
256  return distributions_[0].max();
257 }
258 
259 
262 {
263  return distributions_[0].x(n);
264 }
265 
266 
267 // ************************************************************************* //
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
Random number generator.
Definition: Random.H:58
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
Multiple superimposed normal distributions.
Definition: multiNormal.H:75
multiNormal(const dictionary &dict, Random &rndGen, const label sampleQ)
Construct from a dictionary.
Definition: multiNormal.C:114
virtual scalar min() const
Return the minimum value.
Definition: multiNormal.C:248
virtual scalar sample() const
Sample the distribution.
Definition: multiNormal.C:225
virtual tmp< scalarField > x(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: multiNormal.C:261
virtual ~multiNormal()
Destructor.
Definition: multiNormal.C:219
virtual scalar max() const
Return the maximum value.
Definition: multiNormal.C:254
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:253
scalar sample() const
Sample the distribution.
Definition: unintegrable.C:349
A class for managing temporary objects.
Definition: tmp.H:55
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
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 > &)
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.
IOerror FatalIOError
static const char nl
Definition: Ostream.H:260
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
Random rndGen(label(0))