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::PhiForZeroQ
85 (
86  const scalarField& x
87 ) const
88 {
89  tmp<scalarField> tPhiQ0(new scalarField(x.size(), 0));
90  scalarField& PhiQ0 = tPhiQ0.ref();
91 
92  forAll(distributions_, i)
93  {
94  PhiQ0 +=
95  (cumulativeStrengths_[i + 1] - cumulativeStrengths_[i])
96  *distributions_[i].Phi(0, x);
97  }
98 
99  return tPhiQ0;
100 }
101 
102 
103 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
104 
106 (
107  const unitConversion& units,
108  const dictionary& dict,
109  const label sampleQ,
111 )
112 :
114  (
115  typeName,
116  units,
117  dict,
118  sampleQ,
119  std::move(rndGen)
120  ),
121  cumulativeStrengths_(readCumulativeStrengths(dict))
122 {
123  const scalar min
124  (
125  dict.lookupBackwardsCompatible<scalar>({"min", "minValue"}, units)
126  );
127  const scalar max
128  (
129  dict.lookupBackwardsCompatible<scalar>({"max", "maxValue"}, units)
130  );
131  const scalarList mu
132  (
133  dict.lookupBackwardsCompatible<scalarList>({"mu", "expectation"}, units)
134  );
135  const scalarList sigma
136  (
137  dict.lookup<scalarList>("sigma", units)
138  );
139 
140  if
141  (
142  mu.size() != cumulativeStrengths_.size() - 1
143  || sigma.size() != cumulativeStrengths_.size() - 1
144  )
145  {
147  << type() << ": Differing numbers of means, standard deviations "
148  << "and strengths were given:" << nl << " mu = " << mu
149  << ", sigma = " << sigma << ", strength = "
150  << dict.lookup<scalarList>("strength") << abort(FatalIOError);
151  }
152 
153  distributions_.resize(mu.size());
154  forAll(mu, i)
155  {
156  distributions_.set
157  (
158  i,
159  new normal
160  (
161  0,
162  0,
163  rndGen_.generator(),
164  -1,
165  min,
166  max,
167  mu[i],
168  sigma[i]
169  )
170  );
171  }
172 
173  validateBounds(dict);
174  if (q() != 0) validatePositive(dict);
175  mean();
176 }
177 
178 
180 (
181  const multiNormal& d,
182  const label sampleQ
183 )
184 :
186  cumulativeStrengths_(d.cumulativeStrengths_),
187  distributions_(d.distributions_.size())
188 {
189  forAll(d.distributions_, i)
190  {
191  distributions_.set
192  (
193  i,
194  new normal
195  (
196  0,
197  0,
198  randomGenerator(d.distributions_[i].rndGen_),
199  -1,
200  d.distributions_[i].min_,
201  d.distributions_[i].max_,
202  d.distributions_[i].mu_,
203  d.distributions_[i].sigma_
204  )
205  );
206  }
207 }
208 
209 
210 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
211 
213 {}
214 
215 
216 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
217 
219 {
220  const scalar S = rndGen_.sample01<scalar>();
221 
222  label i = 0;
223  for
224  (
225  ;
226  i < cumulativeStrengths_.size() - 1
227  && cumulativeStrengths_[i + 1] < S;
228  ++ i
229  );
230 
231  const scalar s =
232  (S - cumulativeStrengths_[i])
233  /(cumulativeStrengths_[i + 1] - cumulativeStrengths_[i]);
234 
235  return distributions_[i].sampleForZeroQ(s);
236 }
237 
238 
240 {
241  return distributions_[0].min();
242 }
243 
244 
246 {
247  return distributions_[0].max();
248 }
249 
250 
252 (
253  Ostream& os,
254  const unitConversion& units
255 ) const
256 {
258 
259  writeEntry(os, "min", units, distributions_[0].min());
260  writeEntry(os, "max", units, distributions_[0].max());
261 
262  scalarList mu(distributions_.size()), sigma(distributions_.size());
263  forAll(distributions_, i)
264  {
265  mu[i] = distributions_[i].mu();
266  sigma[i] = distributions_[i].sigma();
267  }
268 
269  writeEntry(os, "mu", units, mu);
270  writeEntry(os, "sigma", units, sigma);
271 }
272 
273 
276 {
277  return distributions_[0].plotX(n);
278 }
279 
280 
281 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
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:239
virtual ~multiNormal()
Destructor.
Definition: multiNormal.C:212
virtual void write(Ostream &os, const unitConversion &units) const
Write stream.
Definition: multiNormal.C:252
virtual scalar sampleForZeroQ() const
Sample the distribution for zero effective size exponent.
Definition: multiNormal.C:218
multiNormal(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen)
Construct from a dictionary.
Definition: multiNormal.C:106
virtual scalar max() const
Return the maximum value.
Definition: multiNormal.C:245
virtual tmp< scalarField > plotX(const label n) const
Return coordinates to plot across the range of the distribution.
Definition: multiNormal.C:275
Normal distribution, scaled so that it spans between a specified minimum and maximum value,...
Definition: normal.H:74
Base class for distributions that have a closed integral form for the cumulative density function (CD...
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(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::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
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
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
static const char nl
Definition: Ostream.H:267
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)