RooAbsOptTestStatistic.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooAbsOptTestStatistic.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 // RooAbsOptTestStatistic is the abstract base class for test
00021 // statistics objects that evaluate a function or PDF at each point of a given
00022 // dataset.  This class provides generic optimizations, such as
00023 // caching and precalculation of constant terms that can be made for
00024 // all such quantities
00025 //
00026 // Implementations should define evaluatePartition(), which calculates the
00027 // value of a (sub)range of the dataset and optionally combinedValue(),
00028 // which combines the values calculated for each partition. If combinedValue()
00029 // is not overloaded, the default implementation will add the partition results
00030 // to obtain the combined result
00031 //
00032 // Support for calculation in partitions is needed to allow multi-core
00033 // parallelized calculation of test statistics
00034 // END_HTML
00035 //
00036 //
00037 
00038 #include "RooFit.h"
00039 
00040 #include "Riostream.h"
00041 #include <string.h>
00042 
00043 
00044 #include "RooAbsOptTestStatistic.h"
00045 #include "RooMsgService.h"
00046 #include "RooAbsPdf.h"
00047 #include "RooAbsData.h"
00048 #include "RooDataHist.h"
00049 #include "RooArgSet.h"
00050 #include "RooRealVar.h"
00051 #include "RooErrorHandler.h"
00052 #include "RooGlobalFunc.h"
00053 #include "RooBinning.h"
00054 #include "RooAbsDataStore.h"
00055 #include "RooCategory.h"
00056 #include "RooDataSet.h"
00057 
00058 ClassImp(RooAbsOptTestStatistic)
00059 ;
00060 
00061 
00062 //_____________________________________________________________________________
00063 RooAbsOptTestStatistic:: RooAbsOptTestStatistic()
00064 {
00065   // Default Constructor
00066 
00067   // Initialize all non-persisted data members
00068   _normSet = 0 ;
00069   _funcCloneSet = 0 ;
00070   _dataClone = 0 ;
00071   _funcClone = 0 ;
00072   _projDeps = 0 ;
00073   _ownData = kTRUE ;
00074 }
00075 
00076 
00077 
00078 //_____________________________________________________________________________
00079 RooAbsOptTestStatistic::RooAbsOptTestStatistic(const char *name, const char *title, RooAbsReal& real, RooAbsData& indata,
00080                                                const RooArgSet& projDeps, const char* rangeName, const char* addCoefRangeName,
00081                                                Int_t nCPU, Bool_t interleave, Bool_t verbose, Bool_t splitCutRange, Bool_t cloneInputData) : 
00082   RooAbsTestStatistic(name,title,real,indata,projDeps,rangeName, addCoefRangeName, nCPU, interleave, verbose, splitCutRange),
00083   _projDeps(0),
00084   _sealed(kFALSE)
00085 {
00086   // Constructor taking function (real), a dataset (data), a set of projected observables (projSet). If 
00087   // rangeName is not null, only events in the dataset inside the range will be used in the test
00088   // statistic calculation. If addCoefRangeName is not null, all RooAddPdf component of 'real' will be
00089   // instructed to fix their fraction definitions to the given named range. If nCPU is greater than
00090   // 1 the test statistic calculation will be paralellized over multiple processes. By default the data
00091   // is split with 'bulk' partitioning (each process calculates a contigious block of fraction 1/nCPU
00092   // of the data). For binned data this approach may be suboptimal as the number of bins with >0 entries
00093   // in each processing block many vary greatly thereby distributing the workload rather unevenly.
00094   // If interleave is set to true, the interleave partitioning strategy is used where each partition
00095   // i takes all bins for which (ibin % ncpu == i) which is more likely to result in an even workload.
00096   // If splitCutRange is true, a different rangeName constructed as rangeName_{catName} will be used
00097   // as range definition for each index state of a RooSimultaneous
00098 
00099 //   cout << "RooAbsOptTestStatistic::ctor(" << GetName() << "," << this << endl ;
00100   //FK: Desperate times, desperate measures. How can RooNLLVar call this ctor with dataClone=kFALSE?
00101 //   cout<<"Setting clonedata to 1"<<endl;
00102   cloneInputData=1;
00103   // Don't do a thing in master mode
00104   if (operMode()!=Slave) {
00105 //     cout << "RooAbsOptTestStatistic::ctor not slave mode, do nothing" << endl ;
00106     _normSet = 0 ;
00107     return ;
00108   }
00109 
00110   RooArgSet obs(*indata.get()) ;
00111   obs.remove(projDeps,kTRUE,kTRUE) ;
00112 
00113   // + ALEX
00114   //   // Check that the FUNC is valid for use with this dataset
00115   //   // Check if there are any unprotected multiple occurrences of dependents
00116   //   if (real.recursiveCheckObservables(&obs)) {
00117   //     coutE(InputArguments) << "RooAbsOptTestStatistic: ERROR in FUNC dependents, abort" << endl ;
00118   //     RooErrorHandler::softAbort() ;
00119   //     return ;
00120   //   }
00121   // - ALEX  
00122 
00123 
00124   // Get list of actual observables of test statistic function
00125   RooArgSet* realDepSet = real.getObservables(&indata) ;
00126 
00127   // Expand list of observables with any observables used in parameterized ranges
00128   TIterator* iter9 = realDepSet->createIterator() ;
00129   RooAbsArg* realDep ;
00130   while((realDep=(RooAbsArg*)iter9->Next())) {
00131     RooAbsRealLValue* realDepRLV = dynamic_cast<RooAbsRealLValue*>(realDep) ;
00132     if (realDepRLV && realDepRLV->isDerived()) {
00133 
00134       RooArgSet tmp ;
00135       realDepRLV->leafNodeServerList(&tmp, 0, kTRUE) ;
00136       realDepSet->add(tmp,kTRUE) ;
00137     }
00138   }
00139   delete iter9 ;
00140 
00141 
00142   // Check if the fit ranges of the dependents in the data and in the FUNC are consistent
00143   const RooArgSet* dataDepSet = indata.get() ;
00144   TIterator* iter = realDepSet->createIterator() ;
00145   RooAbsArg* arg ;
00146   while((arg=(RooAbsArg*)iter->Next())) {
00147     RooRealVar* realReal = dynamic_cast<RooRealVar*>(arg) ;
00148     if (!realReal) continue ;
00149 
00150     
00151     RooRealVar* datReal = dynamic_cast<RooRealVar*>(dataDepSet->find(realReal->GetName())) ;
00152     if (!datReal) continue ;
00153     
00154     if (!realReal->getBinning().lowBoundFunc() && realReal->getMin()<(datReal->getMin()-1e-6)) {
00155       coutE(InputArguments) << "RooAbsOptTestStatistic: ERROR minimum of FUNC observable " << arg->GetName() 
00156                             << "(" << realReal->getMin() << ") is smaller than that of " 
00157                             << arg->GetName() << " in the dataset (" << datReal->getMin() << ")" << endl ;
00158       RooErrorHandler::softAbort() ;
00159       return ;
00160     }
00161     
00162     if (!realReal->getBinning().highBoundFunc() && realReal->getMax()>(datReal->getMax()+1e-6)) {
00163       coutE(InputArguments) << "RooAbsOptTestStatistic: ERROR maximum of FUNC observable " << arg->GetName() 
00164                             << " is larger than that of " << arg->GetName() << " in the dataset" << endl ;
00165       RooErrorHandler::softAbort() ;
00166       return ;
00167     }
00168     
00169   }
00170   
00171   // Copy data and strip entries lost by adjusted fit range, _dataClone ranges will be copied from realDepSet ranges
00172   if (rangeName && strlen(rangeName)) {
00173     if (!cloneInputData) {
00174       coutW(InputArguments) << "RooAbsOptTestStatistic::ctor(" << GetName() 
00175                             << ") WARNING: Must clone input data when a range specification is given, ignoring request to use original input dataset" << endl ; 
00176     }    
00177     _dataClone = ((RooAbsData&)indata).reduce(RooFit::SelectVars(*realDepSet),RooFit::CutRange(rangeName)) ;  
00178 //     cout << "RooAbsOptTestStatistic: reducing dataset to fit in range named " << rangeName << " resulting dataset has " << _dataClone->sumEntries() << " events" << endl ;
00179     _ownData = kTRUE ;
00180   } else {
00181     if (cloneInputData) {
00182       _dataClone = (RooAbsData*) indata.Clone() ;
00183       //reduce(RooFit::SelectVars(*indata.get())) ; //  ((RooAbsData&)data).reduce(RooFit::SelectVars(*realDepSet)) ;  
00184 
00185       _ownData = kTRUE ;
00186     } else {
00187       // coverity[DEADCODE]
00188       _dataClone = &indata ;
00189       _ownData = kFALSE ;
00190     }
00191   }
00192 
00193 
00194 
00195   // Copy any non-shared parameterized range definitions from pdf observables to dataset observables
00196   iter9 = realDepSet->createIterator() ;
00197   while((realDep=(RooAbsRealLValue*)iter9->Next())) {
00198     RooAbsRealLValue* realDepRLV = dynamic_cast<RooAbsRealLValue*>(realDep) ;
00199     if (realDepRLV && !realDepRLV->getBinning().isShareable()) {
00200 
00201       RooRealVar* datReal = dynamic_cast<RooRealVar*>(_dataClone->get()->find(realDepRLV->GetName())) ;
00202       if (datReal) {
00203         datReal->setBinning(realDepRLV->getBinning()) ;
00204         datReal->attachDataSet(*_dataClone) ;
00205       }
00206     }
00207   }
00208   delete iter9 ;
00209 
00210   if (rangeName && strlen(rangeName)) {
00211     
00212     cxcoutI(Fitting) << "RooAbsOptTestStatistic::ctor(" << GetName() << ") constructing test statistic for sub-range named " << rangeName << endl ;
00213 
00214     // Adjust FUNC normalization ranges to requested fitRange, store original ranges for RooAddPdf coefficient interpretation
00215     TIterator* iter2 = _dataClone->get()->createIterator() ;
00216     while((arg=(RooAbsArg*)iter2->Next())) {
00217       RooRealVar* realReal = dynamic_cast<RooRealVar*>(arg) ;
00218       if (realReal) {
00219         if (!(addCoefRangeName && strlen(addCoefRangeName))) {
00220           realReal->setRange(Form("NormalizationRangeFor%s",rangeName),realReal->getMin(),realReal->getMax()) ;
00221         }
00222 
00223         realReal->setRange(realReal->getMin(rangeName),realReal->getMax(rangeName)) ;   
00224       }
00225     }
00226 
00227     // If dataset is binned, activate caching of bins that are invalid because the're outside the
00228     // updated range definition (WVE need to add virtual interface here)
00229     RooDataHist* tmp = dynamic_cast<RooDataHist*>(_dataClone) ;
00230     if (tmp) {
00231       tmp->cacheValidEntries() ;
00232     }
00233 
00234 
00235     // Mark fitted range in original FUNC dependents for future use
00236     if (!_splitRange) {
00237       iter->Reset() ;
00238       while((arg=(RooAbsArg*)iter->Next())) {      
00239         RooRealVar* realReal = dynamic_cast<RooRealVar*>(arg) ;
00240         if (realReal) {
00241           realReal->setStringAttribute("fitrange",Form("fit_%s",GetName())) ;
00242           realReal->setRange(Form("fit_%s",GetName()),realReal->getMin(rangeName),realReal->getMax(rangeName)) ;
00243 
00244           // Keep track of list of fit ranges in string attribute fit range of original p.d.f.
00245           const char* origAttrib = real.getStringAttribute("fitrange") ;          
00246           if (origAttrib) {
00247             real.setStringAttribute("fitrange",Form("%s,fit_%s",origAttrib,GetName())) ;
00248           } else {
00249             real.setStringAttribute("fitrange",Form("fit_%s",GetName())) ;
00250           }
00251         }
00252       }
00253     }
00254     delete iter2 ;
00255   }
00256   delete iter ;
00257 
00258   setEventCount(_dataClone->numEntries()) ;
00259  
00260   // Clone all FUNC compents by copying all branch nodes
00261   RooArgSet tmp("RealBranchNodeList") ;
00262   real.branchNodeServerList(&tmp) ;
00263   _funcCloneSet = (RooArgSet*) tmp.snapshot(kFALSE) ;
00264   
00265   // Find the top level FUNC in the snapshot list
00266   _funcClone = (RooAbsReal*) _funcCloneSet->find(real.GetName()) ;
00267 
00268 
00269   // Fix RooAddPdf coefficients to original normalization range
00270   if (rangeName && strlen(rangeName)) {
00271 
00272     // WVE Remove projected dependents from normalization
00273     _funcClone->fixAddCoefNormalization(*_dataClone->get(),kFALSE) ;
00274     
00275     if (addCoefRangeName && strlen(addCoefRangeName)) {
00276       cxcoutI(Fitting) << "RooAbsOptTestStatistic::ctor(" << GetName() 
00277                        << ") fixing interpretation of coefficients of any RooAddPdf component to range " << addCoefRangeName << endl ;
00278       _funcClone->fixAddCoefRange(addCoefRangeName,kFALSE) ;
00279     } else {
00280         cxcoutI(Fitting) << "RooAbsOptTestStatistic::ctor(" << GetName() 
00281                          << ") fixing interpretation of coefficients of any RooAddPdf to full domain of observables " << endl ;
00282         _funcClone->fixAddCoefRange(Form("NormalizationRangeFor%s",rangeName),kFALSE) ;
00283     }
00284   }
00285 
00286 
00287   // Attach FUNC to data set
00288   _funcClone->attachDataSet(*_dataClone) ;
00289 
00290 
00291   // Store normalization set  
00292   _normSet = (RooArgSet*) indata.get()->snapshot(kFALSE) ;
00293 
00294   // Remove projected dependents from normalization set
00295   if (projDeps.getSize()>0) {
00296 
00297     _projDeps = (RooArgSet*) projDeps.snapshot(kFALSE) ;
00298     
00299     RooArgSet* tobedel = (RooArgSet*) _normSet->selectCommon(*_projDeps) ;
00300     _normSet->remove(*_projDeps,kTRUE,kTRUE) ;
00301 
00302     // Delete owned projected dependent copy in _normSet
00303     TIterator* ii = tobedel->createIterator() ;
00304     RooAbsArg* aa ;
00305     while((aa=(RooAbsArg*)ii->Next())) {
00306       delete aa ;
00307     }
00308     delete ii ;
00309     delete tobedel ;
00310 
00311     // Mark all projected dependents as such
00312     RooArgSet *projDataDeps = (RooArgSet*) _dataClone->get()->selectCommon(*_projDeps) ;
00313     projDataDeps->setAttribAll("projectedDependent") ;
00314     delete projDataDeps ;
00315   } 
00316 
00317 //   cout << "RAOTS: Now evaluating funcClone with _normSet = " << _normSet << " = " << *_normSet << endl ;
00318   _funcClone->getVal(_normSet) ;
00319 
00320   // Add parameters as servers
00321   RooArgSet* params = _funcClone->getParameters(_dataClone) ;
00322   _paramSet.add(*params) ;
00323   delete params ;
00324 
00325   // Mark all projected dependents as such
00326   if (_projDeps) {
00327     RooArgSet *projDataDeps = (RooArgSet*) _dataClone->get()->selectCommon(*_projDeps) ;
00328     projDataDeps->setAttribAll("projectedDependent") ;
00329     delete projDataDeps ;
00330   }
00331 
00332   coutI(Optimization) << "RooAbsOptTestStatistic::ctor(" << GetName() << ") optimizing internal clone of p.d.f for likelihood evaluation." 
00333                         << "Lazy evaluation and associated change tracking will disabled for all nodes that depend on observables" << endl ;
00334 
00335 
00336 
00337   delete realDepSet ;  
00338 
00339   // Redirect pointers of base class to clone 
00340   _func = _funcClone ;
00341   _data = _dataClone ;
00342 
00343 //   cout << "RooAbsOptTestStatistic::ctor call getVal()" << endl ;
00344   _funcClone->getVal(_normSet) ;
00345 
00346   
00347 
00348 //   cout << "RooAbsOptTestStatistic::ctor start caching" << endl ;
00349   optimizeCaching() ;
00350   
00351 }
00352 
00353 
00354 //_____________________________________________________________________________
00355 RooAbsOptTestStatistic::RooAbsOptTestStatistic(const RooAbsOptTestStatistic& other, const char* name) : 
00356   RooAbsTestStatistic(other,name), _sealed(other._sealed), _sealNotice(other._sealNotice)
00357 {
00358   // Copy constructor
00359 //   cout << "RooAbsOptTestStatistic::cctor(" << GetName() << "," << this << endl ;
00360 
00361   // Don't do a thing in master mode
00362   if (operMode()!=Slave) {
00363 //     cout << "RooAbsOptTestStatistic::cctor not slave mode, do nothing" << endl ;
00364     _projDeps = 0 ;
00365     _normSet = other._normSet ? ((RooArgSet*) other._normSet->snapshot()) : 0 ;   
00366     return ;
00367   }
00368 
00369   _funcCloneSet = (RooArgSet*) other._funcCloneSet->snapshot(kFALSE) ;
00370   _funcClone = (RooAbsReal*) _funcCloneSet->find(other._funcClone->GetName()) ;
00371 
00372   // Copy the operMode attribute of all branch nodes
00373   TIterator* iter = _funcCloneSet->createIterator() ;
00374   RooAbsArg* branch ;
00375   while((branch=(RooAbsArg*)iter->Next())) {
00376     branch->setOperMode(other._funcCloneSet->find(branch->GetName())->operMode()) ;
00377   }
00378   delete iter ;
00379 
00380   // WVE Must use clone with cache redirection here
00381   if (other._ownData || other._dataClone->hasFilledCache()) {    
00382     _dataClone = (RooAbsData*) other._dataClone->cacheClone(this,_funcCloneSet) ;
00383     _ownData = kTRUE ;
00384   } else {
00385     _dataClone = other._dataClone ;
00386     _ownData = kFALSE ;
00387     
00388     // Revert any AClean nodes imported from original to ADirty as not optimization is applicable to test statistics with borrowed data
00389     Bool_t wasOpt(kFALSE) ;
00390     TIterator* biter = _funcCloneSet->createIterator() ;
00391     RooAbsArg *branch2 ;
00392     while((branch2=(RooAbsArg*)biter->Next())){
00393       if (branch2->operMode()==RooAbsArg::AClean) {
00394 //      cout << "RooAbsOptTestStatistic::cctor(" << GetName() << " setting branch " << branch2->GetName() << " to ADirty" << endl ;
00395         branch2->setOperMode(RooAbsArg::ADirty) ;
00396         wasOpt=kTRUE ;
00397       }
00398     }
00399     delete biter ;  
00400 
00401     if (wasOpt) {
00402       coutW(Optimization) << "RooAbsOptTestStatistic::cctor(" << GetName() << ") WARNING clone of optimized test statistic with unowned data will not be optimized, "
00403                           << "to retain optimization behavior in cloning, construct test statistic with CloneData(kTRUE)" << endl ;
00404     }
00405   }
00406 
00407   // Attach function clone to dataset
00408   _funcClone->attachDataSet(*_dataClone) ;
00409 
00410   // Explicit attach function clone to current parameter instances
00411   _funcClone->recursiveRedirectServers(_paramSet) ;
00412 
00413   _normSet = (RooArgSet*) other._normSet->snapshot() ;
00414   
00415   if (other._projDeps) {
00416     _projDeps = (RooArgSet*) other._projDeps->snapshot() ;
00417   } else {
00418     _projDeps = 0 ;
00419   }
00420 
00421   _func = _funcClone ;
00422   _data = _dataClone ;
00423 
00424 //   cout << "RooAbsOptTestStatistic::cctor call getVal()" << endl ;
00425   _funcClone->getVal(_normSet) ;
00426 //   cout << "RooAbsOptTestStatistic::cctor start caching" << endl ;
00427   optimizeCaching() ;
00428 }
00429 
00430 
00431 
00432 //_____________________________________________________________________________
00433 RooAbsOptTestStatistic::~RooAbsOptTestStatistic()
00434 {
00435   // Destructor
00436 
00437   if (operMode()==Slave) {
00438     delete _funcCloneSet ;
00439     if (_ownData) {
00440       delete _dataClone ;
00441     } else {
00442       _dataClone->resetCache() ;
00443     }
00444     delete _projDeps ;
00445   } 
00446   delete _normSet ;
00447 }
00448 
00449 
00450 
00451 //_____________________________________________________________________________
00452 Double_t RooAbsOptTestStatistic::combinedValue(RooAbsReal** array, Int_t n) const
00453 {
00454   // Method to combined test statistic results calculated into partitions into
00455   // the global result. This default implementation adds the partition return
00456   // values
00457   
00458   // Default implementation returns sum of components
00459   Double_t sum(0) ;
00460   Int_t i ;
00461   for (i=0 ; i<n ; i++) {
00462     Double_t tmp = array[i]->getVal() ;
00463     if (tmp==0) return 0 ;
00464     sum += tmp ;
00465   }
00466   return sum ;
00467 }
00468 
00469 
00470 
00471 //_____________________________________________________________________________
00472 Bool_t RooAbsOptTestStatistic::redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive) 
00473 {
00474   // Catch server redirect calls and forward to internal clone of function
00475 
00476   RooAbsTestStatistic::redirectServersHook(newServerList,mustReplaceAll,nameChange,isRecursive) ;
00477   if (operMode()!=Slave) return kFALSE ;  
00478   Bool_t ret = _funcClone->recursiveRedirectServers(newServerList,kFALSE,nameChange) ;
00479   return ret ;
00480 }
00481 
00482 
00483 
00484 //_____________________________________________________________________________
00485 void RooAbsOptTestStatistic::printCompactTreeHook(ostream& os, const char* indent) 
00486 {
00487   // Catch print hook function and forward to function clone
00488 
00489   RooAbsTestStatistic::printCompactTreeHook(os,indent) ;
00490   if (operMode()!=Slave) return ;
00491   TString indent2(indent) ;
00492   indent2 += "opt >>" ;
00493   _funcClone->printCompactTree(os,indent2.Data()) ;
00494   os << indent2 << " dataset clone = " << _dataClone << " first obs = " << _dataClone->get()->first() << endl ;
00495 }
00496 
00497 
00498 
00499 //_____________________________________________________________________________
00500 void RooAbsOptTestStatistic::constOptimizeTestStatistic(ConstOpCode opcode) 
00501 {
00502   // Driver function to propagate constant term optimizations in test statistic.
00503   // If code Activate is sent, constant term optimization will be executed.
00504   // If code Deacivate is sent, any existing constant term optimizations will
00505   // be abanoned. If codes ConfigChange or ValueChange are sent, any existing
00506   // constant term optimizations will be redone.
00507 
00508   RooAbsTestStatistic::constOptimizeTestStatistic(opcode);
00509   if (operMode()!=Slave) return ;
00510 
00511   if (_dataClone->hasFilledCache() && _dataClone->store()->cacheOwner()!=this) {
00512     if (opcode==Activate) {
00513       cxcoutW(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() 
00514                             << ") dataset cache is owned by another object, no constant term optimization can be applied" << endl ;
00515     }
00516     return ;
00517   }
00518 
00519 
00520   if (!allowFunctionCache()) {
00521     if (opcode==Activate) {
00522       cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() 
00523                             << ") function caching prohibited by test statistic, no constant term optimization is applied" << endl ;
00524     }
00525     return ;
00526   }
00527 
00528   switch(opcode) {
00529   case Activate:     
00530     cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() 
00531                           << ") optimizing evaluation of test statistic by finding all nodes in p.d.f that depend exclusively"
00532                           << " on observables and constant parameters and precalculating their values" << endl ;
00533     optimizeConstantTerms(kTRUE) ;
00534     break ;
00535 
00536   case DeActivate:  
00537     cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() 
00538                           << ") deactivating optimization of constant terms in test statistic" << endl ;
00539     optimizeConstantTerms(kFALSE) ;
00540     break ;
00541 
00542   case ConfigChange: 
00543     cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() 
00544                           << ") one ore more parameter were changed from constant to floating or vice versa, "
00545                           << "re-evaluating constant term optimization" << endl ;
00546     optimizeConstantTerms(kFALSE) ;
00547     optimizeConstantTerms(kTRUE) ;
00548     break ;
00549 
00550   case ValueChange: 
00551     cxcoutI(Optimization) << "RooAbsOptTestStatistic::constOptimize(" << GetName() 
00552                           << ") the value of one ore more constant parameter were changed re-evaluating constant term optimization" << endl ;
00553     optimizeConstantTerms(kFALSE) ;
00554     optimizeConstantTerms(kTRUE) ;
00555     break ;
00556   }
00557 }
00558 
00559 
00560 
00561 //_____________________________________________________________________________
00562 void RooAbsOptTestStatistic::optimizeCaching() 
00563 {
00564   // This method changes the value caching logic for all nodes that depends on any of the observables
00565   // as defined by the given dataset. When evaluating a test statistic constructed from the RooAbsReal
00566   // with a dataset the observables are guaranteed to change with every call, thus there is no point
00567   // in tracking these changes which result in a net overhead. Thus for observable-dependent nodes, 
00568   // the evaluation mechanism is changed from being dependent on a 'valueDirty' flag to guaranteed evaluation. 
00569   // On the dataset side, the observables objects are modified to no longer send valueDirty messages
00570   // to their client 
00571 
00572 //   cout << "RooAbsOptTestStatistic::optimizeCaching(" << GetName() << "," << this << ")" << endl ;
00573 
00574   // Trigger create of all object caches now in nodes that have deferred object creation
00575   // so that cache contents can be processed immediately
00576   _funcClone->getVal(_normSet) ;
00577 
00578   // Set value caching mode for all nodes that depend on any of the observables to ADirty
00579   _funcClone->optimizeCacheMode(*_dataClone->get()) ;
00580 
00581   // Disable propagation of dirty state flags for observables
00582   _dataClone->setDirtyProp(kFALSE) ;  
00583 
00584   // Disable reading of observables that are not used
00585   _dataClone->optimizeReadingWithCaching(*_funcClone, RooArgSet(),requiredExtraObservables()) ;
00586 }
00587 
00588 
00589 
00590 //_____________________________________________________________________________
00591 void RooAbsOptTestStatistic::optimizeConstantTerms(Bool_t activate)
00592 {
00593   // Driver function to activate global constant term optimization.
00594   // If activated constant terms are found and cached with the dataset
00595   // The operation mode of cached nodes is set to AClean meaning that
00596   // their getVal() call will never result in an evaluate call.
00597   // Finally the branches in the dataset that correspond to observables
00598   // that are exclusively used in constant terms are disabled as
00599   // they serve no more purpose
00600 
00601   if(activate) {
00602     // Trigger create of all object caches now in nodes that have deferred object creation
00603     // so that cache contents can be processed immediately
00604     _funcClone->getVal(_normSet) ;
00605     
00606     // Find all nodes that depend exclusively on constant parameters
00607     RooArgSet cacheableNodes ;
00608 
00609     _funcClone->findConstantNodes(*_dataClone->get(),cacheableNodes) ;
00610 
00611     // Cache constant nodes with dataset 
00612     _dataClone->cacheArgs(this,cacheableNodes,_normSet) ;  
00613     
00614     // Put all cached nodes in AClean value caching mode so that their evaluate() is never called
00615     TIterator* cIter = cacheableNodes.createIterator() ;
00616     RooAbsArg *cacheArg ;
00617     while((cacheArg=(RooAbsArg*)cIter->Next())){
00618       cacheArg->setOperMode(RooAbsArg::AClean) ;
00619     }
00620     delete cIter ;  
00621     
00622     // Disable reading of observables that are no longer used
00623     _dataClone->optimizeReadingWithCaching(*_funcClone, cacheableNodes,requiredExtraObservables()) ;
00624 
00625   } else {
00626     
00627     // Delete the cache
00628     _dataClone->resetCache() ;
00629     
00630     // Reactivate all tree branches
00631     _dataClone->setArgStatus(*_dataClone->get(),kTRUE) ;
00632     
00633     // Reset all nodes to ADirty   
00634     optimizeCaching() ;
00635 
00636     // Disable propagation of dirty state flags for observables
00637     _dataClone->setDirtyProp(kFALSE) ;  
00638     
00639   }
00640 }
00641 
00642 
00643 
00644 //_____________________________________________________________________________
00645 Bool_t RooAbsOptTestStatistic::setData(RooAbsData& indata, Bool_t cloneData) 
00646 { 
00647   // Change dataset that is used to given one. If cloneData is kTRUE, a clone of
00648   // in the input dataset is made.  If the test statistic was constructed with
00649   // a range specification on the data, the cloneData argument is ignore and
00650   // the data is always cloned.
00651 
00652   RooAbsData* origData = _dataClone ;
00653   Bool_t deleteOrigData = _ownData ;
00654 
00655   if (!cloneData && _rangeName.size()>0) {
00656     coutW(InputArguments) << "RooAbsOptTestStatistic::setData(" << GetName() << ") WARNING: test statistic was constructed with range selection on data, "
00657                          << "ignoring request to _not_ clone the input dataset" << endl ; 
00658     cloneData = kTRUE ;
00659   }
00660 
00661   if (cloneData) {
00662     if (_rangeName.size()==0) {
00663       _dataClone = (RooAbsData*) indata.reduce(*indata.get()) ;
00664     } else {
00665       _dataClone = ((RooAbsData&)indata).reduce(RooFit::SelectVars(*indata.get()),RooFit::CutRange(_rangeName.c_str())) ;  
00666     }
00667     _ownData = kTRUE ;
00668   } else {
00669     _dataClone = &indata ;
00670     _ownData = kFALSE ;
00671   }    
00672   // Attach function clone to dataset
00673   _funcClone->attachDataSet(*_dataClone) ;
00674 
00675   _data = _dataClone ;
00676 
00677   if (deleteOrigData) {
00678     delete origData ;
00679   } else {
00680     origData->resetCache() ;
00681   }
00682 
00683   setValueDirty() ;
00684   return kTRUE ;
00685 }
00686 
00687 
00688 
00689 
00690 //_____________________________________________________________________________
00691 RooAbsData& RooAbsOptTestStatistic::data() 
00692 { 
00693   if (_sealed) {
00694     Bool_t notice = (sealNotice() && strlen(sealNotice())) ;
00695     coutW(ObjectHandling) << "RooAbsOptTestStatistic::data(" << GetName() 
00696                           << ") WARNING: object sealed by creator - access to data is not permitted: " 
00697                           << (notice?sealNotice():"<no user notice>") << endl ;
00698     static RooDataSet dummy ("dummy","dummy",RooArgSet()) ;
00699     return dummy ;
00700   }
00701   return *_dataClone ; 
00702 }
00703 
00704 
00705 //_____________________________________________________________________________
00706 const RooAbsData& RooAbsOptTestStatistic::data() const 
00707 { 
00708   if (_sealed) {
00709     Bool_t notice = (sealNotice() && strlen(sealNotice())) ;
00710     coutW(ObjectHandling) << "RooAbsOptTestStatistic::data(" << GetName() 
00711                           << ") WARNING: object sealed by creator - access to data is not permitted: " 
00712                           << (notice?sealNotice():"<no user notice>") << endl ;
00713     static RooDataSet dummy ("dummy","dummy",RooArgSet()) ;
00714     return dummy ;
00715   }
00716   return *_dataClone ; 
00717 }

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