Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Numerics.hh
Go to the documentation of this file.
1// Collection of numerical routines
2
3#ifndef G_NUMERICS_H
4#define G_NUMERICS_H
5
6#include <vector>
7#include <complex>
8
9namespace Garfield {
10
11namespace Numerics {
12
13// Linear algebra routines from CERNLIB
14void Dfact(const int n, std::vector<std::vector<double> >& a,
15 std::vector<int>& ir, int& ifail, double& det, int& jfail);
16void Dfeqn(const int n, std::vector<std::vector<double> >& a,
17 std::vector<int>& ir, std::vector<double>& b);
18void Dfinv(const int n, std::vector<std::vector<double> >& a,
19 std::vector<int>& ir);
20void Deqinv(const int n, std::vector<std::vector<double> >& a, int& ifail,
21 std::vector<double>& b);
22
23void Cfact(const int n, std::vector<std::vector<std::complex<double> > >& a,
24 std::vector<int>& ir, int& ifail, std::complex<double>& det,
25 int& jfail);
26void Cfinv(const int n, std::vector<std::vector<std::complex<double> > >& a,
27 std::vector<int>& ir);
28void Cinv(const int n, std::vector<std::vector<std::complex<double> > >& a,
29 int& ifail);
30
31// Numerical integration using 15-point Gauss-Kronrod algorithm
32double GaussKronrod15(double (*f)(const double), const double a,
33 const double b);
34
35// Modified Bessel functions.
36// Series expansions from Abramowitz and Stegun.
37inline double BesselI0S(const double xx) {
38 return 1. + 3.5156229 * pow(xx / 3.75, 2) + 3.0899424 * pow(xx / 3.75, 4) +
39 1.2067492 * pow(xx / 3.75, 6) + 0.2659732 * pow(xx / 3.75, 8) +
40 0.0360768 * pow(xx / 3.75, 10) + 0.0045813 * pow(xx / 3.75, 12);
41}
42
43inline double BesselI1S(const double xx) {
44 return xx *
45 (0.5 + 0.87890594 * pow(xx / 3.75, 2) +
46 0.51498869 * pow(xx / 3.75, 4) + 0.15084934 * pow(xx / 3.75, 6) +
47 0.02658733 * pow(xx / 3.75, 8) + 0.00301532 * pow(xx / 3.75, 10) +
48 0.00032411 * pow(xx / 3.75, 12));
49}
50
51inline double BesselK0S(const double xx) {
52 return -log(xx / 2.) * BesselI0S(xx) - 0.57721566 +
53 0.42278420 * pow(xx / 2., 2) + 0.23069756 * pow(xx / 2., 4) +
54 0.03488590 * pow(xx / 2., 6) + 0.00262698 * pow(xx / 2., 8) +
55 0.00010750 * pow(xx / 2., 10) + 0.00000740 * pow(xx / 2., 12);
56}
57
58inline double BesselK0L(const double xx) {
59 return (exp(-xx) / sqrt(xx)) *
60 (1.25331414 - 0.07832358 * (2. / xx) + 0.02189568 * pow(2. / xx, 2) -
61 0.01062446 * pow(2. / xx, 3) + 0.00587872 * pow(2. / xx, 4) -
62 0.00251540 * pow(2. / xx, 5) + 0.00053208 * pow(2. / xx, 6));
63}
64
65inline double BesselK1S(const double xx) {
66 return log(xx / 2.) * BesselI1S(xx) +
67 (1. / xx) *
68 (1. + 0.15443144 * pow(xx / 2., 2) - 0.67278579 * pow(xx / 2., 4) -
69 0.18156897 * pow(xx / 2., 6) - 0.01919402 * pow(xx / 2., 8) -
70 0.00110404 * pow(xx / 2., 10) - 0.00004686 * pow(xx / 2., 12));
71}
72
73inline double BesselK1L(const double xx) {
74 return (exp(-xx) / sqrt(xx)) *
75 (1.25331414 + 0.23498619 * (2. / xx) - 0.03655620 * pow(2. / xx, 2) +
76 0.01504268 * pow(2. / xx, 3) - 0.00780353 * pow(2. / xx, 4) +
77 0.00325614 * pow(2. / xx, 5) - 0.00068245 * pow(2. / xx, 6));
78}
79
80double Divdif(const std::vector<double>& f, const std::vector<double>& a,
81 int nn, double x, int mm);
82
83bool Boxin3(std::vector<std::vector<std::vector<double> > >& value,
84 std::vector<double>& xAxis, std::vector<double>& yAxis,
85 std::vector<double>& zAxis, int nx, int ny, int nz, double xx,
86 double yy, double zz, double& f, int iOrder);
87}
88}
89
90#endif
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:336
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:376
double BesselK0L(const double xx)
Definition: Numerics.hh:58
double BesselI0S(const double xx)
Definition: Numerics.hh:37
double BesselI1S(const double xx)
Definition: Numerics.hh:43
void Dfinv(const int n, std::vector< std::vector< double > > &a, std::vector< int > &ir)
Definition: Numerics.cc:138
double BesselK1S(const double xx)
Definition: Numerics.hh:65
void Dfact(const int n, std::vector< std::vector< double > > &a, std::vector< int > &ir, int &ifail, double &det, int &jfail)
Definition: Numerics.cc:10
void Cfact(const int n, std::vector< std::vector< std::complex< double > > > &a, std::vector< int > &ir, int &ifail, std::complex< double > &det, int &jfail)
Definition: Numerics.cc:333
bool Boxin3(std::vector< std::vector< std::vector< double > > > &value, std::vector< double > &xAxis, std::vector< double > &yAxis, std::vector< double > &zAxis, int nx, int ny, int nz, double xx, double yy, double zz, double &f, int iOrder)
Definition: Numerics.cc:771
double GaussKronrod15(double(*f)(const double), const double a, const double b)
Definition: Numerics.cc:599
void Cinv(const int n, std::vector< std::vector< std::complex< double > > > &a, int &ifail)
Definition: Numerics.cc:494
double Divdif(const std::vector< double > &f, const std::vector< double > &a, int nn, double x, int mm)
Definition: Numerics.cc:650
double BesselK0S(const double xx)
Definition: Numerics.hh:51
void Deqinv(const int n, std::vector< std::vector< double > > &a, int &ifail, std::vector< double > &b)
Definition: Numerics.cc:216
void Cfinv(const int n, std::vector< std::vector< std::complex< double > > > &a, std::vector< int > &ir)
Definition: Numerics.cc:424
double BesselK1L(const double xx)
Definition: Numerics.hh:73
void Dfeqn(const int n, std::vector< std::vector< double > > &a, std::vector< int > &ir, std::vector< double > &b)
Definition: Numerics.cc:98