RooRealIntegral.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooRealIntegral.cxx 37380 2010-12-07 21:26:55Z 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 // RooRealIntegral performs hybrid numerical/analytical integrals of RooAbsReal objects
00021 // The class performs none of the actual integration, but only manages the logic
00022 // of what variables can be integrated analytically, accounts for eventual jacobian
00023 // terms and defines what numerical integrations needs to be done to complement the
00024 // analytical integral.
00025 // <p>
00026 // The actual analytical integrations (if any) are done in the PDF themselves, the numerical
00027 // integration is performed in the various implemenations of the RooAbsIntegrator base class.
00028 // END_HTML
00029 //
00030 
00031 #include "RooFit.h"
00032 
00033 #include "TClass.h"
00034 #include "RooMsgService.h"
00035 #include "Riostream.h"
00036 #include "TObjString.h"
00037 #include "TH1.h"
00038 #include "RooRealIntegral.h"
00039 #include "RooArgSet.h"
00040 #include "RooAbsRealLValue.h"
00041 #include "RooAbsCategoryLValue.h"
00042 #include "RooRealBinding.h"
00043 #include "RooRealAnalytic.h"
00044 #include "RooInvTransform.h"
00045 #include "RooSuperCategory.h"
00046 #include "RooNumIntFactory.h"
00047 #include "RooNumIntConfig.h"
00048 #include "RooNameReg.h"
00049 #include "RooExpensiveObjectCache.h"
00050 #include "RooConstVar.h"
00051 #include "RooDouble.h"
00052 
00053 ClassImp(RooRealIntegral) 
00054 ;
00055 
00056 
00057 Int_t RooRealIntegral::_cacheAllNDim(2) ;
00058 
00059 
00060 //_____________________________________________________________________________
00061 RooRealIntegral::RooRealIntegral() : 
00062   _valid(kFALSE),
00063   _funcNormSet(0),
00064   _iconfig(0),
00065   _sumCatIter(0),
00066   _numIntEngine(0),
00067   _numIntegrand(0),
00068   _rangeName(0),
00069   _params(0),
00070   _cacheNum(kFALSE)
00071 {
00072   _facListIter = _facList.createIterator() ;
00073   _jacListIter = _jacList.createIterator() ;
00074 }
00075 
00076 
00077 
00078 //_____________________________________________________________________________
00079 RooRealIntegral::RooRealIntegral(const char *name, const char *title, 
00080                                  const RooAbsReal& function, const RooArgSet& depList,
00081                                  const RooArgSet* funcNormSet, const RooNumIntConfig* config,
00082                                  const char* rangeName) :
00083   RooAbsReal(name,title), 
00084   _valid(kTRUE), 
00085   _sumList("!sumList","Categories to be summed numerically",this,kFALSE,kFALSE), 
00086   _intList("!intList","Variables to be integrated numerically",this,kFALSE,kFALSE), 
00087   _anaList("!anaList","Variables to be integrated analytically",this,kFALSE,kFALSE), 
00088   _jacList("!jacList","Jacobian product term",this,kFALSE,kFALSE), 
00089   _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
00090   _facListIter(_facList.createIterator()),
00091   _jacListIter(_jacList.createIterator()),
00092   _function("!func","Function to be integrated",this,
00093             const_cast<RooAbsReal&>(function),kFALSE,kFALSE), 
00094   _iconfig((RooNumIntConfig*)config),
00095   _sumCat("!sumCat","SuperCategory for summation",this,kFALSE,kFALSE),
00096   _sumCatIter(0),
00097   _mode(0),
00098   _intOperMode(Hybrid), 
00099   _restartNumIntEngine(kFALSE),
00100   _numIntEngine(0), 
00101   _numIntegrand(0),
00102   _rangeName((TNamed*)RooNameReg::ptr(rangeName)),
00103   _params(0),
00104   _cacheNum(kFALSE)
00105 {
00106   // Construct integral of 'function' over observables in 'depList'
00107   // in range 'rangeName'  with normalization observables 'funcNormSet' 
00108   // (for p.d.f.s). In the integral is performed to the maximum extent
00109   // possible the internal (analytical) integrals advertised by function.
00110   // The other integrations are performed numerically. The optional
00111   // config object prescribes how these numeric integrations are configured.
00112   //
00113 
00114   //   A) Check that all dependents are lvalues 
00115   //
00116   //   B) Check if list of dependents can be re-expressed in        
00117   //      lvalues that are higher in the expression tree            
00118   //
00119   //   C) Check for dependents that the PDF insists on integrating  
00120   //      analytically iself                                        
00121   //
00122   //   D) Make list of servers that can be integrated analytically  
00123   //      Add all parameters/dependents as value/shape servers      
00124   //
00125   //   E) Interact with function to make list of objects actually integrated analytically   
00126   //
00127   //   F) Make list of numerical integration variables consisting of:               
00128   //     - Category dependents of RealLValues in analytical integration             
00129   //     - Leaf nodes server lists of function server that are not analytically integrated   
00130   //     - Make Jacobian list for analytically integrated RealLValues            
00131   //
00132   //   G) Split numeric list in integration list and summation list   
00133   //
00134 
00135   oocxcoutI(&function,Integration) << "RooRealIntegral::ctor(" << GetName() << ") Constructing integral of function " 
00136                                      << function.GetName() << " over observables" << depList << " with normalization " 
00137                                      << (funcNormSet?*funcNormSet:RooArgSet()) << " with range identifier " 
00138                                      << (rangeName?rangeName:"<none>") << endl ;
00139 
00140   
00141   // Choose same expensive object cache as integrand
00142   setExpensiveObjectCache(function.expensiveObjectCache()) ;
00143 //   cout << "RRI::ctor(" << GetName() << ") setting expensive object cache to " << &expensiveObjectCache() << " as taken from " << function.GetName() << endl ;
00144 
00145   // Use objects integrator configuration if none is specified
00146   if (!_iconfig) _iconfig = (RooNumIntConfig*) function.getIntegratorConfig() ;
00147 
00148   // Save private copy of funcNormSet, if supplied, excluding factorizing terms
00149   if (funcNormSet) {
00150     _funcNormSet = new RooArgSet ;
00151     TIterator* iter = funcNormSet->createIterator() ;
00152     RooAbsArg* nArg ;  
00153     while ((nArg=(RooAbsArg*)iter->Next())) {
00154       if (function.dependsOn(*nArg)) {
00155         _funcNormSet->addClone(*nArg) ;
00156       }
00157     }
00158     delete iter ;
00159   } else {
00160     _funcNormSet = 0 ;
00161   }
00162 
00163   //_funcNormSet = funcNormSet ? (RooArgSet*)funcNormSet->snapshot(kFALSE) : 0 ;
00164   
00165   // Make internal copy of dependent list
00166   RooArgSet intDepList(depList) ;
00167 
00168   RooAbsArg *arg ;
00169 
00170   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00171   // * A) Check that all dependents are lvalues and filter out any
00172   //      dependents that the PDF doesn't explicitly depend on
00173   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00174   
00175   TIterator* depIter = intDepList.createIterator() ;
00176   while((arg=(RooAbsArg*)depIter->Next())) {
00177     if(!arg->isLValue()) {
00178       coutE(InputArguments) << ClassName() << "::" << GetName() << ": cannot integrate non-lvalue ";
00179       arg->Print("1");
00180       _valid= kFALSE;
00181     }
00182     if (!function.dependsOn(*arg)) {
00183       RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
00184       _facListOwned.addOwned(*argClone) ;
00185       _facList.add(*argClone) ;
00186       addServer(*argClone,kFALSE,kTRUE) ;
00187     }
00188   }
00189 
00190   if (_facList.getSize()>0) {
00191     oocxcoutI(&function,Integration) << function.GetName() << ": Factorizing obserables are " << _facList << endl ;
00192   }
00193     
00194 
00195   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00196   // * B) Check if list of dependents can be re-expressed in       *
00197   // *    lvalues that are higher in the expression tree           *
00198   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * 
00199 
00200 
00201   // Initial fill of list of LValue branches
00202   RooArgSet exclLVBranches("exclLVBranches") ;
00203   RooArgSet branchList,branchListVD ;
00204   function.branchNodeServerList(&branchList) ;
00205 
00206   TIterator* bIter = branchList.createIterator() ;
00207   RooAbsArg* branch ;
00208   while((branch=(RooAbsArg*)bIter->Next())) {
00209     RooAbsRealLValue    *realArgLV = dynamic_cast<RooAbsRealLValue*>(branch) ;
00210     RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(branch) ;
00211     if ((realArgLV && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {
00212       exclLVBranches.add(*branch) ;
00213 //       cout << "exclv branch = " << endl ;
00214 //       branch->printCompactTree() ;
00215     }
00216     if (dependsOnValue(*branch)) {
00217       branchListVD.add(*branch) ;
00218     } else {
00219 //       cout << "value of self does not depend on branch " << branch->GetName() << endl ;
00220     }
00221   }
00222   delete bIter ;
00223   exclLVBranches.remove(depList,kTRUE,kTRUE) ;
00224 //    cout << "exclLVBranches = " << exclLVBranches << endl ;
00225 
00226   // Initial fill of list of LValue leaf servers (put in intDepList)
00227   RooArgSet exclLVServers("exclLVServers") ;
00228   exclLVServers.add(intDepList) ;
00229 
00230 //    cout << "begin exclLVServers = " << exclLVServers << endl ;
00231   
00232   // Obtain mutual exclusive dependence by iterative reduction
00233   TIterator *sIter = exclLVServers.createIterator() ;
00234   bIter = exclLVBranches.createIterator() ;
00235   RooAbsArg *server ;
00236   Bool_t converged(kFALSE) ;
00237   while(!converged) {
00238     converged=kTRUE ;
00239 
00240     // Reduce exclLVServers to only those serving exclusively exclLVBranches
00241     sIter->Reset() ;
00242     while ((server=(RooAbsArg*)sIter->Next())) {
00243       if (!servesExclusively(server,exclLVBranches,branchListVD)) {
00244         exclLVServers.remove(*server) ;
00245 //      cout << "removing " << server->GetName() << " from exclLVServers because servesExclusively(" << server->GetName() << "," << exclLVBranches << ") faile" << endl ;
00246         converged=kFALSE ;
00247       }
00248     }
00249     
00250     // Reduce exclLVBranches to only those depending exclusisvely on exclLVservers
00251     bIter->Reset() ;
00252     while((branch=(RooAbsArg*)bIter->Next())) {
00253       RooArgSet* brDepList = branch->getObservables(&intDepList) ;
00254       RooArgSet bsList(*brDepList,"bsList") ;
00255       delete brDepList ;
00256       bsList.remove(exclLVServers,kTRUE,kTRUE) ;
00257       if (bsList.getSize()>0) {
00258         exclLVBranches.remove(*branch,kTRUE,kTRUE) ;
00259 //      cout << "removing " << branch->GetName() << " from exclLVBranches" << endl ;
00260         converged=kFALSE ;
00261       }
00262     }
00263   }
00264   delete sIter ;
00265   delete bIter ;
00266 
00267 //   cout << "end exclLVServers = " << exclLVServers << endl ;
00268      
00269   // Replace exclusive lvalue branch servers with lvalue branches
00270   if (exclLVServers.getSize()>0) {
00271 //     cout << "activating LVservers " << exclLVServers << " for use in integration " << endl ;
00272     intDepList.remove(exclLVServers) ;
00273     intDepList.add(exclLVBranches) ;
00274   }
00275 
00276      
00277   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00278   // * C) Check for dependents that the PDF insists on integrating *
00279   //      analytically iself                                       *
00280   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00281 
00282   RooArgSet anIntOKDepList ;
00283   depIter->Reset() ;
00284   while((arg=(RooAbsArg*)depIter->Next())) {
00285     if (function.forceAnalyticalInt(*arg)) {
00286       anIntOKDepList.add(*arg) ;
00287     }
00288   }
00289   delete depIter ;
00290   
00291   if (anIntOKDepList.getSize()>0) {
00292     oocxcoutI(&function,Integration) << function.GetName() << ": Observables that function forcibly requires to be integrated internally " << anIntOKDepList << endl ;
00293   }
00294 
00295 
00296   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00297   // * D) Make list of servers that can be integrated analytically *
00298   //      Add all parameters/dependents as value/shape servers     *
00299   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00300   
00301   sIter = function.serverIterator() ;
00302   while((arg=(RooAbsArg*)sIter->Next())) {
00303 
00304     //cout << "considering server" << arg->GetName() << endl ;
00305 
00306     // Dependent or parameter?
00307     if (!arg->dependsOnValue(intDepList)) {
00308 
00309       //cout << " server does not depend on observables, adding server as value server to integral" << endl ;
00310 
00311       if (function.dependsOnValue(*arg)) {
00312         addServer(*arg,kTRUE,kFALSE) ;
00313       }
00314 
00315       continue ;
00316 
00317     } else {
00318 
00319       // Add final dependents of arg as shape servers
00320       RooArgSet argLeafServers ;
00321       arg->leafNodeServerList(&argLeafServers,0,kFALSE) ;
00322 
00323       //arg->printCompactTree() ;
00324       //cout << "leaf nodes of server are " << argLeafServers << " depList = " << depList << endl ;
00325 
00326       // Skip arg if it is neither value or shape server
00327       if (!arg->isValueServer(function) && !arg->isShapeServer(function)) {
00328         //cout << " server is neither value not shape server of function, ignoring" << endl ;
00329         continue ;
00330       }
00331       
00332       TIterator* lIter = argLeafServers.createIterator() ;
00333       RooAbsArg* leaf ;
00334       while((leaf=(RooAbsArg*)lIter->Next())) {
00335 
00336         //cout << " considering leafnode " << leaf->GetName() << " of server " << arg->GetName() << endl ;
00337 
00338         if (depList.find(leaf->GetName()) && function.dependsOnValue(*leaf)) {
00339 
00340           RooAbsRealLValue* leaflv = dynamic_cast<RooAbsRealLValue*>(leaf) ;
00341           if (leaflv && leaflv->getBinning(rangeName).isParameterized()) {
00342             oocxcoutD(&function,Integration) << function.GetName() << " : Observable " << leaf->GetName() << " has parameterized binning, add value dependence of boundary objects rather than shape of leaf" << endl ;
00343             addServer(*leaflv->getBinning(rangeName).lowBoundFunc(),kTRUE,kFALSE) ;
00344             addServer(*leaflv->getBinning(rangeName).highBoundFunc(),kTRUE,kFALSE) ;
00345           } else {
00346             oocxcoutD(&function,Integration) << function.GetName() << ": Adding observable " << leaf->GetName() << " of server " 
00347                                              << arg->GetName() << " as shape dependent" << endl ;
00348             addServer(*leaf,kFALSE,kTRUE) ;
00349           }
00350         } else if (!depList.find(leaf->GetName())) {
00351 
00352           if (function.dependsOnValue(*leaf)) {
00353             oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as value dependent" << endl ;
00354             addServer(*leaf,kTRUE,kFALSE) ;
00355           } else {
00356             oocxcoutD(&function,Integration) << function.GetName() << ": Adding parameter " << leaf->GetName() << " of server " << arg->GetName() << " as shape dependent" << endl ;
00357             addServer(*leaf,kFALSE,kTRUE) ;
00358           }
00359         }       
00360       }
00361       delete lIter ;
00362     }
00363 
00364     // If this dependent arg is self-normalized, stop here
00365     //if (function.selfNormalized()) continue ;
00366 
00367     Bool_t depOK(kFALSE) ;
00368     // Check for integratable AbsRealLValue
00369     if (arg->isDerived()) {
00370       RooAbsRealLValue    *realArgLV = dynamic_cast<RooAbsRealLValue*>(arg) ;
00371       RooAbsCategoryLValue *catArgLV = dynamic_cast<RooAbsCategoryLValue*>(arg) ;
00372 //        cout << "realArgLV = " << realArgLV << " intDepList = " << intDepList << endl ;
00373       if ((realArgLV && intDepList.find(realArgLV->GetName()) && (realArgLV->isJacobianOK(intDepList)!=0)) || catArgLV) {       
00374 
00375         //cout  << " arg " << arg->GetName() << " is derived LValue with valid jacobian" << endl ;
00376 
00377         // Derived LValue with valid jacobian
00378         depOK = kTRUE ;
00379         
00380         // Now, check for overlaps
00381         Bool_t overlapOK = kTRUE ;
00382         RooAbsArg *otherArg ;
00383         TIterator* sIter2 = function.serverIterator() ; 
00384         while((otherArg=(RooAbsArg*)sIter2->Next())) {
00385           // skip comparison with self
00386           if (arg==otherArg) continue ;
00387           if (otherArg->IsA()==RooConstVar::Class()) continue ;
00388           if (arg->overlaps(*otherArg,kTRUE)) {
00389             //cout << "arg " << arg->GetName() << " overlaps with " << otherArg->GetName() << endl ;
00390             //overlapOK=kFALSE ;
00391           }
00392         }       
00393         // coverity[DEADCODE]
00394         if (!overlapOK) depOK=kFALSE ;      
00395 
00396         //cout << "overlap check returns OK=" << (depOK?"T":"F") << endl ;
00397 
00398         delete sIter2 ;
00399       }
00400     } else {
00401       // Fundamental types are always OK
00402       depOK = kTRUE ;
00403     }
00404     
00405     // Add server to list of dependents that are OK for analytical integration
00406     if (depOK) {
00407       anIntOKDepList.add(*arg,kTRUE) ;      
00408       oocxcoutI(&function,Integration) << function.GetName() << ": Observable " << arg->GetName() << " is suitable for analytical integration (if supported by p.d.f)" << endl ;
00409     }
00410   }
00411   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00412   // * E) interact with function to make list of objects actually integrated analytically  *
00413   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00414 
00415   RooArgSet anIntDepList ;
00416 
00417   RooArgSet *anaSet = new RooArgSet( _anaList, Form("UniqueCloneOf_%s",_anaList.GetName()));
00418   _mode = ((RooAbsReal&)_function.arg()).getAnalyticalIntegralWN(anIntOKDepList,*anaSet,_funcNormSet,RooNameReg::str(_rangeName)) ;    
00419   _anaList.removeAll() ;
00420   _anaList.add(*anaSet);    
00421   delete anaSet;
00422 
00423   // Avoid confusion -- if mode is zero no analytical integral is defined regardless of contents of _anaListx
00424   if (_mode==0) {
00425     _anaList.removeAll() ;
00426   }
00427 
00428   if (_mode!=0) {
00429     oocxcoutI(&function,Integration) << function.GetName() << ": Function integrated observables " << _anaList << " internally with code " << _mode << endl ;
00430   }
00431 
00432 
00433   // WVE kludge: synchronize dset for use in analyticalIntegral
00434   function.getVal(_funcNormSet) ;
00435 
00436   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00437   // * F) Make list of numerical integration variables consisting of:            *  
00438   // *   - Category dependents of RealLValues in analytical integration          *  
00439   // *   - Expanded server lists of server that are not analytically integrated  *
00440   // *    Make Jacobian list with analytically integrated RealLValues            *
00441   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00442 
00443   RooArgSet numIntDepList ;
00444 
00445   // Loop over actually analytically integrated dependents
00446   TIterator* aiIter = _anaList.createIterator() ;
00447   while ((arg=(RooAbsArg*)aiIter->Next())) {    
00448 
00449     // Process only derived RealLValues
00450     if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class()) && arg->isDerived() && !arg->isFundamental()) {
00451 
00452       // Add to list of Jacobians to calculate
00453       _jacList.add(*arg) ;
00454 
00455       // Add category dependent of LValueReal used in integration
00456       RooAbsArg *argDep ;
00457       RooArgSet *argDepList = arg->getObservables(&intDepList) ;
00458       TIterator *adIter = argDepList->createIterator() ;
00459       while ((argDep=(RooAbsArg*)adIter->Next())) {
00460         if (argDep->IsA()->InheritsFrom(RooAbsCategoryLValue::Class()) && intDepList.contains(*argDep)) {
00461           numIntDepList.add(*argDep,kTRUE) ;
00462         }
00463       }
00464       delete adIter ;
00465       delete argDepList ;
00466     }
00467   }
00468   delete aiIter ;
00469 
00470   // Loop again over function servers to add remaining numeric integrations
00471   sIter->Reset() ;
00472   while((arg=(RooAbsArg*)sIter->Next())) {
00473 
00474     // Process only servers that are not treated analytically
00475     if (!_anaList.find(arg->GetName()) && arg->dependsOn(intDepList)) {
00476 
00477       // Process only derived RealLValues
00478       if (dynamic_cast<RooAbsLValue*>(arg) && arg->isDerived() && intDepList.contains(*arg)) {
00479         numIntDepList.add(*arg,kTRUE) ; 
00480       } else {
00481         
00482         // Expand server in final dependents 
00483         RooArgSet *argDeps = arg->getObservables(&intDepList) ;
00484 
00485         // Add final dependents, that are not forcibly integrated analytically, 
00486         // to numerical integration list      
00487         TIterator* iter = argDeps->createIterator() ;
00488         RooAbsArg* dep ;
00489         while((dep=(RooAbsArg*)iter->Next())) {
00490           if (!_anaList.find(dep->GetName())) {
00491             numIntDepList.add(*dep,kTRUE) ;
00492           }
00493         }      
00494         delete iter ;
00495         delete argDeps ; 
00496       }
00497 
00498     }
00499   }
00500   delete sIter ;
00501 
00502   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00503   // * G) Split numeric list in integration list and summation list  *
00504   // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
00505 
00506   // Split numeric integration list in summation and integration lists
00507   TIterator* numIter=numIntDepList.createIterator() ;
00508   while ((arg=(RooAbsArg*)numIter->Next())) {
00509 
00510     if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
00511       _intList.add(*arg) ;
00512     } else if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
00513       _sumList.add(*arg) ;
00514     }
00515   }
00516   delete numIter ;
00517 
00518   if (_anaList.getSize()>0) {
00519     oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _anaList << " are analytically integrated with code " << _mode << endl ;
00520   }
00521   if (_intList.getSize()>0) {
00522     oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _intList << " are numerically integrated" << endl ;
00523   }
00524   if (_sumList.getSize()>0) {
00525     oocxcoutI(&function,Integration) << function.GetName() << ": Observables " << _sumList << " are numerically summed" << endl ;
00526   }
00527   
00528 
00529   // Determine operating mode
00530   if (numIntDepList.getSize()>0) {
00531     // Numerical and optional Analytical integration
00532     _intOperMode = Hybrid ;
00533   } else if (_anaList.getSize()>0) {
00534     // Purely analytical integration
00535     _intOperMode = Analytic ;    
00536   } else {
00537     // No integration performed
00538     _intOperMode = PassThrough ;
00539   }
00540 
00541   // Determine auto-dirty status  
00542   autoSelectDirtyMode() ;
00543 
00544   // Create value caches for _intList and _sumList
00545   _intList.snapshot(_saveInt) ;
00546   _sumList.snapshot(_saveSum) ;
00547 
00548   
00549   if (_sumList.getSize()>0) {
00550     RooSuperCategory *sumCat = new RooSuperCategory(Form("%s_sumCat",GetName()),"sumCat",_sumList) ;
00551     _sumCatIter = sumCat->typeIterator() ;    
00552     _sumCat.addOwned(*sumCat) ;
00553   }
00554 
00555 }
00556 
00557 
00558 
00559 //_____________________________________________________________________________
00560 void RooRealIntegral::autoSelectDirtyMode() 
00561 {
00562   // Set appropriate cache operation mode for integral depending on cache operation
00563   // mode of server objects
00564 
00565   //cout << "RooRealIntegral::autoSelectDirtyMode(" << GetName() << ")" << endl ;
00566 
00567   // If any of our servers are is forcedDirty or a projectedDependent, then we need to be ADirty
00568   TIterator* siter = serverIterator() ;  
00569   RooAbsArg* server ;
00570   while((server=(RooAbsArg*)siter->Next())){
00571     RooArgSet leafSet ;
00572     server->leafNodeServerList(&leafSet) ;
00573     TIterator* liter = leafSet.createIterator() ;
00574     RooAbsArg* leaf ;
00575     while((leaf=(RooAbsArg*)liter->Next())) {
00576       if (leaf->operMode()==ADirty && leaf->isValueServer(*this)) {      
00577         //cout << "RooRealIntegral::autoSelectDirtyMode(" << GetName() << ") selecting ADirty mode because value server leaf " 
00578         //     << leaf->GetName() << " is also " << endl ;
00579         setOperMode(ADirty) ;
00580         break ;
00581       }
00582       if (leaf->getAttribute("projectedDependent")) {
00583         //cout << "RooRealIntegral::autoSelectDirtyMode(" << GetName() << ") selecting ADirty mode because leaf " 
00584         //    << leaf->GetName() << " is projectedDependent " << endl ;
00585         setOperMode(ADirty) ;
00586         break ;
00587       }
00588     }
00589     delete liter ;
00590   }
00591   delete siter ;
00592 }
00593 
00594 
00595 
00596 //_____________________________________________________________________________
00597 Bool_t RooRealIntegral::servesExclusively(const RooAbsArg* server,const RooArgSet& exclLVBranches, const RooArgSet& allBranches) const
00598 {  
00599   // Utility function that returns true if 'object server' is a server
00600   // to exactly one of the RooAbsArgs in 'exclLVBranches'
00601 
00602   // Determine if given server serves exclusively exactly one of the given nodes in exclLVBranches
00603 
00604   // Special case, no LV servers available
00605   if (exclLVBranches.getSize()==0) return kFALSE ;
00606 
00607   // If server has no clients and is not an LValue itself, return false
00608    if (server->_clientList.GetSize()==0 && exclLVBranches.find(server->GetName())) {
00609      return kFALSE ;
00610    }
00611 
00612    // WVE must check for value relations only here!!!!
00613 
00614 
00615 //    cout << "servesExclusively: does " << server->GetName() << " serve only one of " << exclLVBranches << endl ;
00616 
00617    // Loop over all clients
00618    Int_t numLVServ(0) ;
00619    RooAbsArg* client ;
00620    TIterator* cIter = server->valueClientIterator() ;
00621    while((client=(RooAbsArg*)cIter->Next())) {
00622 //      cout << "now checking value client " << client->GetName() << " of server " << server->GetName() << endl ;
00623      // If client is not an LValue, recurse
00624      if (!(exclLVBranches.find(client->GetName())==client)) {
00625 //        cout << " client " << client->GetName() << "is not an lvalue" << endl ;
00626        if (allBranches.find(client->GetName())==client) {
00627 //       cout << " ... recursing call" << endl ;
00628          if (!servesExclusively(client,exclLVBranches,allBranches)) {
00629          // Client is a non-LValue that doesn't have an exclusive LValue server
00630          delete cIter ;
00631 //       cout << "client " << client->GetName() << " is a non-lvalue that doesn't have an exclusive lvalue server" << endl ;
00632          return kFALSE ;         
00633          }
00634        }
00635      } else {
00636        // Client is an LValue       
00637 //        cout << "client " << client->GetName() << " of server " << server->GetName() << " is an LValue " << endl ;
00638        numLVServ++ ;
00639      }
00640    }
00641 
00642    delete cIter ;
00643 //    cout << "numLVserv = " << numLVServ << endl ;
00644    return (numLVServ==1) ;
00645 }
00646 
00647 
00648 
00649 
00650 //_____________________________________________________________________________
00651 Bool_t RooRealIntegral::initNumIntegrator() const
00652 {
00653   // (Re)Initialize numerical integration engine if necessary. Return kTRUE if
00654   // successful, or otherwise kFALSE.
00655 
00656   // if we already have an engine, check if it still works for the present limits.
00657   if(0 != _numIntEngine) {
00658     if(_numIntEngine->isValid() && _numIntEngine->checkLimits() && !_restartNumIntEngine ) return kTRUE;
00659     // otherwise, cleanup the old engine
00660     delete _numIntEngine ;
00661     _numIntEngine= 0;
00662     if(0 != _numIntegrand) {
00663       delete _numIntegrand;
00664       _numIntegrand= 0;
00665     }
00666   }
00667 
00668   // All done if there are no arguments to integrate numerically
00669   if(0 == _intList.getSize()) return kTRUE;
00670   
00671   // Bind the appropriate analytic integral (specified by _mode) of our RooRealVar object to
00672   // those of its arguments that will be integrated out numerically.
00673   if(_mode != 0) {
00674     _numIntegrand= new RooRealAnalytic(_function.arg(),_intList,_mode,_funcNormSet,_rangeName);
00675   }
00676   else {
00677     _numIntegrand= new RooRealBinding(_function.arg(),_intList,_funcNormSet,kFALSE,_rangeName);
00678   }
00679   if(0 == _numIntegrand || !_numIntegrand->isValid()) {
00680     coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrand." << endl;
00681     return kFALSE;
00682   }
00683 
00684   // Create appropriate numeric integrator using factory
00685   _numIntEngine = RooNumIntFactory::instance().createIntegrator(*_numIntegrand,*_iconfig) ;
00686 
00687   if(0 == _numIntEngine || !_numIntEngine->isValid()) {
00688     coutE(Integration) << ClassName() << "::" << GetName() << ": failed to create valid integrator." << endl;
00689     return kFALSE;
00690   }
00691 
00692   cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") using numeric integrator " 
00693                           << _numIntEngine->IsA()->GetName() << " to calculate Int" << _intList << endl ;
00694 
00695   if (_intList.getSize()>3) {
00696     cxcoutI(NumIntegration) << "RooRealIntegral::init(" << GetName() << ") evaluation requires " << _intList.getSize() << "-D numeric integration step. Evaluation may be slow, sufficient numeric precision for fitting & minimization is not guaranteed" << endl ;
00697   }
00698 
00699   _restartNumIntEngine = kFALSE ;
00700   return kTRUE;
00701 }
00702 
00703 
00704 
00705 //_____________________________________________________________________________
00706 RooRealIntegral::RooRealIntegral(const RooRealIntegral& other, const char* name) : 
00707   RooAbsReal(other,name), 
00708   _valid(other._valid),
00709   _sumList("!sumList",this,other._sumList),
00710   _intList("!intList",this,other._intList), 
00711   _anaList("!anaList",this,other._anaList),
00712   _jacList("!jacList",this,other._jacList),
00713   _facList("!facList","Variables independent of function",this,kFALSE,kTRUE),
00714   _facListIter(_facList.createIterator()),
00715   _jacListIter(_jacList.createIterator()),
00716   _function("!func",this,other._function), 
00717   _iconfig(other._iconfig),
00718   _sumCat("!sumCat",this,other._sumCat),
00719   _sumCatIter(0),
00720   _mode(other._mode),
00721   _intOperMode(other._intOperMode), 
00722   _restartNumIntEngine(kFALSE),
00723   _numIntEngine(0), 
00724   _numIntegrand(0),
00725   _rangeName(other._rangeName),
00726   _params(0),
00727   _cacheNum(kFALSE)
00728 {
00729   // Copy constructor
00730 
00731  _funcNormSet = other._funcNormSet ? (RooArgSet*)other._funcNormSet->snapshot(kFALSE) : 0 ;
00732 
00733  other._facListIter->Reset() ;
00734  RooAbsArg* arg ;
00735  while((arg=(RooAbsArg*)other._facListIter->Next())) {
00736    RooAbsArg* argClone = (RooAbsArg*) arg->Clone() ;
00737    _facListOwned.addOwned(*argClone) ;
00738    _facList.add(*argClone) ;
00739    addServer(*argClone,kFALSE,kTRUE) ;
00740  }
00741 
00742  other._intList.snapshot(_saveInt) ;
00743  other._sumList.snapshot(_saveSum) ;
00744 
00745 }
00746 
00747 
00748 
00749 //_____________________________________________________________________________
00750 RooRealIntegral::~RooRealIntegral()
00751   // Destructor
00752 {
00753   if (_numIntEngine) delete _numIntEngine ;
00754   if (_numIntegrand) delete _numIntegrand ;
00755   if (_funcNormSet) delete _funcNormSet ;
00756   delete _facListIter ;
00757   delete _jacListIter ;
00758   if (_sumCatIter)  delete _sumCatIter ;
00759 }
00760 
00761 
00762 
00763 
00764 
00765 //_____________________________________________________________________________
00766 RooAbsReal* RooRealIntegral::createIntegral(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const 
00767 {
00768   // Special handling of integral of integral, return RooRealIntegral that represents integral over all dimensions in one pass
00769   RooArgSet isetAll(iset) ;
00770   isetAll.add(_sumList) ;
00771   isetAll.add(_intList) ;
00772   isetAll.add(_anaList) ;
00773   isetAll.add(_facList) ;
00774 
00775   const RooArgSet* newNormSet(0) ;
00776   RooArgSet* tmp(0) ;
00777   if (nset && !_funcNormSet) {
00778     newNormSet = nset ;
00779   } else if (!nset && _funcNormSet) {
00780     newNormSet = _funcNormSet ;
00781   } else if (nset && _funcNormSet) {
00782     tmp = new RooArgSet ;
00783     tmp->add(*nset) ;
00784     tmp->add(*_funcNormSet,kTRUE) ;
00785     newNormSet = tmp ;
00786   } 
00787   RooAbsReal* ret =  _function.arg().createIntegral(isetAll,newNormSet,cfg,rangeName) ;
00788 
00789   if (tmp) {
00790     delete tmp ;
00791   }
00792 
00793   return ret ;
00794 }
00795 
00796 
00797 
00798 
00799 //_____________________________________________________________________________
00800 Double_t RooRealIntegral::evaluate() const 
00801 {  
00802   // Perform the integration and return the result
00803 
00804   Double_t retVal(0) ;
00805   switch (_intOperMode) {    
00806     
00807   case Hybrid: 
00808     {      
00809       // Cache numeric integrals in >1d expensive object cache
00810       RooDouble* cacheVal(0) ;
00811       if ((_cacheNum && _intList.getSize()>0) || _intList.getSize()>=_cacheAllNDim) {
00812         cacheVal = (RooDouble*) expensiveObjectCache().retrieveObject(GetName(),RooDouble::Class(),parameters())  ;
00813       }
00814 
00815       if (cacheVal) {
00816         retVal = *cacheVal ;
00817 //      cout << "using cached value of integral" << GetName() << endl ;
00818       } else {
00819 
00820 
00821         // Find any function dependents that are AClean 
00822         // and switch them temporarily to ADirty
00823         setACleanADirty(kTRUE) ;
00824         
00825         // try to initialize our numerical integration engine
00826         if(!(_valid= initNumIntegrator())) {
00827           coutE(Integration) << ClassName() << "::" << GetName()
00828                              << ":evaluate: cannot initialize numerical integrator" << endl;
00829           return 0;
00830         }
00831         
00832         // Save current integral dependent values 
00833         _saveInt = _intList ;
00834         _saveSum = _sumList ;
00835         
00836         // Evaluate sum/integral
00837         retVal = sum() ;
00838         
00839         // Restore integral dependent values
00840         _intList=_saveInt ;
00841         _sumList=_saveSum ;
00842         
00843 
00844         // Cache numeric integrals in >1d expensive object cache
00845         if ((_cacheNum && _intList.getSize()>0) || _intList.getSize()>=_cacheAllNDim) {
00846           RooDouble* val = new RooDouble(retVal) ;
00847           expensiveObjectCache().registerObject(_function.arg().GetName(),GetName(),*val,parameters())  ;
00848 //        cout << "### caching value of integral" << GetName() << " in " << &expensiveObjectCache() << endl ;
00849         }
00850         
00851         setACleanADirty(kFALSE) ;
00852       }
00853       break ;
00854     }
00855   case Analytic:
00856     {
00857       retVal =  ((RooAbsReal&)_function.arg()).analyticalIntegralWN(_mode,_funcNormSet,RooNameReg::str(_rangeName)) / jacobianProduct() ;
00858       cxcoutD(Tracing) << "RooRealIntegral::evaluate_analytic(" << GetName() 
00859                        << ")func = " << _function.arg().IsA()->GetName() << "::" << _function.arg().GetName()
00860                        << " raw = " << retVal << " _funcNormSet = " << (_funcNormSet?*_funcNormSet:RooArgSet()) << endl ;
00861 
00862       
00863       break ;
00864     }
00865 
00866   case PassThrough:
00867     {
00868       //setDirtyInhibit(kTRUE) ;
00869       retVal= _function.arg().getVal(_funcNormSet) ;      
00870       //setDirtyInhibit(kFALSE) ;
00871       break ;
00872     }
00873   }
00874   
00875 
00876   // Multiply answer with integration ranges of factorized variables
00877   if (_facList.getSize()>0) {
00878     RooAbsArg *arg ;
00879     _facListIter->Reset() ;
00880     while((arg=(RooAbsArg*)_facListIter->Next())) {
00881       // Multiply by fit range for 'real' dependents
00882       if (arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
00883         RooAbsRealLValue* argLV = (RooAbsRealLValue*)arg ;
00884         retVal *= (argLV->getMax() - argLV->getMin()) ;
00885       }
00886       // Multiply by number of states for category dependents
00887       if (arg->IsA()->InheritsFrom(RooAbsCategoryLValue::Class())) {
00888         RooAbsCategoryLValue* argLV = (RooAbsCategoryLValue*)arg ;
00889         retVal *= argLV->numTypes() ;
00890       }    
00891     } 
00892   }
00893 
00894 
00895   if (dologD(Tracing)) {
00896     cxcoutD(Tracing) << "RooRealIntegral::evaluate(" << GetName() << ") anaInt = " << _anaList << " numInt = " << _intList << _sumList << " mode = " ;
00897     switch(_mode) {
00898     case Hybrid: ccoutD(Tracing) << "Hybrid" ; break ;
00899     case Analytic: ccoutD(Tracing) << "Analytic" ; break ;
00900     case PassThrough: ccoutD(Tracing) << "PassThrough" ; break ;
00901     }
00902 
00903     ccxcoutD(Tracing) << "raw*fact = " << retVal << endl ;
00904   }
00905 
00906   //   cout << "RooRealIntegral::evaluate(" << GetName() << ") value = " << retVal << endl ;
00907 
00908   return retVal ;
00909 }
00910 
00911 
00912 
00913 //_____________________________________________________________________________
00914 Double_t RooRealIntegral::jacobianProduct() const 
00915 {
00916   // Return product of jacobian terms originating from analytical integration
00917 
00918   if (_jacList.getSize()==0) {
00919     return 1 ;
00920   }
00921 
00922   Double_t jacProd(1) ;
00923   _jacListIter->Reset() ;
00924   RooAbsRealLValue* arg ;
00925   while ((arg=(RooAbsRealLValue*)_jacListIter->Next())) {
00926     jacProd *= arg->jacobian() ;
00927   }
00928 
00929   // Take fabs() here: if jacobian is negative, min and max are swapped and analytical integral
00930   // will be positive, so must multiply with positive jacobian.
00931   return fabs(jacProd) ;
00932 }
00933 
00934 
00935 
00936 //_____________________________________________________________________________
00937 Double_t RooRealIntegral::sum() const
00938 {  
00939   // Perform summation of list of category dependents to be integrated
00940  
00941   if (_sumList.getSize()!=0) {
00942  
00943     // Add integrals for all permutations of categories summed over
00944     Double_t total(0) ;
00945 
00946     _sumCatIter->Reset() ;
00947     RooCatType* type ;
00948     RooSuperCategory* sumCat = (RooSuperCategory*) _sumCat.first() ;
00949     while((type=(RooCatType*)_sumCatIter->Next())) {
00950       sumCat->setIndex(type->getVal()) ;
00951       if (!_rangeName || sumCat->inRange(RooNameReg::str(_rangeName))) {
00952         total += integrate() / jacobianProduct() ;
00953       }
00954     }
00955 
00956     return total ;
00957 
00958   } else {
00959 
00960     // Simply return integral 
00961     Double_t ret = integrate() / jacobianProduct() ;
00962     return ret ;
00963   }
00964 }
00965 
00966 
00967 
00968 
00969 //_____________________________________________________________________________
00970 Double_t RooRealIntegral::integrate() const
00971 {
00972   // Perform hybrid numerical/analytical integration over all real-valued dependents
00973 
00974   if (!_numIntEngine) {
00975     // Trivial case, fully analytical integration
00976     return ((RooAbsReal&)_function.arg()).analyticalIntegralWN(_mode,_funcNormSet,RooNameReg::str(_rangeName)) ;
00977   }
00978   else {
00979     return _numIntEngine->calculate()  ;
00980   }
00981 }
00982 
00983 
00984 
00985 //_____________________________________________________________________________
00986 Bool_t RooRealIntegral::redirectServersHook(const RooAbsCollection& /*newServerList*/, 
00987                                             Bool_t /*mustReplaceAll*/, Bool_t /*nameChange*/, Bool_t /*isRecursive*/) 
00988 {
00989   // Intercept server redirects and reconfigure internal object accordingly
00990 
00991   _restartNumIntEngine = kTRUE ;
00992 
00993   autoSelectDirtyMode() ;
00994 
00995   // Update contents value caches for _intList and _sumList
00996   _saveInt.removeAll() ;
00997   _saveSum.removeAll() ;
00998   _intList.snapshot(_saveInt) ;
00999   _sumList.snapshot(_saveSum) ;
01000 
01001   // Delete parameters cache if we have one
01002   if (_params) {
01003     delete _params ;
01004     _params = 0 ;
01005   }
01006 
01007   return kFALSE ;
01008 }
01009 
01010 
01011 
01012 //_____________________________________________________________________________
01013 const RooArgSet& RooRealIntegral::parameters() const
01014 {
01015   if (!_params) {
01016     _params = new RooArgSet("params") ;
01017     
01018     TIterator* siter = serverIterator() ;
01019     RooArgSet params ;
01020     RooAbsArg* server ;
01021     while((server = (RooAbsArg*)siter->Next())) {
01022       if (server->isValueServer(*this)) _params->add(*server) ;
01023     }
01024     delete siter ;
01025   }
01026 
01027   return *_params ;
01028 }
01029 
01030 
01031 
01032 //_____________________________________________________________________________
01033 void RooRealIntegral::operModeHook()
01034 {
01035   // Dummy
01036   if (_operMode==ADirty) {    
01037 //     cout << "RooRealIntegral::operModeHook(" << GetName() << " warning: mode set to ADirty" << endl ;
01038 //     if (TString(GetName()).Contains("FULL")) {
01039 //       cout << "blah" << endl ;
01040 //     }
01041   }
01042 }
01043 
01044 
01045 
01046 //_____________________________________________________________________________
01047 Bool_t RooRealIntegral::isValidReal(Double_t /*value*/, Bool_t /*printError*/) const 
01048 {
01049   // Check if current value is valid
01050   return kTRUE ;
01051 }
01052 
01053 
01054 
01055 //_____________________________________________________________________________
01056 void RooRealIntegral::printMetaArgs(ostream& os) const
01057 {
01058   // Customized printing of arguments of a RooRealIntegral to more intuitively reflect the contents of the
01059   // integration operation
01060 
01061 
01062   if (intVars().getSize()!=0) {
01063     os << "Int " ;
01064   }
01065   os << _function.arg().GetName() ;
01066   if (_funcNormSet) {
01067     os << "_Norm" ;
01068     os << *_funcNormSet ;
01069     os << " " ;
01070   }
01071  
01072   // List internally integrated observables and factorizing observables as analytically integrated
01073   RooArgSet tmp(_anaList) ;
01074   tmp.add(_facList) ;
01075   if (tmp.getSize()>0) {
01076     os << "d[Ana]" ;
01077     os << tmp ;
01078     os << " " ;
01079   }
01080   
01081   // List numerically integrated and summed observables as numerically integrated
01082   RooArgSet tmp2(_intList) ;
01083   tmp2.add(_sumList) ;
01084   if (tmp2.getSize()>0) {
01085     os << " d[Num]" ;
01086     os << tmp2 ;
01087     os << " " ;
01088   }
01089 }
01090 
01091 
01092 
01093 //_____________________________________________________________________________
01094 void RooRealIntegral::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
01095 {
01096   // Print the state of this object to the specified output stream.
01097 
01098   RooAbsReal::printMultiline(os,contents,verbose,indent) ;
01099   os << indent << "--- RooRealIntegral ---" << endl; 
01100   os << indent << "  Integrates ";
01101   _function.arg().printStream(os,kName|kArgs,kSingleLine,indent);
01102   TString deeper(indent);
01103   deeper.Append("  ");
01104   os << indent << "  operating mode is " 
01105      << (_intOperMode==Hybrid?"Hybrid":(_intOperMode==Analytic?"Analytic":"PassThrough")) << endl ;
01106   os << indent << "  Summed discrete args are " << _sumList << endl ;
01107   os << indent << "  Numerically integrated args are " << _intList << endl;
01108   os << indent << "  Analytically integrated args using mode " << _mode << " are " << _anaList << endl ;
01109   os << indent << "  Arguments included in Jacobian are " << _jacList << endl ;
01110   os << indent << "  Factorized arguments are " << _facList << endl ;
01111   os << indent << "  Function normalization set " ;
01112   if (_funcNormSet) 
01113     _funcNormSet->Print("1") ; 
01114   else
01115     os << "<none>" ;
01116   
01117   os << endl ;
01118 } 
01119 
01120 
01121 

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