ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / src / lagrangian / distributionModels / exponential / exponential.C
blob72df563826687e789fc5efbb7a0c822e6eabd169
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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 "exponential.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 namespace Foam
33     namespace distributionModels
34     {
35         defineTypeNameAndDebug(exponential, 0);
36         addToRunTimeSelectionTable(distributionModel, exponential, dictionary);
37     }
40 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
42 Foam::distributionModels::exponential::exponential
44     const dictionary& dict,
45     cachedRandom& rndGen
48     distributionModel(typeName, dict, rndGen),
49     minValue_(readScalar(distributionModelDict_.lookup("minValue"))),
50     maxValue_(readScalar(distributionModelDict_.lookup("maxValue"))),
51     lambda_(readScalar(distributionModelDict_.lookup("lambda")))
53     check();
57 Foam::distributionModels::exponential::exponential(const exponential& p)
59     distributionModel(p),
60     minValue_(p.minValue_),
61     maxValue_(p.maxValue_),
62     lambda_(p.lambda_)
66 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
68 Foam::distributionModels::exponential::~exponential()
72 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
74 Foam::scalar Foam::distributionModels::exponential::sample() const
76     scalar y = rndGen_.sample01<scalar>();
77     scalar K = exp(-lambda_*maxValue_) - exp(-lambda_*minValue_);
78     return -(1.0/lambda_)*log(exp(-lambda_*minValue_) + y*K);
82 Foam::scalar Foam::distributionModels::exponential::minValue() const
84     return minValue_;
88 Foam::scalar Foam::distributionModels::exponential::maxValue() const
90     return maxValue_;
94 // ************************************************************************* //