exponential.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-2018 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 "exponential.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace distributionModels
34 {
37 }
38 }
39 
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
41 
43 (
44  const dictionary& dict,
45  Random& rndGen
46 )
47 :
48  distributionModel(typeName, dict, rndGen),
49  minValue_(readScalar(distributionModelDict_.lookup("minValue"))),
50  maxValue_(readScalar(distributionModelDict_.lookup("maxValue"))),
51  lambda_(readScalar(distributionModelDict_.lookup("lambda")))
52 {
53  check();
54 }
55 
56 
58 :
60  minValue_(p.minValue_),
61  maxValue_(p.maxValue_),
62  lambda_(p.lambda_)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
67 
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
75 {
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);
79 }
80 
81 
83 {
84  return minValue_;
85 }
86 
87 
89 {
90  return maxValue_;
91 }
92 
93 
95 {
96  return 1.0/lambda_;
97 }
98 
99 
100 // ************************************************************************* //
virtual scalar minValue() const
Return the minimum value.
Definition: exponential.C:82
dimensionedScalar log(const dimensionedScalar &ds)
Random & rndGen_
Reference to the random number generator.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual scalar maxValue() const
Return the maximum value.
Definition: exponential.C:88
virtual scalar sample() const
Sample the distributionModel.
Definition: exponential.C:74
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
virtual ~exponential()
Destructor.
Definition: exponential.C:68
scalar y
exponential distribution model
Definition: exponential.H:50
Type sample01()
Advance the state and return a sample of a given type from a.
Definition: RandomI.H:70
dimensionedScalar exp(const dimensionedScalar &ds)
virtual scalar meanValue() const
Return the mean value.
Definition: exponential.C:94
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
Random number generator.
Definition: Random.H:57
defineTypeNameAndDebug(exponential, 0)
A library of runtime-selectable distribution models.
addToRunTimeSelectionTable(distributionModel, exponential, dictionary)
volScalarField & p
exponential(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: exponential.C:43
Namespace for OpenFOAM.