FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

normal.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "normal.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <OpenFOAM/mathematicalConstants.H>
00029 
00030 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034     namespace pdfs
00035     {
00036         defineTypeNameAndDebug(normal, 0);
00037         addToRunTimeSelectionTable(pdf, normal, dictionary);
00038     }
00039 }
00040 
00041 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00042 
00043 Foam::pdfs::normal::normal(const dictionary& dict, Random& rndGen)
00044 :
00045     pdf(typeName, dict, rndGen),
00046     minValue_(readScalar(pdfDict_.lookup("minValue"))),
00047     maxValue_(readScalar(pdfDict_.lookup("maxValue"))),
00048     expectation_(readScalar(pdfDict_.lookup("expectation"))),
00049     variance_(readScalar(pdfDict_.lookup("variance"))),
00050     a_(0.147)
00051 {
00052     if (minValue_ < 0)
00053     {
00054         FatalErrorIn("normal::normal(const dictionary&, Random&)")
00055             << "Minimum value must be greater than zero. "
00056             << "Supplied minValue = " << minValue_
00057             << abort(FatalError);
00058     }
00059 
00060     if (maxValue_ < minValue_)
00061     {
00062         FatalErrorIn("normal::normal(const dictionary&, Random&)")
00063             << "Maximum value is smaller than the minimum value:"
00064             << "    maxValue = " << maxValue_ << ", minValue = " << minValue_
00065             << abort(FatalError);
00066     }
00067 }
00068 
00069 
00070 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00071 
00072 Foam::pdfs::normal::~normal()
00073 {}
00074 
00075 
00076 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00077 
00078 Foam::scalar Foam::pdfs::normal::sample() const
00079 {
00080 
00081     scalar a = erf((minValue_ - expectation_)/variance_);
00082     scalar b = erf((maxValue_ - expectation_)/variance_);
00083 
00084     scalar y = rndGen_.scalar01();
00085     scalar x = erfInv(y*(b - a) + a)*variance_ + expectation_;
00086 
00087     // Note: numerical approximation of the inverse function yields slight
00088     //       inaccuracies
00089 
00090     x = min(max(x, minValue_), maxValue_);
00091 
00092     return x;
00093 }
00094 
00095 
00096 Foam::scalar Foam::pdfs::normal::minValue() const
00097 {
00098     return minValue_;
00099 }
00100 
00101 
00102 Foam::scalar Foam::pdfs::normal::maxValue() const
00103 {
00104     return maxValue_;
00105 }
00106 
00107 
00108 Foam::scalar Foam::pdfs::normal::erfInv(const scalar y) const
00109 {
00110     scalar k = 2.0/(mathematicalConstant::pi*a_) +  0.5*log(1.0 - y*y);
00111     scalar h = log(1.0 - y*y)/a_;
00112     scalar x = sqrt(-k + sqrt(k*k - h));
00113     if (y < 0.0)
00114     {
00115         x *= -1.0;
00116     }
00117     return x;
00118 }
00119 
00120 
00121 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines