multiNormal.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "multiNormal.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace distributionModels
36 {
37  defineTypeNameAndDebug(multiNormal, 0);
38  addToRunTimeSelectionTable(distributionModel, multiNormal, dictionary);
39 }
40 }
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  Random& rndGen
48 )
49 :
50  distributionModel(typeName, dict, rndGen),
51  minValue_(distributionModelDict_.get<scalar>("minValue")),
52  maxValue_(distributionModelDict_.get<scalar>("maxValue")),
53  range_(maxValue_ - minValue_),
54  expectation_(distributionModelDict_.lookup("expectation")),
55  variance_(distributionModelDict_.lookup("variance")),
56  strength_(distributionModelDict_.lookup("strength"))
57 {
58  check();
59 
60  scalar sMax = 0;
61  label n = strength_.size();
62  for (label i=0; i<n; i++)
63  {
64  scalar x = expectation_[i];
65  scalar s = strength_[i];
66  for (label j=0; j<n; j++)
67  {
68  if (i!=j)
69  {
70  scalar x2 = (x - expectation_[j])/variance_[j];
71  scalar y = exp(-0.5*x2*x2);
72  s += strength_[j]*y;
73  }
74  }
75 
76  sMax = max(sMax, s);
77  }
78 
79  for (label i=0; i<n; i++)
80  {
81  strength_[i] /= sMax;
82  }
83 }
84 
85 
87 :
89  minValue_(p.minValue_),
90  maxValue_(p.maxValue_),
91  range_(p.range_),
92  expectation_(p.expectation_),
93  variance_(p.variance_),
94  strength_(p.strength_)
95 {}
96 
97 
98 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
99 
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
107 {
108  scalar y = 0;
109  scalar x = 0;
110  label n = expectation_.size();
111  bool success = false;
112 
113  while (!success)
114  {
115  x = minValue_ + range_*rndGen_.sample01<scalar>();
116  y = rndGen_.sample01<scalar>();
117  scalar p = 0.0;
118 
119  for (label i=0; i<n; i++)
120  {
121  scalar nu = expectation_[i];
122  scalar sigma = variance_[i];
123  scalar s = strength_[i];
124  scalar v = (x - nu)/sigma;
125  p += s*exp(-0.5*v*v);
126  }
127 
128  if (y<p)
129  {
130  success = true;
131  }
132  }
133 
134  return x;
135 }
136 
137 
139 {
140  return minValue_;
141 }
142 
143 
145 {
146  return maxValue_;
147 }
148 
149 
151 {
152  scalar mean = 0.0;
153  forAll(strength_, i)
154  {
155  mean += strength_[i]*expectation_[i];
156  }
157 
158  return mean;
159 }
160 
161 
162 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
multiNormal.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::distributionModels::multiNormal::~multiNormal
virtual ~multiNormal()
Destructor.
Definition: multiNormal.C:100
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
success
bool success
Definition: LISASMDCalcMethod1.H:16
Foam::distributionModels::multiNormal::maxValue
virtual scalar maxValue() const
Return the maximum value.
Definition: multiNormal.C:144
Foam::distributionModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::distributionModel
A library of runtime-selectable distribution models.
Definition: distributionModel.H:70
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::distributionModels::multiNormal
A multiNormal distribution model.
Definition: multiNormal.H:57
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distributionModels::multiNormal::minValue
virtual scalar minValue() const
Return the minimum value.
Definition: multiNormal.C:138
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::distributionModels::multiNormal::meanValue
virtual scalar meanValue() const
Return the mean value.
Definition: multiNormal.C:150
rndGen
Random rndGen
Definition: createFields.H:23
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::distributionModels::multiNormal::sample
virtual scalar sample() const
Sample the distributionModel.
Definition: multiNormal.C:106
Foam::distributionModels::multiNormal::multiNormal
multiNormal(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: multiNormal.C:45
y
scalar y
Definition: LISASMDCalcMethod1.H:14