00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
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
00066
00067
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
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102 cloneInputData=1;
00103
00104 if (operMode()!=Slave) {
00105
00106 _normSet = 0 ;
00107 return ;
00108 }
00109
00110 RooArgSet obs(*indata.get()) ;
00111 obs.remove(projDeps,kTRUE,kTRUE) ;
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125 RooArgSet* realDepSet = real.getObservables(&indata) ;
00126
00127
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
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
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
00179 _ownData = kTRUE ;
00180 } else {
00181 if (cloneInputData) {
00182 _dataClone = (RooAbsData*) indata.Clone() ;
00183
00184
00185 _ownData = kTRUE ;
00186 } else {
00187
00188 _dataClone = &indata ;
00189 _ownData = kFALSE ;
00190 }
00191 }
00192
00193
00194
00195
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
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
00228
00229 RooDataHist* tmp = dynamic_cast<RooDataHist*>(_dataClone) ;
00230 if (tmp) {
00231 tmp->cacheValidEntries() ;
00232 }
00233
00234
00235
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
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
00261 RooArgSet tmp("RealBranchNodeList") ;
00262 real.branchNodeServerList(&tmp) ;
00263 _funcCloneSet = (RooArgSet*) tmp.snapshot(kFALSE) ;
00264
00265
00266 _funcClone = (RooAbsReal*) _funcCloneSet->find(real.GetName()) ;
00267
00268
00269
00270 if (rangeName && strlen(rangeName)) {
00271
00272
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
00288 _funcClone->attachDataSet(*_dataClone) ;
00289
00290
00291
00292 _normSet = (RooArgSet*) indata.get()->snapshot(kFALSE) ;
00293
00294
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
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
00312 RooArgSet *projDataDeps = (RooArgSet*) _dataClone->get()->selectCommon(*_projDeps) ;
00313 projDataDeps->setAttribAll("projectedDependent") ;
00314 delete projDataDeps ;
00315 }
00316
00317
00318 _funcClone->getVal(_normSet) ;
00319
00320
00321 RooArgSet* params = _funcClone->getParameters(_dataClone) ;
00322 _paramSet.add(*params) ;
00323 delete params ;
00324
00325
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
00340 _func = _funcClone ;
00341 _data = _dataClone ;
00342
00343
00344 _funcClone->getVal(_normSet) ;
00345
00346
00347
00348
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
00359
00360
00361
00362 if (operMode()!=Slave) {
00363
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
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
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
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
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
00408 _funcClone->attachDataSet(*_dataClone) ;
00409
00410
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
00425 _funcClone->getVal(_normSet) ;
00426
00427 optimizeCaching() ;
00428 }
00429
00430
00431
00432
00433 RooAbsOptTestStatistic::~RooAbsOptTestStatistic()
00434 {
00435
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
00455
00456
00457
00458
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
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
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
00503
00504
00505
00506
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
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576 _funcClone->getVal(_normSet) ;
00577
00578
00579 _funcClone->optimizeCacheMode(*_dataClone->get()) ;
00580
00581
00582 _dataClone->setDirtyProp(kFALSE) ;
00583
00584
00585 _dataClone->optimizeReadingWithCaching(*_funcClone, RooArgSet(),requiredExtraObservables()) ;
00586 }
00587
00588
00589
00590
00591 void RooAbsOptTestStatistic::optimizeConstantTerms(Bool_t activate)
00592 {
00593
00594
00595
00596
00597
00598
00599
00600
00601 if(activate) {
00602
00603
00604 _funcClone->getVal(_normSet) ;
00605
00606
00607 RooArgSet cacheableNodes ;
00608
00609 _funcClone->findConstantNodes(*_dataClone->get(),cacheableNodes) ;
00610
00611
00612 _dataClone->cacheArgs(this,cacheableNodes,_normSet) ;
00613
00614
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
00623 _dataClone->optimizeReadingWithCaching(*_funcClone, cacheableNodes,requiredExtraObservables()) ;
00624
00625 } else {
00626
00627
00628 _dataClone->resetCache() ;
00629
00630
00631 _dataClone->setArgStatus(*_dataClone->get(),kTRUE) ;
00632
00633
00634 optimizeCaching() ;
00635
00636
00637 _dataClone->setDirtyProp(kFALSE) ;
00638
00639 }
00640 }
00641
00642
00643
00644
00645 Bool_t RooAbsOptTestStatistic::setData(RooAbsData& indata, Bool_t cloneData)
00646 {
00647
00648
00649
00650
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
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 }