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  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "exponential.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace distributionModels
37 {
38  defineTypeNameAndDebug(exponential, 0);
40 }
41 }
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const dictionary& dict,
48  Random& rndGen
49 )
50 :
51  distributionModel(typeName, dict, rndGen),
52  lambda_(distributionModelDict_.get<scalar>("lambda"))
53 {
54  if (lambda_ < VSMALL)
55  {
57  << "Rate parameter cannot be equal to or less than zero:" << nl
58  << " lambda = " << lambda_
60  }
61 
62  check();
63 }
64 
65 
67 :
69  lambda_(p.lambda_)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  const scalar u = rndGen_.sample01<scalar>();
78  const scalar qMin = exp(-lambda_*minValue_);
79  const scalar qMax = exp(-lambda_*maxValue_);
80  return -(scalar(1)/lambda_)*log(qMin + u*(qMax - qMin));
81 }
82 
83 
85 {
86  return scalar(1)/lambda_;
87 }
88 
89 
90 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
addToRunTimeSelectionTable(distributionModel, binned, dictionary)
virtual scalar sample() const
Sample the distribution.
Definition: exponential.C:68
Macros for easy insertion into run-time selection tables.
Particle-size distribution model wherein random samples are drawn from the doubly-truncated exponenti...
Definition: exponential.H:185
defineTypeNameAndDebug(binned, 0)
dimensionedScalar exp(const dimensionedScalar &ds)
virtual void check() const
Check that the distribution model is valid.
virtual scalar meanValue() const
Return the theoretical mean of the distribution.
Definition: exponential.C:77
Random number generator.
Definition: Random.H:55
A library of runtime-selectable doubly-truncated probability distribution models. Returns random samp...
volScalarField & p
exponential(const dictionary &dict, Random &rndGen)
Construct from components.
Definition: exponential.C:39
Namespace for OpenFOAM.