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
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392
00393
00394
00395
00396
00397
00398
00399
00400
00401
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418 #include "RooFit.h"
00419
00420 #include <string.h>
00421 #include <string.h>
00422
00423 #ifndef _WIN32
00424 #include <strings.h>
00425 #else
00426
00427
00428 static char *strtok_r(char *s1, const char *s2, char **lasts)
00429 {
00430 char *ret;
00431
00432 if (s1 == NULL)
00433 s1 = *lasts;
00434 while(*s1 && strchr(s2, *s1))
00435 ++s1;
00436 if(*s1 == '\0')
00437 return NULL;
00438 ret = s1;
00439 while(*s1 && !strchr(s2, *s1))
00440 ++s1;
00441 if(*s1)
00442 *s1++ = '\0';
00443 *lasts = s1;
00444 return ret;
00445 }
00446
00447 #endif
00448
00449 #include "Riostream.h"
00450 #include "RooSimPdfBuilder.h"
00451
00452 #include "RooRealVar.h"
00453 #include "RooFormulaVar.h"
00454 #include "RooAbsCategory.h"
00455 #include "RooCategory.h"
00456 #include "RooStringVar.h"
00457 #include "RooMappedCategory.h"
00458 #include "RooRealIntegral.h"
00459 #include "RooDataSet.h"
00460 #include "RooArgSet.h"
00461 #include "RooPlot.h"
00462 #include "RooAddPdf.h"
00463 #include "RooLinearVar.h"
00464 #include "RooTruthModel.h"
00465 #include "RooAddModel.h"
00466 #include "RooProdPdf.h"
00467 #include "RooCustomizer.h"
00468 #include "RooThresholdCategory.h"
00469 #include "RooMultiCategory.h"
00470 #include "RooSuperCategory.h"
00471 #include "RooSimultaneous.h"
00472 #include "RooTrace.h"
00473 #include "RooFitResult.h"
00474 #include "RooDataHist.h"
00475 #include "RooGenericPdf.h"
00476 #include "RooMsgService.h"
00477
00478 using namespace std ;
00479
00480 ClassImp(RooSimPdfBuilder)
00481 ;
00482
00483
00484
00485
00486 RooSimPdfBuilder::RooSimPdfBuilder(const RooArgSet& protoPdfSet) :
00487 _protoPdfSet(protoPdfSet)
00488 {
00489 _compSplitCatSet.setHashTableSize(1000) ;
00490 _splitNodeList.setHashTableSize(10000) ;
00491 _splitNodeListOwned.setHashTableSize(10000) ;
00492 }
00493
00494
00495
00496
00497
00498 RooArgSet* RooSimPdfBuilder::createProtoBuildConfig()
00499 {
00500
00501 RooArgSet* buildConfig = new RooArgSet ;
00502 buildConfig->addOwned(* new RooStringVar("physModels","List and mapping of physics models to include in build","",4096)) ;
00503 buildConfig->addOwned(* new RooStringVar("splitCats","List of categories used for splitting","",1024)) ;
00504
00505 TIterator* iter = _protoPdfSet.createIterator() ;
00506 RooAbsPdf* proto ;
00507 while ((proto=(RooAbsPdf*)iter->Next())) {
00508 buildConfig->addOwned(* new RooStringVar(proto->GetName(),proto->GetName(),"",4096)) ;
00509 }
00510 delete iter ;
00511
00512 return buildConfig ;
00513 }
00514
00515
00516
00517
00518 void RooSimPdfBuilder::addSpecializations(const RooArgSet& specSet)
00519 {
00520 _splitNodeList.add(specSet) ;
00521 }
00522
00523
00524
00525
00526 RooSimultaneous* RooSimPdfBuilder::buildPdf(const RooArgSet& buildConfig, const RooArgSet& dependents,
00527 const RooArgSet* auxSplitCats, Bool_t verbose)
00528 {
00529
00530 const char* spaceChars = " \t" ;
00531
00532
00533 Int_t buflen = strlen(((RooStringVar*)buildConfig.find("physModels"))->getVal())+1 ;
00534 char *buf = new char[buflen] ;
00535
00536 strlcpy(buf,((RooStringVar*)buildConfig.find("physModels"))->getVal(),buflen) ;
00537 RooAbsCategoryLValue* physCat(0) ;
00538 if (strstr(buf," : ")) {
00539 const char* physCatName = strtok(buf,spaceChars) ;
00540 physCat = dynamic_cast<RooAbsCategoryLValue*>(dependents.find(physCatName)) ;
00541 if (!physCat) {
00542 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR physics index category " << physCatName
00543 << " not found in dataset variables" << endl ;
00544 delete[] buf ;
00545 return 0 ;
00546 }
00547 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: category indexing physics model: " << physCatName << endl ;
00548 }
00549
00550
00551 char *physName ;
00552 RooArgSet physModelSet ;
00553 if (physCat) {
00554
00555 strtok(0,spaceChars) ;
00556 physName = strtok(0,spaceChars) ;
00557 } else {
00558 physName = strtok(buf,spaceChars) ;
00559 }
00560
00561 if (!physName) {
00562 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR: No models specified, nothing to do!" << endl ;
00563 delete[] buf ;
00564 return 0 ;
00565 }
00566
00567 Bool_t first(kTRUE) ;
00568 RooArgSet stateMap ;
00569 while(physName) {
00570
00571 char *stateName(0) ;
00572
00573
00574 if (strchr(physName,'=')) {
00575
00576 if (!physCat) {
00577 coutW(ObjectHandling) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification "
00578 << "<physCatState>=<pdfProtoName> association is meaningless" << endl ;
00579 }
00580 stateName = physName ;
00581 physName = strchr(stateName,'=') ;
00582 if (physName) {
00583 *(physName++) = 0 ;
00584 }
00585 } else {
00586 stateName = physName ;
00587 }
00588
00589 RooAbsPdf* physModel = (RooAbsPdf*) (physName ? _protoPdfSet.find(physName) : 0 );
00590 if (!physModel) {
00591 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested physics model "
00592 << (physName?physName:"(null)") << " is not defined" << endl ;
00593 delete[] buf ;
00594 return 0 ;
00595 }
00596
00597
00598 if (stateMap.find(stateName)) {
00599 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: multiple PDFs specified for state "
00600 << stateName << ", only first will be used" << endl ;
00601 continue ;
00602 }
00603
00604
00605 physModelSet.add(*physModel,kTRUE) ;
00606
00607
00608 stateMap.addOwned(* new RooStringVar(stateName,stateName,physName)) ;
00609
00610
00611 physName = strtok(0,spaceChars) ;
00612 if (first) {
00613 first = kFALSE ;
00614 } else if (physCat==0) {
00615 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification, only the first model will be used" << endl ;
00616 break ;
00617 }
00618 }
00619 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of physics models " << physModelSet << endl ;
00620
00621
00622
00623
00624 TList splitStateList ;
00625 RooArgSet splitCatSet ;
00626
00627 delete[] buf ;
00628 buflen = strlen(((RooStringVar*)buildConfig.find("splitCats"))->getVal())+1 ;
00629 buf = new char[buflen] ;
00630 strlcpy(buf,((RooStringVar*)buildConfig.find("splitCats"))->getVal(),buflen) ;
00631
00632 char *catName = strtok(buf,spaceChars) ;
00633 char *stateList(0) ;
00634 while(catName) {
00635
00636
00637 char* tokenPtr(0) ;
00638 if (strchr(catName,'(')) {
00639
00640 catName = strtok_r(catName,"(",&tokenPtr) ;
00641 stateList = strtok_r(0,")",&tokenPtr) ;
00642
00643 } else {
00644 stateList = 0 ;
00645 }
00646
00647 RooCategory* splitCat = catName ? dynamic_cast<RooCategory*>(dependents.find(catName)) : 0 ;
00648 if (!splitCat) {
00649 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested split category " << (catName?catName:"(null)")
00650 << " is not a RooCategory in the dataset" << endl ;
00651 delete[] buf ;
00652 return 0 ;
00653 }
00654 splitCatSet.add(*splitCat) ;
00655
00656
00657 if (stateList) {
00658 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: splitting of category " << catName
00659 << " restricted to states (" << stateList << ")" << endl ;
00660
00661
00662 TList* slist = new TList ;
00663 slist->SetName(catName) ;
00664 splitStateList.Add(slist) ;
00665
00666 char* stateLabel = strtok_r(stateList,",",&tokenPtr) ;
00667
00668 while(stateLabel) {
00669
00670 const RooCatType* type = splitCat->lookupType(stateLabel) ;
00671 if (!type) {
00672 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitCat " << splitCat->GetName()
00673 << " doesn't have a state named " << stateLabel << endl ;
00674 splitStateList.Delete() ;
00675 delete[] buf ;
00676 return 0 ;
00677 }
00678 slist->Add((TObject*)type) ;
00679
00680 stateLabel = strtok_r(0,",",&tokenPtr) ;
00681 }
00682 }
00683
00684 catName = strtok(0,spaceChars) ;
00685 }
00686 if (physCat) splitCatSet.add(*physCat) ;
00687 RooSuperCategory masterSplitCat("masterSplitCat","Master splitting category",splitCatSet) ;
00688
00689 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of splitting categories " << splitCatSet << endl ;
00690
00691
00692 RooArgSet auxSplitSet ;
00693 RooArgSet* auxSplitCloneSet(0) ;
00694 if (auxSplitCats) {
00695
00696 auxSplitCloneSet = (RooArgSet*) auxSplitCats->snapshot(kTRUE) ;
00697 if (!auxSplitCloneSet) {
00698 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf(" << GetName() << ") Couldn't deep-clone set auxiliary splitcats, abort." << endl ;
00699 delete[] buf ;
00700 return 0 ;
00701 }
00702
00703 TIterator* iter = auxSplitCats->createIterator() ;
00704 RooAbsArg* arg ;
00705 while((arg=(RooAbsArg*)iter->Next())) {
00706
00707 RooAbsArg* aux = auxSplitCats->find(arg->GetName()) ;
00708
00709
00710 if (splitCatSet.find(aux->GetName())) {
00711 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: dataset contains a fundamental splitting category " << endl
00712 << " with the same name as an auxiliary split function (" << aux->GetName() << "). " << endl
00713 << " Auxiliary split function will be ignored" << endl ;
00714 continue ;
00715 }
00716
00717
00718 RooArgSet* parSet = aux->getParameters(splitCatSet) ;
00719 if (parSet->getSize()>0) {
00720 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: ignoring auxiliary category " << aux->GetName()
00721 << " because it has servers that are not listed in splitCatSet: " << *parSet << endl ;
00722 delete parSet ;
00723 continue ;
00724 }
00725
00726
00727 aux->recursiveRedirectServers(splitCatSet) ;
00728
00729
00730 auxSplitSet.add(*aux) ;
00731 }
00732 delete iter ;
00733
00734 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of auxiliary splitting categories " << auxSplitSet << endl ;
00735 }
00736
00737
00738 TList* customizerList = new TList ;
00739
00740
00741 TIterator* physIter = physModelSet.createIterator() ;
00742 RooAbsPdf* physModel ;
00743 while((physModel=(RooAbsPdf*)physIter->Next())) {
00744 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: processing physics model " << physModel->GetName() << endl ;
00745
00746 RooCustomizer* physCustomizer = new RooCustomizer(*physModel,masterSplitCat,_splitNodeList) ;
00747 customizerList->Add(physCustomizer) ;
00748
00749
00750 RooStringVar* ruleStr = (RooStringVar*) buildConfig.find(physModel->GetName()) ;
00751 if (ruleStr) {
00752
00753 delete[] buf ;
00754 buflen = strlen(ruleStr->getVal())+1 ;
00755 buf = new char[buflen] ;
00756
00757 strlcpy(buf,ruleStr->getVal(),buflen) ;
00758 char *tokenPtr(0) ;
00759
00760 char* token = strtok_r(buf,spaceChars,&tokenPtr) ;
00761
00762 enum Mode { SplitCat, Colon, ParamList } ;
00763 Mode mode(SplitCat) ;
00764
00765 char* splitCatName ;
00766 RooAbsCategory* splitCat(0) ;
00767
00768 while(token) {
00769
00770 switch (mode) {
00771 case SplitCat:
00772 {
00773 splitCatName = token ;
00774
00775 if (strchr(splitCatName,',')) {
00776
00777
00778
00779 splitCat = (RooAbsCategory*) _compSplitCatSet.find(splitCatName) ;
00780 TString origCompCatName(splitCatName) ;
00781 if (!splitCat) {
00782
00783
00784 char *tokptr = 0;
00785 char *catName2 = strtok_r(token,",",&tokptr) ;
00786
00787 RooArgSet compCatSet ;
00788 while(catName2) {
00789 RooAbsArg* cat = splitCatSet.find(catName2) ;
00790
00791
00792 if (!cat) {
00793 cat = (RooAbsCategory*) auxSplitSet.find(catName2) ;
00794 }
00795
00796 if (!cat) {
00797 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << catName2
00798 << " not found in the primary or auxilary splitcat list" << endl ;
00799 customizerList->Delete() ;
00800 delete customizerList ;
00801
00802 splitStateList.Delete() ;
00803 delete[] buf ;
00804 return 0 ;
00805 }
00806 compCatSet.add(*cat) ;
00807
00808 catName2 = strtok_r(0,",",&tokptr) ;
00809 }
00810
00811
00812
00813
00814 TIterator* iter = compCatSet.createIterator() ;
00815 RooAbsArg* arg ;
00816 while((arg=(RooAbsArg*)iter->Next())) {
00817 RooArgSet tmp(compCatSet) ;
00818 tmp.remove(*arg) ;
00819 if (arg->dependsOnValue(tmp)) {
00820 coutE(InputArguments) << "RooSimPdfBuilder::buildPDF: ERROR: Ill defined split: auxiliary splitting category " << arg->GetName()
00821 << " used in composite split " << compCatSet << " depends on one or more of the other splitting categories in the composite split" << endl ;
00822
00823
00824 customizerList->Delete() ;
00825 delete customizerList ;
00826 splitStateList.Delete() ;
00827 delete[] buf ;
00828 return 0 ;
00829 }
00830 }
00831 delete iter ;
00832
00833 splitCat = new RooMultiCategory(origCompCatName,origCompCatName,compCatSet) ;
00834 _compSplitCatSet.addOwned(*splitCat) ;
00835
00836 }
00837 } else {
00838
00839
00840
00841 splitCat = (RooAbsCategory*) splitCatSet.find(splitCatName) ;
00842
00843
00844 if (!splitCat) {
00845 splitCat = (RooAbsCategory*) auxSplitSet.find(splitCatName) ;
00846 }
00847
00848 if (!splitCat) {
00849 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitting category "
00850 << splitCatName << " not found in the primary or auxiliary splitcat list" << endl ;
00851 customizerList->Delete() ;
00852 delete customizerList ;
00853 splitStateList.Delete() ;
00854 delete[] buf ;
00855 return 0 ;
00856 }
00857 }
00858
00859 mode = Colon ;
00860 break ;
00861 }
00862 case Colon:
00863 {
00864 if (strcmp(token,":")) {
00865 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected ':' after "
00866 << splitCat << ", found " << token << endl ;
00867 customizerList->Delete() ;
00868 delete customizerList ;
00869 splitStateList.Delete() ;
00870 delete[] buf ;
00871 return 0 ;
00872 }
00873 mode = ParamList ;
00874 break ;
00875 }
00876 case ParamList:
00877 {
00878
00879 RooArgSet splitParamList ;
00880 RooArgSet* paramList = physModel->getParameters(dependents) ;
00881
00882
00883 RooArgSet* compList = physModel->getComponents() ;
00884 paramList->add(*compList) ;
00885 delete compList ;
00886
00887 Bool_t lastCharIsComma = (token[strlen(token)-1]==',') ;
00888
00889 char *tokptr = 0 ;
00890 char *paramName = strtok_r(token,",",&tokptr) ;
00891
00892
00893 char *remainderState = 0 ;
00894 char *tokptr2 = 0 ;
00895 if (paramName && strtok_r(paramName,"[",&tokptr2)) {
00896 remainderState = strtok_r(0,"]",&tokptr2) ;
00897 }
00898
00899 while(paramName) {
00900
00901
00902 if (remainderState) {
00903 if (!splitCat->lookupType(remainderState)) {
00904 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split of parameter "
00905 << paramName << " has invalid remainder state name: " << remainderState << endl ;
00906 delete paramList ;
00907 customizerList->Delete() ;
00908 delete customizerList ;
00909 splitStateList.Delete() ;
00910 delete[] buf ;
00911 return 0 ;
00912 }
00913 }
00914
00915 RooAbsArg* param = paramList->find(paramName) ;
00916 if (!param) {
00917 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName
00918 << " is not a parameter of physics model " << physModel->GetName() << endl ;
00919 delete paramList ;
00920 customizerList->Delete() ;
00921 delete customizerList ;
00922 splitStateList.Delete() ;
00923 delete[] buf ;
00924 return 0 ;
00925 }
00926 splitParamList.add(*param) ;
00927
00928
00929 if (remainderState) {
00930
00931
00932 if (!dynamic_cast<RooAbsReal*>(param)) {
00933 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split requested of non-real valued parameter "
00934 << param->GetName() << endl ;
00935 delete paramList ;
00936 customizerList->Delete() ;
00937 delete customizerList ;
00938 splitStateList.Delete() ;
00939 delete[] buf ;
00940 return 0 ;
00941 }
00942
00943
00944 TList* remStateSplitList = static_cast<TList*>(splitStateList.FindObject(splitCat->GetName())) ;
00945
00946
00947 if (remStateSplitList && !remStateSplitList->FindObject(remainderState)) {
00948 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName
00949 << " remainder state " << remainderState << " in parameter split "
00950 << param->GetName() << " is not actually being built" << endl ;
00951 delete paramList ;
00952 customizerList->Delete() ;
00953 delete customizerList ;
00954 splitStateList.Delete() ;
00955 delete[] buf ;
00956 return 0 ;
00957 }
00958
00959 TIterator* iter = splitCat->typeIterator() ;
00960 RooCatType* type ;
00961 RooArgList fracLeafList ;
00962 TString formExpr("1") ;
00963 Int_t i(0) ;
00964
00965 while((type=(RooCatType*)iter->Next())) {
00966
00967
00968 if (!TString(type->GetName()).CompareTo(remainderState)) continue ;
00969
00970
00971 if (remStateSplitList && !remStateSplitList->FindObject(type->GetName())) {
00972 continue ;
00973 }
00974
00975
00976 TString splitLeafName(param->GetName()) ;
00977 splitLeafName.Append("_") ;
00978 splitLeafName.Append(type->GetName()) ;
00979
00980
00981 RooAbsArg* splitLeaf = _splitNodeList.find(splitLeafName) ;
00982 if (!splitLeaf) {
00983
00984 splitLeaf = (RooAbsArg*) param->clone(splitLeafName) ;
00985 _splitNodeList.add(*splitLeaf) ;
00986 _splitNodeListOwned.addOwned(*splitLeaf) ;
00987 }
00988 fracLeafList.add(*splitLeaf) ;
00989 formExpr.Append(Form("-@%d",i++)) ;
00990 }
00991 delete iter ;
00992
00993
00994 TString remLeafName(param->GetName()) ;
00995 remLeafName.Append("_") ;
00996 remLeafName.Append(remainderState) ;
00997
00998
00999 if (!_splitNodeList.find(remLeafName)) {
01000 RooAbsArg* remLeaf = new RooFormulaVar(remLeafName,formExpr,fracLeafList) ;
01001 _splitNodeList.add(*remLeaf) ;
01002 _splitNodeListOwned.addOwned(*remLeaf) ;
01003 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: creating remainder fraction formula for " << remainderState
01004 << " specialization of split parameter " << param->GetName() << " " << formExpr << endl ;
01005 }
01006 }
01007
01008
01009 paramName = strtok_r(0,",",&tokptr) ;
01010 if (paramName && strtok_r(paramName,"[",&tokptr2)) {
01011 remainderState = strtok_r(0,"]",&tokptr2) ;
01012 }
01013 }
01014
01015
01016 physCustomizer->splitArgs(splitParamList,*splitCat) ;
01017
01018 delete paramList ;
01019
01020 if (!lastCharIsComma) mode = SplitCat ;
01021 break ;
01022 }
01023 }
01024
01025 token = strtok_r(0,spaceChars,&tokenPtr) ;
01026
01027 }
01028 if (mode!=SplitCat) {
01029 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected "
01030 << (mode==Colon?":":"parameter list") << " after " << (token?token:"(null)") << endl ;
01031 }
01032
01033
01034 } else {
01035 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: no splitting rules for " << physModel->GetName() << endl ;
01036 }
01037 }
01038
01039 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: configured customizers for all physics models" << endl ;
01040 if (oodologI((TObject*)0,ObjectHandling)) {
01041 customizerList->Print() ;
01042 }
01043
01044
01045 RooArgSet fitCatList ;
01046 if (physCat) fitCatList.add(*physCat) ;
01047 fitCatList.add(splitCatSet) ;
01048 TIterator* fclIter = fitCatList.createIterator() ;
01049 RooSuperCategory *fitCat = new RooSuperCategory("fitCat","fitCat",fitCatList) ;
01050
01051
01052 RooSimultaneous* simPdf = new RooSimultaneous("simPdf","simPdf",*fitCat) ;
01053
01054
01055 TIterator* fcIter = fitCat->typeIterator() ;
01056
01057 RooCatType* fcState ;
01058 while((fcState=(RooCatType*)fcIter->Next())) {
01059
01060 fitCat->setLabel(fcState->GetName()) ;
01061
01062
01063 fclIter->Reset() ;
01064 RooAbsCategory* splitCat ;
01065 Bool_t select(kTRUE) ;
01066 while((splitCat=(RooAbsCategory*)fclIter->Next())) {
01067
01068 TList* slist = (TList*) splitStateList.FindObject(splitCat->GetName()) ;
01069 if (!slist) continue ;
01070 RooCatType* type = (RooCatType*) slist->FindObject(splitCat->getLabel()) ;
01071 if (!type) {
01072 select = kFALSE ;
01073 }
01074 }
01075 if (!select) continue ;
01076
01077
01078
01079 RooCustomizer* physCustomizer ;
01080 if (physCat) {
01081 RooStringVar* physNameVar = (RooStringVar*) stateMap.find(physCat->getLabel()) ;
01082 if (!physNameVar) continue ;
01083 physCustomizer = (RooCustomizer*) customizerList->FindObject(physNameVar->getVal());
01084 } else {
01085 physCustomizer = (RooCustomizer*) customizerList->First() ;
01086 }
01087
01088 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: Customizing physics model " << physCustomizer->GetName()
01089 << " for mode " << fcState->GetName() << endl ;
01090
01091
01092 RooAbsPdf* fcPdf = (RooAbsPdf*) physCustomizer->build(masterSplitCat.getLabel(),verbose) ;
01093 simPdf->addPdf(*fcPdf,fcState->GetName()) ;
01094 }
01095 delete fcIter ;
01096
01097
01098 _retiredCustomizerList.AddAll(customizerList) ;
01099 delete customizerList ;
01100
01101 delete fclIter ;
01102 splitStateList.Delete() ;
01103
01104 if (auxSplitCloneSet) delete auxSplitCloneSet ;
01105 delete physIter ;
01106
01107 delete[] buf ;
01108 _simPdfList.push_back(simPdf) ;
01109 _fitCatList.push_back(fitCat) ;
01110 return simPdf ;
01111 }
01112
01113
01114
01115
01116
01117
01118 RooSimPdfBuilder::~RooSimPdfBuilder()
01119 {
01120 _retiredCustomizerList.Delete() ;
01121
01122 std::list<RooSimultaneous*>::iterator iter = _simPdfList.begin() ;
01123 while(iter != _simPdfList.end()) {
01124 delete *iter ;
01125 ++iter ;
01126 }
01127
01128 std::list<RooSuperCategory*>::iterator iter2 = _fitCatList.begin() ;
01129 while(iter2 != _fitCatList.end()) {
01130 delete *iter2 ;
01131 ++iter2 ;
01132 }
01133
01134 }
01135
01136