00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include "Math/GaussIntegrator.h"
00012 #include <cmath>
00013
00014 namespace ROOT {
00015 namespace Math {
00016
00017
00018 bool GaussIntegrator::fgAbsValue = false;
00019
00020 GaussIntegrator::GaussIntegrator(double eps)
00021 {
00022
00023
00024 fEpsilon = eps;
00025 fLastResult = fLastError = 0;
00026 fUsedOnce = false;
00027 fFunction = 0;
00028 }
00029
00030 GaussIntegrator::~GaussIntegrator()
00031 {
00032
00033 }
00034
00035 void GaussIntegrator::AbsValue(bool flag)
00036 { fgAbsValue = flag; }
00037
00038 double GaussIntegrator::Integral(double a, double b) {
00039 return DoIntegral(a, b, fFunction);
00040 }
00041
00042 double GaussIntegrator::Integral () {
00043 IntegrandTransform it(this->fFunction);
00044 return DoIntegral(0., 1., it.Clone());
00045 }
00046
00047 double GaussIntegrator::IntegralUp (double a) {
00048 IntegrandTransform it(a, IntegrandTransform::kPlus, this->fFunction);
00049 return DoIntegral(0., 1., it.Clone());
00050 }
00051
00052 double GaussIntegrator::IntegralLow (double b) {
00053 IntegrandTransform it(b, IntegrandTransform::kMinus, this->fFunction);
00054 return DoIntegral(0., 1., it.Clone());
00055 }
00056
00057 double GaussIntegrator::DoIntegral(double a, double b, const IGenFunction* function)
00058 {
00059
00060
00061 const double kHF = 0.5;
00062 const double kCST = 5./1000;
00063
00064 double x[12] = { 0.96028985649753623, 0.79666647741362674,
00065 0.52553240991632899, 0.18343464249564980,
00066 0.98940093499164993, 0.94457502307323258,
00067 0.86563120238783174, 0.75540440835500303,
00068 0.61787624440264375, 0.45801677765722739,
00069 0.28160355077925891, 0.09501250983763744};
00070
00071 double w[12] = { 0.10122853629037626, 0.22238103445337447,
00072 0.31370664587788729, 0.36268378337836198,
00073 0.02715245941175409, 0.06225352393864789,
00074 0.09515851168249278, 0.12462897125553387,
00075 0.14959598881657673, 0.16915651939500254,
00076 0.18260341504492359, 0.18945061045506850};
00077
00078 double h, aconst, bb, aa, c1, c2, u, s8, s16, f1, f2;
00079 double xx[1];
00080 int i;
00081
00082 if ( fFunction == 0 )
00083 {
00084 MATH_ERROR_MSG("ROOT::Math::GausIntegratorOneDim", "A function must be set first!");
00085 return 0.0;
00086 }
00087
00088 h = 0;
00089 fUsedOnce = true;
00090 if (b == a) return h;
00091 aconst = kCST/std::abs(b-a);
00092 bb = a;
00093 CASE1:
00094 aa = bb;
00095 bb = b;
00096 CASE2:
00097 c1 = kHF*(bb+aa);
00098 c2 = kHF*(bb-aa);
00099 s8 = 0;
00100 for (i=0;i<4;i++) {
00101 u = c2*x[i];
00102 xx[0] = c1+u;
00103 f1 = (*function)(xx);
00104 if (fgAbsValue) f1 = std::abs(f1);
00105 xx[0] = c1-u;
00106 f2 = (*function) (xx);
00107 if (fgAbsValue) f2 = std::abs(f2);
00108 s8 += w[i]*(f1 + f2);
00109 }
00110 s16 = 0;
00111 for (i=4;i<12;i++) {
00112 u = c2*x[i];
00113 xx[0] = c1+u;
00114 f1 = (*function) (xx);
00115 if (fgAbsValue) f1 = std::abs(f1);
00116 xx[0] = c1-u;
00117 f2 = (*function) (xx);
00118 if (fgAbsValue) f2 = std::abs(f2);
00119 s16 += w[i]*(f1 + f2);
00120 }
00121 s16 = c2*s16;
00122 if (std::abs(s16-c2*s8) <= fEpsilon*(1. + std::abs(s16))) {
00123 h += s16;
00124 if(bb != b) goto CASE1;
00125 } else {
00126 bb = c1;
00127 if(1. + aconst*std::abs(c2) != 1) goto CASE2;
00128 h = s8;
00129 }
00130
00131 fLastResult = h;
00132 fLastError = std::abs(s16-c2*s8);
00133
00134 return h;
00135 }
00136
00137
00138 void GaussIntegrator::SetRelTolerance (double eps)
00139 { fEpsilon = eps; }
00140
00141 void GaussIntegrator::SetAbsTolerance (double)
00142 { MATH_WARN_MSG("ROOT::Math::GaussIntegrator", "There is no Absolute Tolerance!"); }
00143
00144 double GaussIntegrator::Result () const
00145 {
00146
00147
00148 if (!fUsedOnce)
00149 MATH_ERROR_MSG("ROOT::Math::GaussIntegrator", "You must calculate the result at least once!");
00150
00151 return fLastResult;
00152 }
00153
00154 double GaussIntegrator::Error() const
00155 { return fLastError; }
00156
00157 int GaussIntegrator::Status() const
00158 { return (fUsedOnce) ? 0 : -1; }
00159
00160 void GaussIntegrator::SetFunction (const IGenFunction & function)
00161 {
00162
00163 fFunction = &function;
00164
00165 fUsedOnce = false;
00166 }
00167
00168 double GaussIntegrator::Integral (const std::vector< double > &)
00169 {
00170
00171 MATH_WARN_MSG("ROOT::Math::GaussIntegrator", "This method is not implemented in this class !");
00172 return -1.0;
00173 }
00174
00175 double GaussIntegrator::IntegralCauchy (double , double , double )
00176 {
00177
00178 MATH_WARN_MSG("ROOT::Math::GaussIntegrator", "This method is not implemented in this class !");
00179 return -1.0;
00180 }
00181
00182 void GaussIntegrator::SetOptions(const ROOT::Math::IntegratorOneDimOptions & opt) {
00183
00184 SetRelTolerance(opt.RelTolerance() );
00185 }
00186
00187 ROOT::Math::IntegratorOneDimOptions GaussIntegrator::Options() const {
00188
00189 ROOT::Math::IntegratorOneDimOptions opt;
00190 opt.SetIntegrator("Gauss");
00191
00192 opt.SetRelTolerance(fEpsilon);
00193 opt.SetWKSize(0);
00194 opt.SetNPoints(0);
00195 return opt;
00196 }
00197
00198
00199
00200
00201
00202 IntegrandTransform::IntegrandTransform(const IGenFunction* integrand)
00203 : fSign(kPlus), fIntegrand(integrand), fBoundary(0.), fInfiniteInterval(true) {
00204 }
00205
00206 IntegrandTransform::IntegrandTransform(const double boundary, ESemiInfinitySign sign, const IGenFunction* integrand)
00207 : fSign(sign), fIntegrand(integrand), fBoundary(boundary), fInfiniteInterval(false) {
00208 }
00209
00210 double IntegrandTransform::DoEval(double x) const {
00211 double result = DoEval(x, fBoundary, fSign);
00212 return (result += (fInfiniteInterval ? DoEval(x, 0., -1) : 0.));
00213 }
00214
00215 double IntegrandTransform::DoEval(double x, double boundary, int sign) const {
00216 double mappedX = 1. / x - 1.;
00217 return (*fIntegrand)(boundary + sign * mappedX) * std::pow(mappedX + 1., 2);;
00218 }
00219
00220 double IntegrandTransform::operator()(double x) const {
00221 return DoEval(x);
00222 }
00223
00224 IGenFunction* IntegrandTransform::Clone() const {
00225 return (fInfiniteInterval ? new IntegrandTransform(fIntegrand) : new IntegrandTransform(fBoundary, fSign, fIntegrand));
00226 }
00227
00228
00229 }
00230 }