BUG: pointHitSort: define operator<
[OpenFOAM-1.7.x.git] / src / thermophysicalModels / pdfs / normal / normal.C
blobb62c65554590211312355297a8969f4502b9320e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "normal.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "mathematicalConstants.H"
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 namespace Foam
34     namespace pdfs
35     {
36         defineTypeNameAndDebug(normal, 0);
37         addToRunTimeSelectionTable(pdf, normal, dictionary);
38     }
41 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
43 Foam::pdfs::normal::normal(const dictionary& dict, Random& rndGen)
45     pdf(typeName, dict, rndGen),
46     minValue_(readScalar(pdfDict_.lookup("minValue"))),
47     maxValue_(readScalar(pdfDict_.lookup("maxValue"))),
48     expectation_(readScalar(pdfDict_.lookup("expectation"))),
49     variance_(readScalar(pdfDict_.lookup("variance"))),
50     a_(0.147)
52     if (minValue_ < 0)
53     {
54         FatalErrorIn("normal::normal(const dictionary&, Random&)")
55             << "Minimum value must be greater than zero. "
56             << "Supplied minValue = " << minValue_
57             << abort(FatalError);
58     }
60     if (maxValue_ < minValue_)
61     {
62         FatalErrorIn("normal::normal(const dictionary&, Random&)")
63             << "Maximum value is smaller than the minimum value:"
64             << "    maxValue = " << maxValue_ << ", minValue = " << minValue_
65             << abort(FatalError);
66     }
70 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
72 Foam::pdfs::normal::~normal()
76 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
78 Foam::scalar Foam::pdfs::normal::sample() const
81     scalar a = erf((minValue_ - expectation_)/variance_);
82     scalar b = erf((maxValue_ - expectation_)/variance_);
84     scalar y = rndGen_.scalar01();
85     scalar x = erfInv(y*(b - a) + a)*variance_ + expectation_;
87     // Note: numerical approximation of the inverse function yields slight
88     //       inaccuracies
90     x = min(max(x, minValue_), maxValue_);
92     return x;
96 Foam::scalar Foam::pdfs::normal::minValue() const
98     return minValue_;
102 Foam::scalar Foam::pdfs::normal::maxValue() const
104     return maxValue_;
108 Foam::scalar Foam::pdfs::normal::erfInv(const scalar y) const
110     scalar k = 2.0/(mathematicalConstant::pi*a_) +  0.5*log(1.0 - y*y);
111     scalar h = log(1.0 - y*y)/a_;
112     scalar x = sqrt(-k + sqrt(k*k - h));
113     if (y < 0.0)
114     {
115         x *= -1.0;
116     }
117     return x;
121 // ************************************************************************* //