CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/HiggsAnalysis/CombinedLimit/src/VerticalInterpPdf.cc

Go to the documentation of this file.
00001 #include "../interface/VerticalInterpPdf.h"
00002 
00003 #include "RooFit.h"
00004 #include "Riostream.h"
00005 
00006 #include "TIterator.h"
00007 #include "TList.h"
00008 #include "RooRealProxy.h"
00009 #include "RooPlot.h"
00010 #include "RooRealVar.h"
00011 #include "RooAddGenContext.h"
00012 #include "RooRealConstant.h"
00013 #include "RooRealIntegral.h"
00014 #include "RooMsgService.h"
00015 
00016 
00017 
00018 ClassImp(VerticalInterpPdf)
00019 
00020 
00021 //_____________________________________________________________________________
00022 VerticalInterpPdf::VerticalInterpPdf() 
00023 {
00024   // Default constructor
00025   // coverity[UNINIT_CTOR]
00026   _funcIter  = _funcList.createIterator() ;
00027   _coefIter  = _coefList.createIterator() ;
00028   _quadraticRegion = 0;
00029 }
00030 
00031 
00032 //_____________________________________________________________________________
00033 VerticalInterpPdf::VerticalInterpPdf(const char *name, const char *title, const RooArgList& inFuncList, const RooArgList& inCoefList, Double_t quadraticRegion, Int_t quadraticAlgo) :
00034   RooAbsPdf(name,title),
00035   _normIntMgr(this,10),
00036   _funcList("!funcList","List of functions",this),
00037   _coefList("!coefList","List of coefficients",this),
00038   _quadraticRegion(quadraticRegion),
00039   _quadraticAlgo(quadraticAlgo)
00040 { 
00041 
00042   if (inFuncList.getSize()!=2*inCoefList.getSize()+1) {
00043     coutE(InputArguments) << "VerticalInterpPdf::VerticalInterpPdf(" << GetName() 
00044                           << ") number of pdfs and coefficients inconsistent, must have Nfunc=1+2*Ncoef" << endl ;
00045     assert(0);
00046   }
00047 
00048   TIterator* funcIter = inFuncList.createIterator() ;
00049   RooAbsArg* func;
00050   while((func = (RooAbsArg*)funcIter->Next())) {
00051     if (!dynamic_cast<RooAbsReal*>(func)) {
00052       coutE(InputArguments) << "ERROR: VerticalInterpPdf::VerticalInterpPdf(" << GetName() << ") function  " << func->GetName() << " is not of type RooAbsReal" << endl;
00053       assert(0);
00054     }
00055     _funcList.add(*func) ;
00056   }
00057   delete funcIter;
00058 
00059   TIterator* coefIter = inCoefList.createIterator() ;
00060   RooAbsArg* coef;
00061   while((coef = (RooAbsArg*)coefIter->Next())) {
00062     if (!dynamic_cast<RooAbsReal*>(coef)) {
00063       coutE(InputArguments) << "ERROR: VerticalInterpPdf::VerticalInterpPdf(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal" << endl;
00064       assert(0);
00065     }
00066     _coefList.add(*coef) ;    
00067   }
00068   delete coefIter;
00069 
00070   _funcIter  = _funcList.createIterator() ;
00071   _coefIter = _coefList.createIterator() ;
00072 
00073   if (_quadraticAlgo == -1) { 
00074     // multiplicative morphing: no way to do analytical integrals.
00075     _forceNumInt = kTRUE; 
00076   } else if (_quadraticAlgo >= 100) {
00077       _quadraticAlgo -= 100;
00078       _forceNumInt = kTRUE; 
00079   }
00080 }
00081 
00082 
00083 
00084 
00085 //_____________________________________________________________________________
00086 VerticalInterpPdf::VerticalInterpPdf(const VerticalInterpPdf& other, const char* name) :
00087   RooAbsPdf(other,name),
00088   _normIntMgr(other._normIntMgr,this),
00089   _funcList("!funcList",this,other._funcList),
00090   _coefList("!coefList",this,other._coefList),
00091   _quadraticRegion(other._quadraticRegion),
00092   _quadraticAlgo(other._quadraticAlgo)
00093 {
00094   // Copy constructor
00095 
00096   _funcIter  = _funcList.createIterator() ;
00097   _coefIter = _coefList.createIterator() ;
00098 }
00099 
00100 
00101 
00102 //_____________________________________________________________________________
00103 VerticalInterpPdf::~VerticalInterpPdf()
00104 {
00105   // Destructor
00106   delete _funcIter ;
00107   delete _coefIter ;
00108 }
00109 
00110 //_____________________________________________________________________________
00111 Double_t VerticalInterpPdf::evaluate() const 
00112 {
00113   // Calculate the current value
00114   Double_t value(0) ;
00115 
00116   // Do running sum of coef/func pairs, calculate lastCoef.
00117   _funcIter->Reset() ;
00118   _coefIter->Reset() ;
00119   RooAbsReal* coef ;
00120   RooAbsReal* func = (RooAbsReal*)_funcIter->Next();
00121 
00122   Double_t central = func->getVal();
00123   value = central;
00124 
00125   if (_quadraticAlgo >= 0) {
00126       // additive interpolation
00127       while((coef=(RooAbsReal*)_coefIter->Next())) {
00128           Double_t coefVal = coef->getVal() ;
00129           RooAbsReal* funcUp = (RooAbsReal*)_funcIter->Next() ;
00130           RooAbsReal* funcDn = (RooAbsReal*)_funcIter->Next() ;
00131           value += interpolate(coefVal, central, funcUp, funcDn);
00132       }
00133   } else {
00134       // multiplicative interpolation
00135       while((coef=(RooAbsReal*)_coefIter->Next())) {
00136           Double_t coefVal = coef->getVal() ;
00137           RooAbsReal* funcUp = (RooAbsReal*)_funcIter->Next() ;
00138           RooAbsReal* funcDn = (RooAbsReal*)_funcIter->Next() ;
00139           value *= interpolate(coefVal, central, funcUp, funcDn);
00140       }
00141   }
00142    
00143   return value > 0 ? value : 1E-9 ;
00144 }
00145 
00146 
00147 
00148 
00149 //_____________________________________________________________________________
00150 Bool_t VerticalInterpPdf::checkObservables(const RooArgSet* nset) const 
00151 {
00152   // Check if FUNC is valid for given normalization set.
00153   // Coeffient and FUNC must be non-overlapping, but func-coefficient 
00154   // pairs may overlap each other
00155   //
00156   // In the present implementation, coefficients may not be observables or derive
00157   // from observables
00158   if (_quadraticAlgo == -1) return false; // multiplicative morphing. we don't care.
00159   
00160 
00161   _coefIter->Reset() ;
00162   RooAbsReal* coef ;
00163   while((coef=(RooAbsReal*)_coefIter->Next())) {
00164     if (coef->dependsOn(*nset)) {
00165       coutE(InputArguments) << "RooRealPdf::checkObservables(" << GetName() << "): ERROR coefficient " << coef->GetName() 
00166                             << " depends on one or more of the following observables" ; nset->Print("1") ;
00167       return true;
00168     }
00169   }
00170 
00171   _funcIter->Reset() ;
00172   RooAbsReal* func ;
00173   while((func = (RooAbsReal*)_funcIter->Next())) { 
00174     if (func->observableOverlaps(nset,*coef)) {
00175       coutE(InputArguments) << "VerticalInterpPdf::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName() 
00176                             << " and FUNC " << func->GetName() << " have one or more observables in common" << endl ;
00177       return true;
00178     }
00179   }
00180   
00181   return false;
00182 }
00183 
00184 
00185 
00186 
00187 //_____________________________________________________________________________
00188 Int_t VerticalInterpPdf::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
00189                                              const RooArgSet* normSet2, const char* /*rangeName*/) const 
00190 {
00191   // Advertise that all integrals can be handled internally.
00192 
00193   // Handle trivial no-integration scenario
00194   if (allVars.getSize()==0) return 0 ;
00195   if (_forceNumInt) return 0 ;
00196 
00197   // Select subset of allVars that are actual dependents
00198   analVars.add(allVars) ;
00199   RooArgSet* normSet = normSet2 ? getObservables(normSet2) : 0 ;
00200 
00201 
00202   // Check if this configuration was created before
00203   Int_t sterileIdx(-1) ;
00204   CacheElem* cache = (CacheElem*) _normIntMgr.getObj(normSet,&analVars,&sterileIdx,(const char *)0);
00205   if (cache) {
00206     return _normIntMgr.lastIndex()+1 ;
00207   }
00208   
00209   // Create new cache element
00210   cache = new CacheElem ;
00211 
00212   // Make list of function projection and normalization integrals 
00213   _funcIter->Reset() ;
00214   RooAbsReal *func ;
00215   while((func=(RooAbsReal*)_funcIter->Next())) {
00216     RooAbsReal* funcInt = func->createIntegral(analVars) ;
00217     cache->_funcIntList.addOwned(*funcInt) ;
00218     if (normSet && normSet->getSize()>0) {
00219       RooAbsReal* funcNorm = func->createIntegral(*normSet) ;
00220       cache->_funcNormList.addOwned(*funcNorm) ;
00221     }
00222   }
00223 
00224   // Store cache element
00225   Int_t code = _normIntMgr.setObj(normSet,&analVars,(RooAbsCacheElement*)cache,0) ;
00226 
00227   if (normSet) {
00228     delete normSet ;
00229   }
00230 
00231   return code+1 ; 
00232 }
00233 
00234 
00235 
00236 
00237 //_____________________________________________________________________________
00238 Double_t VerticalInterpPdf::analyticalIntegralWN(Int_t code, const RooArgSet* normSet2, const char* /*rangeName*/) const 
00239 {
00240   // Implement analytical integrations by deferring integration of component
00241   // functions to integrators of components
00242 
00243   // Handle trivial passthrough scenario
00244   if (code==0) return getVal(normSet2) ;
00245 
00246   RooAbsReal *coef;
00247   Double_t value = 0;
00248 
00249   // WVE needs adaptation for rangeName feature
00250   CacheElem* cache = (CacheElem*) _normIntMgr.getObjByIndex(code-1) ;
00251 
00252   TIterator* funcIntIter = cache->_funcIntList.createIterator() ;
00253   RooAbsReal *funcInt = (RooAbsReal *) funcIntIter->Next();
00254   Double_t central = funcInt->getVal();
00255   value += central;
00256 
00257   _coefIter->Reset() ;
00258   while((coef=(RooAbsReal*)_coefIter->Next())) {
00259     Double_t coefVal = coef->getVal(normSet2) ;
00260     RooAbsReal * funcIntUp = (RooAbsReal*)funcIntIter->Next() ;
00261     RooAbsReal * funcIntDn = (RooAbsReal*)funcIntIter->Next() ;
00262     value += interpolate(coefVal, central, funcIntUp, funcIntDn);
00263   }
00264   
00265   delete funcIntIter ;
00266   
00267   Double_t normVal(1) ;
00268   if (normSet2) {
00269     normVal = 0 ;
00270 
00271     TIterator* funcNormIter = cache->_funcNormList.createIterator() ;
00272 
00273     RooAbsReal* funcNorm = (RooAbsReal*) funcNormIter->Next();
00274     central = funcNorm->getVal(normSet2) ;
00275 
00276     _coefIter->Reset() ;
00277     while((coef=(RooAbsReal*)_coefIter->Next())) {
00278       RooAbsReal *funcNormUp = (RooAbsReal*)funcNormIter->Next() ;
00279       RooAbsReal *funcNormDn = (RooAbsReal*)funcNormIter->Next() ;
00280       Double_t coefVal = coef->getVal(normSet2) ;
00281       normVal += interpolate(coefVal, central, funcNormUp, funcNormDn);
00282     }
00283     
00284     delete funcNormIter ;      
00285   }
00286 
00287   return ( value > 0 ? value : 1E-9 ) / normVal;
00288 }
00289 
00290 Double_t VerticalInterpPdf::interpolate(Double_t coeff, Double_t central, RooAbsReal *fUp, RooAbsReal *fDn) const  
00291 {
00292     if (_quadraticAlgo == -1) {
00293         Double_t kappa = (coeff > 0 ? fUp->getVal()/central : central/fDn->getVal());
00294         return pow(kappa,fabs(coeff));
00295     }
00296 
00297     if (fabs(coeff) >= _quadraticRegion) {
00298         return coeff * (coeff > 0 ? fUp->getVal() - central : central - fDn->getVal());
00299     } else {
00300         // quadratic interpolation coefficients between the three
00301         if (_quadraticAlgo != 1) {
00302             // quadratic interpolation null in zero and continuos at boundaries, but not differentiable at boundaries
00303             // conditions:
00304             //   c_up (+_quadraticRegion) = +_quadraticRegion
00305             //   c_cen(+_quadraticRegion) = -_quadraticRegion
00306             //   c_dn (+_quadraticRegion) = 0
00307             //   c_up (-_quadraticRegion) = 0 
00308             //   c_cen(-_quadraticRegion) = -_quadraticRegion
00309             //   c_dn (-_quadraticRegion) = +_quadraticRegion
00310             //   c_up(0) = c_dn(0) = c_cen(0) = 0
00311             Double_t c_up  = + coeff * (_quadraticRegion + coeff) / (2 * _quadraticRegion);
00312             Double_t c_dn  = - coeff * (_quadraticRegion - coeff) / (2 * _quadraticRegion);
00313             Double_t c_cen = - coeff * coeff / _quadraticRegion;
00314             return c_up * fUp->getVal() + c_dn * fDn->getVal() + c_cen * central;
00315         } else { //if (_quadraticAlgo == 1) { 
00316             // quadratic interpolation that is everywhere differentiable, but it's not null in zero
00317             // conditions on the function
00318             //   c_up (+_quadraticRegion) = +_quadraticRegion
00319             //   c_cen(+_quadraticRegion) = -_quadraticRegion
00320             //   c_dn (+_quadraticRegion) = 0
00321             //   c_up (-_quadraticRegion) = 0 
00322             //   c_cen(-_quadraticRegion) = -_quadraticRegion
00323             //   c_dn (-_quadraticRegion) = +_quadraticRegion
00324             // conditions on the derivatives
00325             //   c_up '(+_quadraticRegion) = +1
00326             //   c_cen'(+_quadraticRegion) = -1
00327             //   c_dn '(+_quadraticRegion) = 0
00328             //   c_up '(-_quadraticRegion) = 0
00329             //   c_cen'(-_quadraticRegion) = +1
00330             //   c_dn '(-_quadraticRegion) = -1
00331             Double_t c_up  = (_quadraticRegion + coeff) * (_quadraticRegion + coeff) / (4 * _quadraticRegion);
00332             Double_t c_dn  = (_quadraticRegion - coeff) * (_quadraticRegion - coeff) / (4 * _quadraticRegion);
00333             Double_t c_cen = - c_up - c_dn;
00334             return c_up * fUp->getVal() + c_dn * fDn->getVal() + c_cen * central;
00335         } 
00336     }
00337 
00338 }