#include <G4PolynomialPDF.hh>
|
| G4PolynomialPDF (size_t n=0, const double *coeffs=nullptr, G4double x1=0, G4double x2=1) |
|
| ~G4PolynomialPDF () |
|
void | SetNCoefficients (size_t n) |
|
size_t | GetNCoefficients () const |
|
void | SetCoefficients (const std::vector< G4double > &v) |
|
G4double | GetCoefficient (size_t i) const |
|
void | SetCoefficient (size_t i, G4double value) |
|
void | SetCoefficients (size_t n, const G4double *coeffs) |
|
void | SetDomain (G4double x1, G4double x2) |
|
void | Normalize () |
|
G4double | Evaluate (G4double x, G4int ddxPower=0) |
|
G4double | GetRandomX () |
|
void | SetTolerance (G4double tolerance) |
|
G4double | GetX (G4double p, G4double x1, G4double x2, G4int ddxPower=0, G4double guess=1.e99, G4bool bisect=true) |
|
G4double | EvalInverseCDF (G4double p) |
|
G4double | Bisect (G4double p, G4double x1, G4double x2) |
|
void | Dump () |
|
Definition at line 49 of file G4PolynomialPDF.hh.
◆ G4PolynomialPDF()
G4PolynomialPDF::G4PolynomialPDF |
( |
size_t |
n = 0 , |
|
|
const double * |
coeffs = nullptr , |
|
|
G4double |
x1 = 0 , |
|
|
G4double |
x2 = 1 |
|
) |
| |
Definition at line 44 of file G4PolynomialPDF.cc.
void SetCoefficients(const std::vector< G4double > &v)
void SetNCoefficients(size_t n)
◆ ~G4PolynomialPDF()
G4PolynomialPDF::~G4PolynomialPDF |
( |
| ) |
|
|
inline |
◆ Bisect()
Definition at line 334 of file G4PolynomialPDF.cc.
G4double GetX(G4double p, G4double x1, G4double x2, G4int ddxPower=0, G4double guess=1.e99, G4bool bisect=true)
G4double Bisect(G4double p, G4double x1, G4double x2)
G4double Evaluate(G4double x, G4int ddxPower=0)
◆ Dump()
void G4PolynomialPDF::Dump |
( |
| ) |
|
Definition at line 343 of file G4PolynomialPDF.cc.
345 G4cout <<
"G4PolynomialPDF::Dump() - PDF(x) = ";
350 if(i>1)
G4cout <<
"^" << i;
353 G4cout <<
"G4PolynomialPDF::Dump() - Interval: " <<
fX1 <<
" <= x < "
size_t GetNCoefficients() const
G4GLOB_DLL std::ostream G4cout
G4double GetCoefficient(size_t i) const
◆ EvalInverseCDF()
Definition at line 98 of file G4PolynomialPDF.hh.
G4double GetX(G4double p, G4double x1, G4double x2, G4int ddxPower=0, G4double guess=1.e99, G4bool bisect=true)
◆ Evaluate()
Evaluate f(x) ddxPower = -1: f = CDF ddxPower = 0: f = PDF ddxPower = 1: f = (d/dx) PDF ddxPower = 2: f = (d2/dx2) PDF
Definition at line 109 of file G4PolynomialPDF.cc.
116 if(ddxPower < -1 || ddxPower > 2) {
117 G4cout <<
"G4PolynomialPDF::GetX(): ddxPower " << ddxPower
118 <<
" not implemented" <<
G4endl;
131 else if(ddxPower == 1) {
134 else if(ddxPower == 2) {
size_t GetNCoefficients() const
G4GLOB_DLL std::ostream G4cout
G4double GetCoefficient(size_t i) const
◆ GetCoefficient()
G4double G4PolynomialPDF::GetCoefficient |
( |
size_t |
i | ) |
const |
|
inline |
◆ GetNCoefficients()
size_t G4PolynomialPDF::GetNCoefficients |
( |
| ) |
const |
|
inline |
◆ GetRandomX()
G4double G4PolynomialPDF::GetRandomX |
( |
| ) |
|
Definition at line 181 of file G4PolynomialPDF.cc.
186 G4cout <<
"G4PolynomialPDF::GetRandomX(): PDF has negative values, returning 0..." <<
G4endl;
G4GLOB_DLL std::ostream G4cout
G4bool HasNegativeMinimum(G4double x1, G4double x2)
G4double EvalInverseCDF(G4double p)
◆ GetX()
Find a value of X between x1 and x2 at which f(x) = p. ddxPower = -1: f = CDF ddxPower = 0: f = PDF ddxPower = 1: f = (d/dx) PDF Uses the Newton-Raphson method to find the zero of f(x) - p. If not found in range, returns the nearest boundary
Definition at line 194 of file G4PolynomialPDF.cc.
206 G4cout <<
"G4PolynomialPDF::GetX(): no PDF defined!" <<
G4endl;
209 if(ddxPower < -1 || ddxPower > 1) {
210 G4cout <<
"G4PolynomialPDF::GetX(): ddxPower " << ddxPower
211 <<
" not implemented" <<
G4endl;
214 if(ddxPower == -1 && (p<0 || p>1)) {
215 G4cout <<
"G4PolynomialPDF::GetX(): p is out of range" <<
G4endl;
220 if(
x2 <=
x1 || x1 < fX1 || x2 >
fX2) {
221 G4cout <<
"G4PolynomialPDF::GetX(): domain must have fX1 <= x1 < x2 <= fX2. " 222 <<
"You sent x1 = " <<
x1 <<
", x2 = " <<
x2 <<
"." <<
G4endl;
236 if(ddxPower == 1) slope *= 2.;
241 else if(value >
x2) {
256 if(ddxPower == 1) b *= 2.;
258 if(ddxPower == 1) a *= 3;
259 else if(ddxPower == -1) a *= 0.5;
260 double sqrtFactor = b*b - 4.*a*
c;
261 if(sqrtFactor < 0)
return x2;
262 sqrtFactor = sqrt(sqrtFactor)/2./fabs(a);
263 G4double valueMinus = -b/2./a - sqrtFactor;
264 if(valueMinus >=
x1 && valueMinus <=
x2)
return valueMinus;
265 else if(valueMinus >
x2)
return x2;
266 G4double valuePlus = -b/2./a + sqrtFactor;
267 if(valuePlus >=
x1 && valuePlus <=
x2)
return valuePlus;
268 else if(valuePlus <
x1)
return x2;
269 return (
x1-valueMinus <= valuePlus-
x2) ?
x1 :
x2;
274 if(guess < x1 || guess >
x2) guess = (x2+
x1)*0.5;
276 size_t iterations = 0;
289 else if(ddxPower == 0) {
299 if(f == 0)
return guess;
301 G4cout <<
"G4PolynomialPDF::GetX(): got f != 0 but slope = 0 for ddxPower = " 305 lastChange = - f/dfdx;
307 if(guess + lastChange <
x1) {
308 lastChange =
x1 - guess;
309 }
else if(guess + lastChange > x2) {
310 lastChange = x2 - guess;
314 lastChange /= (fX2-
fX1);
317 if(iterations > 50) {
319 G4cout <<
"G4PolynomialPDF::GetX(): got stuck searching for " << p
320 <<
" between " <<
x1 <<
" and " << x2 <<
" with ddxPower = " 322 <<
". Last guess was " << guess <<
"." <<
G4endl;
324 if(ddxPower==-1 && bisect) {
G4double Bisect(G4double p, G4double x1, G4double x2)
size_t GetNCoefficients() const
G4GLOB_DLL std::ostream G4cout
G4double GetCoefficient(size_t i) const
◆ HasNegativeMinimum()
Definition at line 142 of file G4PolynomialPDF.cc.
147 if(x1 < fX1 || x2 >
fX2 ||
x2 <
x1) {
148 G4cout <<
"G4PolynomialPDF::HasNegativeMinimum(): Invalid range " 166 if(xMin <
x1) xMin =
x1;
167 if(xMin >
x2) xMin =
x2;
G4double GetX(G4double p, G4double x1, G4double x2, G4int ddxPower=0, G4double guess=1.e99, G4bool bisect=true)
size_t GetNCoefficients() const
G4GLOB_DLL std::ostream G4cout
G4bool HasNegativeMinimum(G4double x1, G4double x2)
G4double Evaluate(G4double x, G4int ddxPower=0)
G4double GetCoefficient(size_t i) const
◆ Normalize()
void G4PolynomialPDF::Normalize |
( |
| ) |
|
Normalize PDF to 1 over domain fX1 to fX2. Double-check that the highest-order coefficient is non-zero.
Definition at line 81 of file G4PolynomialPDF.cc.
98 G4cout <<
"G4PolynomialPDF::Normalize() PDF has non-positive area: "
std::vector< G4double > fCoefficients
size_t GetNCoefficients() const
G4GLOB_DLL std::ostream G4cout
void SetCoefficient(size_t i, G4double value)
G4double GetCoefficient(size_t i) const
◆ SetCoefficient()
void G4PolynomialPDF::SetCoefficient |
( |
size_t |
i, |
|
|
G4double |
value |
|
) |
| |
◆ SetCoefficients() [1/2]
void G4PolynomialPDF::SetCoefficients |
( |
const std::vector< G4double > & |
v | ) |
|
|
inline |
◆ SetCoefficients() [2/2]
void G4PolynomialPDF::SetCoefficients |
( |
size_t |
n, |
|
|
const G4double * |
coeffs |
|
) |
| |
Definition at line 59 of file G4PolynomialPDF.cc.
void SetNCoefficients(size_t n)
size_t GetNCoefficients() const
void SetCoefficient(size_t i, G4double value)
◆ SetDomain()
Definition at line 69 of file G4PolynomialPDF.cc.
72 G4cout <<
"G4PolynomialPDF::SetDomain(): Invalide domain! " 73 <<
"(x1 = " <<
x1 <<
", x2 = " <<
x2 <<
")." <<
G4endl;
G4GLOB_DLL std::ostream G4cout
◆ SetNCoefficients()
void G4PolynomialPDF::SetNCoefficients |
( |
size_t |
n | ) |
|
|
inline |
◆ SetTolerance()
void G4PolynomialPDF::SetTolerance |
( |
G4double |
tolerance | ) |
|
|
inline |
◆ fChanged
G4bool G4PolynomialPDF::fChanged |
|
protected |
◆ fCoefficients
std::vector<G4double> G4PolynomialPDF::fCoefficients |
|
protected |
◆ fTolerance
◆ fX1
◆ fX2
The documentation for this class was generated from the following files: