exponential.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 "exponential.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace distributionModels
36 {
37  defineTypeNameAndDebug(exponential, 0);
38  addToRunTimeSelectionTable(distributionModel, exponential, 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  lambda_(distributionModelDict_.get<scalar>("lambda"))
54 {
55  check();
56 }
57 
58 
60 :
62  minValue_(p.minValue_),
63  maxValue_(p.maxValue_),
64  lambda_(p.lambda_)
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
69 
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  scalar y = rndGen_.sample01<scalar>();
79  scalar K = exp(-lambda_*maxValue_) - exp(-lambda_*minValue_);
80  return -(1.0/lambda_)*log(exp(-lambda_*minValue_) + y*K);
81 }
82 
83 
85 {
86  return minValue_;
87 }
88 
89 
91 {
92  return maxValue_;
93 }
94 
95 
97 {
98  return 1.0/lambda_;
99 }
100 
101 
102 // ************************************************************************* //
Foam::Random
Random number generator.
Definition: Random.H:59
exponential.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::distributionModels::exponential::sample
virtual scalar sample() const
Sample the distributionModel.
Definition: exponential.C:76
Foam::distributionModels::exponential::~exponential
virtual ~exponential()
Destructor.
Definition: exponential.C:70
Foam::distributionModels::exponential::exponential
exponential(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: exponential.C:45
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
Foam::distributionModels::exponential::minValue
virtual scalar minValue() const
Return the minimum value.
Definition: exponential.C:84
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::distributionModels::exponential::maxValue
virtual scalar maxValue() const
Return the maximum value.
Definition: exponential.C:90
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
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distributionModels::exponential::meanValue
virtual scalar meanValue() const
Return the mean value.
Definition: exponential.C:96
Foam::distributionModels::defineTypeNameAndDebug
defineTypeNameAndDebug(binned, 0)
rndGen
Random rndGen
Definition: createFields.H:23
Foam::distributionModels::exponential
exponential distribution model
Definition: exponential.H:52
y
scalar y
Definition: LISASMDCalcMethod1.H:14