RooGaussKronrodIntegrator1D.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooGaussKronrodIntegrator1D.cxx 28259 2009-04-16 16:21:16Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 //
00019 // BEGIN_HTML
00020 // RooGaussKronrodIntegrator1D implements the Gauss-Kronrod integration algorithm.
00021 //
00022 // An Gaussian quadrature method for numerical integration in which
00023 // error is estimation based on evaluation at special points known as
00024 // "Kronrod points."  By suitably picking these points, abscissas from
00025 // previous iterations can be reused as part of the new set of points,
00026 // whereas usual Gaussian quadrature would require recomputation of
00027 // all abscissas at each iteration.
00028 //
00029 // This class automatically handles (-inf,+inf) integrals by dividing
00030 // the integration in three regions (-inf,-1), (-1,1), (1,inf) and
00031 // calculating the 1st and 3rd term using a x -> 1/x coordinate
00032 // transformation
00033 //
00034 // This class embeds the Gauss-Kronrod integrator from the GNU
00035 // Scientific Library version 1.5 and applies the 10-, 21-, 43- and
00036 // 87-point rule in succession until the required target precision is
00037 // reached
00038 // END_HTML
00039 //
00040 
00041 
00042 
00043 #include "RooFit.h"
00044 
00045 #include <assert.h>
00046 #include <math.h>
00047 #include <float.h>
00048 #include <stdlib.h>
00049 #include "Riostream.h"
00050 #include "TMath.h"
00051 #include "RooGaussKronrodIntegrator1D.h"
00052 #include "RooArgSet.h"
00053 #include "RooRealVar.h"
00054 #include "RooNumber.h"
00055 #include "RooNumIntFactory.h"
00056 #include "RooIntegratorBinding.h"
00057 #include "RooMsgService.h"
00058 
00059 
00060 ClassImp(RooGaussKronrodIntegrator1D)
00061 ;
00062 
00063 
00064 // --- From GSL_MATH.h -------------------------------------------
00065 struct gsl_function_struct
00066 {
00067   double (* function) (double x, void * params);
00068   void * params;
00069 };
00070 typedef struct gsl_function_struct gsl_function ;
00071 #define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params)
00072 //----From GSL_INTEGRATION.h ---------------------------------------
00073 int gsl_integration_qng (const gsl_function * f,
00074                          double a, double b,
00075                          double epsabs, double epsrel,
00076                          double *result, double *abserr,
00077                          size_t * neval);
00078 //-------------------------------------------------------------------
00079 
00080 
00081 
00082 //_____________________________________________________________________________
00083 void RooGaussKronrodIntegrator1D::registerIntegrator(RooNumIntFactory& fact)
00084 {
00085   // Register RooGaussKronrodIntegrator1D, its parameters and capabilities with RooNumIntConfig
00086 
00087   fact.storeProtoIntegrator(new RooGaussKronrodIntegrator1D(),RooArgSet()) ;
00088 }
00089 
00090 
00091 
00092 //_____________________________________________________________________________
00093 RooGaussKronrodIntegrator1D::RooGaussKronrodIntegrator1D() : _x(0)
00094 {
00095   // Default constructor
00096 }
00097 
00098 
00099 
00100 //_____________________________________________________________________________
00101 RooGaussKronrodIntegrator1D::RooGaussKronrodIntegrator1D(const RooAbsFunc& function, const RooNumIntConfig& config) :
00102   RooAbsIntegrator(function),
00103   _epsAbs(config.epsRel()),
00104   _epsRel(config.epsAbs())
00105 {
00106   // Construct integral on 'function' using given configuration object. The integration
00107   // range is taken from the definition in the function binding
00108 
00109   _useIntegrandLimits= kTRUE;
00110   _valid= initialize();
00111 } 
00112 
00113 
00114 
00115 //_____________________________________________________________________________
00116 RooGaussKronrodIntegrator1D::RooGaussKronrodIntegrator1D(const RooAbsFunc& function, 
00117                                                          Double_t xmin, Double_t xmax, const RooNumIntConfig& config) :
00118   RooAbsIntegrator(function),
00119   _epsAbs(config.epsRel()),
00120   _epsRel(config.epsAbs()),
00121   _xmin(xmin),
00122   _xmax(xmax)
00123 {
00124   // Construct integral on 'function' using given configuration object in the given range
00125 
00126   _useIntegrandLimits= kFALSE;
00127   _valid= initialize();
00128 } 
00129 
00130 
00131 
00132 //_____________________________________________________________________________
00133 RooAbsIntegrator* RooGaussKronrodIntegrator1D::clone(const RooAbsFunc& function, const RooNumIntConfig& config) const
00134 {
00135   // Clone integrator with given function and configuration. Needed for RooNumIntFactory
00136 
00137   return new RooGaussKronrodIntegrator1D(function,config) ;
00138 }
00139 
00140 
00141 
00142 //_____________________________________________________________________________
00143 Bool_t RooGaussKronrodIntegrator1D::initialize()
00144 {
00145   // Perform one-time initialization of integrator
00146 
00147   // Allocate coordinate buffer size after number of function dimensions
00148   _x = new Double_t[_function->getDimension()] ;
00149 
00150   return checkLimits();
00151 }
00152 
00153 
00154 
00155 //_____________________________________________________________________________
00156 RooGaussKronrodIntegrator1D::~RooGaussKronrodIntegrator1D()
00157 {
00158   // Destructor
00159 
00160   if (_x) {
00161     delete[] _x ;
00162   }
00163 }
00164 
00165 
00166 
00167 //_____________________________________________________________________________
00168 Bool_t RooGaussKronrodIntegrator1D::setLimits(Double_t* xmin, Double_t* xmax) 
00169 {
00170   // Change our integration limits. Return kTRUE if the new limits are
00171   // ok, or otherwise kFALSE. Always returns kFALSE and does nothing
00172   // if this object was constructed to always use our integrand's limits.
00173 
00174   if(_useIntegrandLimits) {
00175     oocoutE((TObject*)0,Eval) << "RooGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
00176     return kFALSE;
00177   }
00178   _xmin= *xmin;
00179   _xmax= *xmax;
00180   return checkLimits();
00181 }
00182 
00183 
00184 
00185 //_____________________________________________________________________________
00186 Bool_t RooGaussKronrodIntegrator1D::checkLimits() const 
00187 {
00188   // Check that our integration range is finite and otherwise return kFALSE.
00189   // Update the limits from the integrand if requested.
00190 
00191   if(_useIntegrandLimits) {
00192     assert(0 != integrand() && integrand()->isValid());
00193     _xmin= integrand()->getMinLimit(0);
00194     _xmax= integrand()->getMaxLimit(0);
00195   }
00196   return kTRUE ;
00197 }
00198 
00199 
00200 
00201 double RooGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data) 
00202 {
00203   RooGaussKronrodIntegrator1D* instance = (RooGaussKronrodIntegrator1D*) data ;
00204   return instance->integrand(instance->xvec(x)) ;
00205 }
00206 
00207 
00208 
00209 //_____________________________________________________________________________
00210 Double_t RooGaussKronrodIntegrator1D::integral(const Double_t *yvec) 
00211 {
00212   // Calculate and return integral
00213 
00214   assert(isValid());
00215 
00216   // Copy yvec to xvec if provided
00217   if (yvec) {
00218     UInt_t i ; for (i=0 ; i<_function->getDimension()-1 ; i++) {
00219       _x[i+1] = yvec[i] ;
00220     }
00221   }
00222 
00223   // Setup glue function
00224   gsl_function F;
00225   F.function = &RooGaussKronrodIntegrator1D_GSL_GlueFunction ;
00226   F.params = this ;
00227 
00228   // Return values
00229   double result, error;
00230   size_t neval = 0 ;
00231 
00232   // Call GSL implementation of integeator
00233   gsl_integration_qng (&F, _xmin, _xmax, _epsAbs, _epsRel, &result, &error, &neval); 
00234 
00235   return result;
00236 }
00237 
00238 
00239 // ----------------------------------------------------------------------------
00240 // ---------- Code below imported from GSL version 1.5 ------------------------
00241 // ----------------------------------------------------------------------------
00242 
00243 /*
00244  * 
00245  * Copyright (C) 1996, 1997, 1998, 1999, 2000 Brian Gough
00246  * 
00247  * This program is free software; you can redistribute it and/or modify
00248  * it under the terms of the GNU General Public License as published by
00249  * the Free Software Foundation; either version 2 of the License, or (at
00250  * your option) any later version.
00251  * 
00252  * This program is distributed in the hope that it will be useful, but
00253  * WITHOUT ANY WARRANTY; without even the implied warranty of
00254  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00255  * General Public License for more details.
00256  * 
00257  * You should have received a copy of the GNU General Public License
00258  * along with this program; if not, write to the Free Software
00259  * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
00260  */
00261 
00262 #define GSL_SUCCESS 0
00263 #define GSL_EBADTOL 13  /* user specified an invalid tolerance */
00264 #define GSL_ETOL    14  /* failed to reach the specified tolerance */
00265 #define GSL_ERROR(a,b) oocoutE((TObject*)0,Eval) << "RooGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ;
00266 #define GSL_DBL_MIN        2.2250738585072014e-308
00267 #define GSL_DBL_EPSILON    2.2204460492503131e-16
00268 
00269 
00270 // INCLUDED BELOW #include "qng.c"
00271 
00272 int gsl_integration_qng (const gsl_function * f,
00273                          double a, double b,
00274                          double epsabs, double epsrel,
00275                          double *result, double *abserr,
00276                          size_t * neval);
00277 
00278 
00279 
00280 // INCLUDED BELOW #include "err.c"
00281 static double rescale_error (double err, const double result_abs, const double result_asc) ;
00282 
00283 static double
00284 rescale_error (double err, const double result_abs, const double result_asc)
00285 {
00286   err = fabs(err) ;
00287 
00288   if (result_asc != 0 && err != 0)
00289       {
00290         double scale = TMath::Power((200 * err / result_asc), 1.5) ;
00291         
00292         if (scale < 1)
00293           {
00294             err = result_asc * scale ;
00295           }
00296         else 
00297           {
00298             err = result_asc ;
00299           }
00300       }
00301   if (result_abs > GSL_DBL_MIN / (50 * GSL_DBL_EPSILON))
00302     {
00303       double min_err = 50 * GSL_DBL_EPSILON * result_abs ;
00304 
00305       if (min_err > err) 
00306         {
00307           err = min_err ;
00308         }
00309     }
00310   
00311   return err ;
00312 }
00313 
00314 
00315 // INCLUDED BELOW #include "qng.h"
00316 /* Gauss-Kronrod-Patterson quadrature coefficients for use in
00317    quadpack routine qng. These coefficients were calculated with
00318    101 decimal digit arithmetic by L. W. Fullerton, Bell Labs, Nov
00319    1981. */
00320 
00321 /* x1, abscissae common to the 10-, 21-, 43- and 87-point rule */
00322 static const double x1[5] = {
00323   0.973906528517171720077964012084452,
00324   0.865063366688984510732096688423493,
00325   0.679409568299024406234327365114874,
00326   0.433395394129247190799265943165784,
00327   0.148874338981631210884826001129720
00328 } ;
00329 
00330 /* w10, weights of the 10-point formula */
00331 static const double w10[5] = {
00332   0.066671344308688137593568809893332,
00333   0.149451349150580593145776339657697,
00334   0.219086362515982043995534934228163,
00335   0.269266719309996355091226921569469,
00336   0.295524224714752870173892994651338
00337 } ;
00338 
00339 /* x2, abscissae common to the 21-, 43- and 87-point rule */
00340 static const double x2[5] = {
00341   0.995657163025808080735527280689003,
00342   0.930157491355708226001207180059508,
00343   0.780817726586416897063717578345042,
00344   0.562757134668604683339000099272694,
00345   0.294392862701460198131126603103866
00346 } ;
00347 
00348 /* w21a, weights of the 21-point formula for abscissae x1 */
00349 static const double w21a[5] = {
00350   0.032558162307964727478818972459390,
00351   0.075039674810919952767043140916190,
00352   0.109387158802297641899210590325805,
00353   0.134709217311473325928054001771707,
00354   0.147739104901338491374841515972068
00355 } ;
00356 
00357 /* w21b, weights of the 21-point formula for abscissae x2 */
00358 static const double w21b[6] = {
00359   0.011694638867371874278064396062192,
00360   0.054755896574351996031381300244580,
00361   0.093125454583697605535065465083366,
00362   0.123491976262065851077958109831074,
00363   0.142775938577060080797094273138717,
00364   0.149445554002916905664936468389821
00365 } ;
00366 
00367 /* x3, abscissae common to the 43- and 87-point rule */
00368 static const double x3[11] = {
00369   0.999333360901932081394099323919911,
00370   0.987433402908088869795961478381209,
00371   0.954807934814266299257919200290473,
00372   0.900148695748328293625099494069092,
00373   0.825198314983114150847066732588520,
00374   0.732148388989304982612354848755461,
00375   0.622847970537725238641159120344323,
00376   0.499479574071056499952214885499755,
00377   0.364901661346580768043989548502644,
00378   0.222254919776601296498260928066212,
00379   0.074650617461383322043914435796506
00380 } ;
00381 
00382 /* w43a, weights of the 43-point formula for abscissae x1, x3 */
00383 static const double w43a[10] = {
00384   0.016296734289666564924281974617663,
00385   0.037522876120869501461613795898115,
00386   0.054694902058255442147212685465005,
00387   0.067355414609478086075553166302174,
00388   0.073870199632393953432140695251367,
00389   0.005768556059769796184184327908655,
00390   0.027371890593248842081276069289151,
00391   0.046560826910428830743339154433824,
00392   0.061744995201442564496240336030883,
00393   0.071387267268693397768559114425516
00394 } ;
00395 
00396 /* w43b, weights of the 43-point formula for abscissae x3 */
00397 static const double w43b[12] = {
00398   0.001844477640212414100389106552965,
00399   0.010798689585891651740465406741293,
00400   0.021895363867795428102523123075149,
00401   0.032597463975345689443882222526137,
00402   0.042163137935191811847627924327955,
00403   0.050741939600184577780189020092084,
00404   0.058379395542619248375475369330206,
00405   0.064746404951445885544689259517511,
00406   0.069566197912356484528633315038405,
00407   0.072824441471833208150939535192842,
00408   0.074507751014175118273571813842889,
00409   0.074722147517403005594425168280423
00410 } ;
00411 
00412 /* x4, abscissae of the 87-point rule */
00413 static const double x4[22] = {
00414   0.999902977262729234490529830591582,
00415   0.997989895986678745427496322365960,
00416   0.992175497860687222808523352251425,
00417   0.981358163572712773571916941623894,
00418   0.965057623858384619128284110607926,
00419   0.943167613133670596816416634507426,
00420   0.915806414685507209591826430720050,
00421   0.883221657771316501372117548744163,
00422   0.845710748462415666605902011504855,
00423   0.803557658035230982788739474980964,
00424   0.757005730685495558328942793432020,
00425   0.706273209787321819824094274740840,
00426   0.651589466501177922534422205016736,
00427   0.593223374057961088875273770349144,
00428   0.531493605970831932285268948562671,
00429   0.466763623042022844871966781659270,
00430   0.399424847859218804732101665817923,
00431   0.329874877106188288265053371824597,
00432   0.258503559202161551802280975429025,
00433   0.185695396568346652015917141167606,
00434   0.111842213179907468172398359241362,
00435   0.037352123394619870814998165437704
00436 } ;
00437 
00438 /* w87a, weights of the 87-point formula for abscissae x1, x2, x3 */
00439 static const double w87a[21] = {
00440   0.008148377384149172900002878448190,
00441   0.018761438201562822243935059003794,
00442   0.027347451050052286161582829741283,
00443   0.033677707311637930046581056957588,
00444   0.036935099820427907614589586742499,
00445   0.002884872430211530501334156248695,
00446   0.013685946022712701888950035273128,
00447   0.023280413502888311123409291030404,
00448   0.030872497611713358675466394126442,
00449   0.035693633639418770719351355457044,
00450   0.000915283345202241360843392549948,
00451   0.005399280219300471367738743391053,
00452   0.010947679601118931134327826856808,
00453   0.016298731696787335262665703223280,
00454   0.021081568889203835112433060188190,
00455   0.025370969769253827243467999831710,
00456   0.029189697756475752501446154084920,
00457   0.032373202467202789685788194889595,
00458   0.034783098950365142750781997949596,
00459   0.036412220731351787562801163687577,
00460   0.037253875503047708539592001191226
00461 } ;
00462 
00463 /* w87b, weights of the 87-point formula for abscissae x4    */
00464 static const double w87b[23] = {
00465   0.000274145563762072350016527092881,
00466   0.001807124155057942948341311753254,
00467   0.004096869282759164864458070683480,
00468   0.006758290051847378699816577897424,
00469   0.009549957672201646536053581325377,
00470   0.012329447652244853694626639963780,
00471   0.015010447346388952376697286041943,
00472   0.017548967986243191099665352925900,
00473   0.019938037786440888202278192730714,
00474   0.022194935961012286796332102959499,
00475   0.024339147126000805470360647041454,
00476   0.026374505414839207241503786552615,
00477   0.028286910788771200659968002987960,
00478   0.030052581128092695322521110347341,
00479   0.031646751371439929404586051078883,
00480   0.033050413419978503290785944862689,
00481   0.034255099704226061787082821046821,
00482   0.035262412660156681033782717998428,
00483   0.036076989622888701185500318003895,
00484   0.036698604498456094498018047441094,
00485   0.037120549269832576114119958413599,
00486   0.037334228751935040321235449094698,
00487   0.037361073762679023410321241766599
00488 } ;
00489 
00490 
00491 int
00492 gsl_integration_qng (const gsl_function *f,
00493                      double a, double b,
00494                      double epsabs, double epsrel,
00495                      double * result, double * abserr, size_t * neval)
00496 {
00497   double fv1[5], fv2[5], fv3[5], fv4[5];
00498   double savfun[21];  /* array of function values which have been computed */
00499   double res10, res21, res43, res87;    /* 10, 21, 43 and 87 point results */
00500   double result_kronrod, err ; 
00501   double resabs; /* approximation to the integral of abs(f) */
00502   double resasc; /* approximation to the integral of abs(f-i/(b-a)) */
00503 
00504   const double half_length =  0.5 * (b - a);
00505   const double abs_half_length = fabs (half_length);
00506   const double center = 0.5 * (b + a);
00507   const double f_center = GSL_FN_EVAL(f, center);
00508 
00509   int k ;
00510 
00511   if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28))
00512     {
00513       * result = 0;
00514       * abserr = 0;
00515       * neval = 0;
00516       GSL_ERROR ("tolerance cannot be acheived with given epsabs and epsrel",
00517                  GSL_EBADTOL);
00518     };
00519 
00520   /* Compute the integral using the 10- and 21-point formula. */
00521 
00522   res10 = 0;
00523   res21 = w21b[5] * f_center;
00524   resabs = w21b[5] * fabs (f_center);
00525 
00526   for (k = 0; k < 5; k++)
00527     {
00528       const double abscissa = half_length * x1[k];
00529       const double fval1 = GSL_FN_EVAL(f, center + abscissa);
00530       const double fval2 = GSL_FN_EVAL(f, center - abscissa);
00531       const double fval = fval1 + fval2;
00532       res10 += w10[k] * fval;
00533       res21 += w21a[k] * fval;
00534       resabs += w21a[k] * (fabs (fval1) + fabs (fval2));
00535       savfun[k] = fval;
00536       fv1[k] = fval1;
00537       fv2[k] = fval2;
00538     }
00539 
00540   for (k = 0; k < 5; k++)
00541     {
00542       const double abscissa = half_length * x2[k];
00543       const double fval1 = GSL_FN_EVAL(f, center + abscissa);
00544       const double fval2 = GSL_FN_EVAL(f, center - abscissa);
00545       const double fval = fval1 + fval2;
00546       res21 += w21b[k] * fval;
00547       resabs += w21b[k] * (fabs (fval1) + fabs (fval2));
00548       savfun[k + 5] = fval;
00549       fv3[k] = fval1;
00550       fv4[k] = fval2;
00551     }
00552 
00553   resabs *= abs_half_length ;
00554 
00555   { 
00556     const double mean = 0.5 * res21;
00557   
00558     resasc = w21b[5] * fabs (f_center - mean);
00559     
00560     for (k = 0; k < 5; k++)
00561       {
00562         resasc +=
00563           (w21a[k] * (fabs (fv1[k] - mean) + fabs (fv2[k] - mean))
00564           + w21b[k] * (fabs (fv3[k] - mean) + fabs (fv4[k] - mean)));
00565       }
00566     resasc *= abs_half_length ;
00567   }
00568 
00569   result_kronrod = res21 * half_length;
00570   
00571   err = rescale_error ((res21 - res10) * half_length, resabs, resasc) ;
00572 
00573   /*   test for convergence. */
00574 
00575   if (err < epsabs || err < epsrel * fabs (result_kronrod))
00576     {
00577       * result = result_kronrod ;
00578       * abserr = err ;
00579       * neval = 21;
00580       return GSL_SUCCESS;
00581     }
00582 
00583   /* compute the integral using the 43-point formula. */
00584 
00585   res43 = w43b[11] * f_center;
00586 
00587   for (k = 0; k < 10; k++)
00588     {
00589       res43 += savfun[k] * w43a[k];
00590     }
00591 
00592   for (k = 0; k < 11; k++)
00593     {
00594       const double abscissa = half_length * x3[k];
00595       const double fval = (GSL_FN_EVAL(f, center + abscissa) 
00596                            + GSL_FN_EVAL(f, center - abscissa));
00597       res43 += fval * w43b[k];
00598       savfun[k + 10] = fval;
00599     }
00600 
00601   /*  test for convergence */
00602 
00603   result_kronrod = res43 * half_length;
00604   err = rescale_error ((res43 - res21) * half_length, resabs, resasc);
00605 
00606   if (err < epsabs || err < epsrel * fabs (result_kronrod))
00607     {
00608       * result = result_kronrod ;
00609       * abserr = err ;
00610       * neval = 43;
00611       return GSL_SUCCESS;
00612     }
00613 
00614   /* compute the integral using the 87-point formula. */
00615 
00616   res87 = w87b[22] * f_center;
00617 
00618   for (k = 0; k < 21; k++)
00619     {
00620       res87 += savfun[k] * w87a[k];
00621     }
00622 
00623   for (k = 0; k < 22; k++)
00624     {
00625       const double abscissa = half_length * x4[k];
00626       res87 += w87b[k] * (GSL_FN_EVAL(f, center + abscissa) 
00627                           + GSL_FN_EVAL(f, center - abscissa));
00628     }
00629 
00630   /*  test for convergence */
00631 
00632   result_kronrod = res87 * half_length ;
00633   
00634   err = rescale_error ((res87 - res43) * half_length, resabs, resasc);
00635   
00636   if (err < epsabs || err < epsrel * fabs (result_kronrod))
00637     {
00638       * result = result_kronrod ;
00639       * abserr = err ;
00640       * neval = 87;
00641       return GSL_SUCCESS;
00642     }
00643 
00644   /* failed to converge */
00645 
00646   * result = result_kronrod ;
00647   * abserr = err ;
00648   * neval = 87;
00649 
00650   // GSL_ERROR("failed to reach tolerance with highest-order rule", GSL_ETOL) ;
00651   return GSL_ETOL ;
00652 }

Generated on Tue Jul 5 15:06:40 2011 for ROOT_528-00b_version by  doxygen 1.5.1