RooParametricStepFunction.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitBabar                                                      *
00004  * @(#)root/roofit:$Id: RooParametricStepFunction.cxx 24286 2008-06-16 15:47:04Z wouter $
00005  * Authors:                                                                  *
00006  *    Aaron Roodman, Stanford Linear Accelerator Center, Stanford University *
00007  *                                                                           *
00008  * Copyright (c) 2004, Stanford University. All rights reserved.        *
00009  *           
00010  * Redistribution and use in source and binary forms,                        *
00011  * with or without modification, are permitted according to the terms        *
00012  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00013  *****************************************************************************/
00014 
00015 //////////////////////////////////////////////////////////////////////////////
00016 //
00017 // The Parametric Step Function PDF is a binned distribution whose parameters 
00018 // are the heights of each bin.  This PDF was first used in BaBar's B0->pi0pi0
00019 // paper BABAR Collaboration (B. Aubert et al.) Phys.Rev.Lett.91:241801,2003.
00020 // 
00021 // This PDF may be used to describe oddly shaped distributions.  It differs
00022 // from a RooKeysPdf or a RooHistPdf in that a RooParametricStepFunction
00023 // has free parameters.  In particular, any statistical uncertainty in 
00024 // sample used to model this PDF may be understood with these free parameters;
00025 // this is not possible with non-parametric PDFs.
00026 //
00027 // The RooParametricStepFunction has Nbins-1 free parameters. Note that
00028 // the limits of the dependent varaible must match the low and hi bin limits.
00029 //
00030 // An example of usage is:
00031 //
00032 // Int_t nbins(10);
00033 // TArrayD limits(nbins+1);
00034 // limits[0] = 0.0; //etc...
00035 // RooArgList* list = new RooArgList("list");
00036 // RooRealVar* binHeight0 = new RooRealVar("binHeight0","bin 0 Value",0.1,0.0,1.0);
00037 // list->add(binHeight0); // up to binHeight8, ie. 9 parameters
00038 //
00039 // RooParametricStepFunction  aPdf = ("aPdf","PSF",*x,
00040 //                                    *list,limits,nbins);
00041 //
00042                                 
00043 
00044 #include "RooFit.h"
00045 
00046 #include "Riostream.h"
00047 #include "TArrayD.h"
00048 #include <math.h>
00049 
00050 #include "RooParametricStepFunction.h"
00051 #include "RooAbsReal.h"
00052 #include "RooRealVar.h"
00053 #include "RooArgList.h"
00054 
00055 ClassImp(RooParametricStepFunction)
00056 ;
00057 
00058 
00059 //_____________________________________________________________________________
00060 RooParametricStepFunction::RooParametricStepFunction(const char* name, const char* title, 
00061                              RooAbsReal& x, const RooArgList& coefList, TArrayD& limits, Int_t nBins) :
00062   RooAbsPdf(name, title),
00063   _x("x", "Dependent", this, x),
00064   _coefList("coefList","List of coefficients",this),
00065   _nBins(nBins)
00066 {
00067   // Constructor
00068   _coefIter = _coefList.createIterator() ;
00069 
00070   // Check lowest order
00071   if (_nBins<0) {
00072     cout << "RooParametricStepFunction::ctor(" << GetName() 
00073          << ") WARNING: nBins must be >=0, setting value to 0" << endl ;
00074     _nBins=0 ;
00075   }
00076 
00077   TIterator* coefIter = coefList.createIterator() ;
00078   RooAbsArg* coef ;
00079   while((coef = (RooAbsArg*)coefIter->Next())) {
00080     if (!dynamic_cast<RooAbsReal*>(coef)) {
00081       cout << "RooParametricStepFunction::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName() 
00082            << " is not of type RooAbsReal" << endl ;
00083       assert(0) ;
00084     }
00085     _coefList.add(*coef) ;
00086   }
00087   delete coefIter ;
00088 
00089   // Bin limits  
00090   limits.Copy(_limits);
00091 
00092 }
00093 
00094 
00095 
00096 //_____________________________________________________________________________
00097 RooParametricStepFunction::RooParametricStepFunction(const RooParametricStepFunction& other, const char* name) :
00098   RooAbsPdf(other, name), 
00099   _x("x", this, other._x), 
00100   _coefList("coefList",this,other._coefList),
00101   _nBins(other._nBins)
00102 {
00103   // Copy constructor
00104   _coefIter = _coefList.createIterator();
00105   (other._limits).Copy(_limits);
00106 }
00107 
00108 
00109 
00110 //_____________________________________________________________________________
00111 RooParametricStepFunction::~RooParametricStepFunction()
00112 {
00113   // Destructor
00114   delete _coefIter ;
00115 }
00116 
00117 
00118 
00119 //_____________________________________________________________________________
00120 Int_t RooParametricStepFunction::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const 
00121 {
00122   if (matchArgs(allVars, analVars, _x)) return 1;
00123   return 0;
00124 }
00125 
00126 
00127 
00128 //_____________________________________________________________________________
00129 Double_t RooParametricStepFunction::analyticalIntegral(Int_t code, const char* rangeName) const 
00130 {
00131   assert(code==1) ;
00132   
00133   // Case without range is trivial: p.d.f is by construction normalized 
00134   if (!rangeName) {
00135     return 1.0 ;
00136   }
00137 
00138   // Case with ranges, calculate integral explicitly 
00139   Double_t xmin = _x.min(rangeName) ;
00140   Double_t xmax = _x.max(rangeName) ;
00141 
00142   Double_t sum=0 ;
00143   Int_t i ;
00144   for (i=1 ; i<=_nBins ; i++) {
00145     Double_t binVal = (i<_nBins) ? (static_cast<RooRealVar*>(_coefList.at(i-1))->getVal()) : lastBinValue() ;      
00146     if (_limits[i-1]>=xmin && _limits[i]<=xmax) {
00147       // Bin fully in the integration domain
00148       sum += (_limits[i]-_limits[i-1])*binVal ;
00149     } else if (_limits[i-1]<xmin && _limits[i]>xmax) {
00150       // Domain is fully contained in this bin
00151       sum += (xmax-xmin)*binVal ;
00152       // Exit here, this is the last bin to be processed by construction
00153       return sum ;
00154     } else if (_limits[i-1]<xmin && _limits[i]<=xmax && _limits[i]>xmin) {
00155       // Lower domain boundary is in bin
00156       sum +=  (_limits[i]-xmin)*binVal ;
00157     } else if (_limits[i-1]>=xmin && _limits[i]>xmax && _limits[i-1]<xmax) {
00158       sum +=  (xmax-_limits[i-1])*binVal ;      
00159       // Upper domain boundary is in bin
00160       // Exit here, this is the last bin to be processed by construction
00161       return sum ;
00162     }
00163   }
00164 
00165   return sum;  
00166 }
00167 
00168 
00169 
00170 //_____________________________________________________________________________
00171 Double_t RooParametricStepFunction::lastBinValue() const
00172 {
00173   Double_t sum(0.);
00174   Double_t binSize(0.);
00175   for (Int_t j=1;j<_nBins;j++){
00176     RooRealVar* tmp = (RooRealVar*) _coefList.at(j-1);
00177     binSize = _limits[j] - _limits[j-1];
00178     sum = sum + tmp->getVal()*binSize;
00179   }
00180   binSize = _limits[_nBins] - _limits[_nBins-1];
00181   return (1.0 - sum)/binSize;
00182 }
00183 
00184 
00185 
00186 //_____________________________________________________________________________
00187 Double_t RooParametricStepFunction::evaluate() const 
00188 {
00189   Double_t xval(0.);
00190   xval = _x;
00191   Double_t value(0.);
00192   if (_x >= _limits[0] && _x < _limits[_nBins]){
00193 
00194     for (Int_t i=1;i<=_nBins;i++){
00195       if (_x < _limits[i]){
00196         // in Bin i-1 (starting with Bin 0)
00197         if (i<_nBins) {
00198           // not in last Bin
00199           RooRealVar* tmp = (RooRealVar*) _coefList.at(i-1);
00200           value =  tmp->getVal();
00201           break;
00202         } else {
00203           // in last Bin
00204           Double_t sum(0.);
00205           Double_t binSize(0.);
00206           for (Int_t j=1;j<_nBins;j++){
00207             RooRealVar* tmp = (RooRealVar*) _coefList.at(j-1);
00208             binSize = _limits[j] - _limits[j-1];
00209             sum = sum + tmp->getVal()*binSize;
00210           }
00211           binSize = _limits[_nBins] - _limits[_nBins-1];
00212           value = (1.0 - sum)/binSize;
00213           if (value<=0.0){
00214             value = 0.000001;
00215             //      cout << "RooParametricStepFunction: sum of values gt 1.0 -- beware!!" <<endl;
00216           }
00217           break;
00218         }
00219       } 
00220     }
00221 
00222   } 
00223   return value;
00224 
00225 }
00226 
00227 
00228 //_____________________________________________________________________________
00229 Int_t RooParametricStepFunction::getnBins(){
00230   return _nBins;
00231 }
00232 
00233 
00234 //_____________________________________________________________________________
00235 Double_t* RooParametricStepFunction::getLimits(){
00236   Double_t* limoutput = _limits.GetArray();
00237   return limoutput;
00238 }

Generated on Tue Jul 5 14:55:24 2011 for ROOT_528-00b_version by  doxygen 1.5.1