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 #include <iostream>
00027 #include <fstream>
00028
00029 #include "RooStats/HLFactory.h"
00030 #include "TFile.h"
00031 #include "TObject.h"
00032 #include "TObjArray.h"
00033 #include "TObjString.h"
00034
00035 #include "RooSimultaneous.h"
00036
00037 ClassImp(RooStats::HLFactory) ;
00038
00039
00040 using namespace RooStats;
00041 using namespace RooFit;
00042
00043
00044 HLFactory::HLFactory(const char *name,
00045 const char *fileName,
00046 bool isVerbose):
00047 TNamed(name,name),
00048 fComboCat(0),
00049 fComboBkgPdf(0),
00050 fComboSigBkgPdf(0),
00051 fComboDataset(0),
00052 fCombinationDone(false),
00053 fVerbose(isVerbose),
00054 fInclusionLevel(0),
00055 fOwnWs(true){
00056
00057
00058
00059
00060 TString wsName(name);
00061 wsName+="_ws";
00062 fWs = new RooWorkspace(wsName,true);
00063
00064 fSigBkgPdfNames.SetOwner();
00065 fBkgPdfNames.SetOwner();
00066 fDatasetsNames.SetOwner();
00067
00068
00069 fReadFile(fileName);
00070 }
00071
00072
00073 HLFactory::HLFactory(const char* name,
00074 RooWorkspace* externalWs,
00075 bool isVerbose):
00076 TNamed(name,name),
00077 fComboCat(0),
00078 fComboBkgPdf(0),
00079 fComboSigBkgPdf(0),
00080 fComboDataset(0),
00081 fCombinationDone(false),
00082 fVerbose(isVerbose),
00083 fInclusionLevel(0),
00084 fOwnWs(false){
00085
00086
00087 fWs=externalWs;
00088 fSigBkgPdfNames.SetOwner();
00089 fBkgPdfNames.SetOwner();
00090 fDatasetsNames.SetOwner();
00091
00092 }
00093
00094
00095 HLFactory::HLFactory():
00096 TNamed("hlfactory","hlfactory"),
00097 fComboCat(0),
00098 fComboBkgPdf(0),
00099 fComboSigBkgPdf(0),
00100 fComboDataset(0),
00101 fCombinationDone(false),
00102 fVerbose(false),
00103 fInclusionLevel(0),
00104 fOwnWs(true){
00105
00106 fWs = new RooWorkspace("hlfactory_ws",true);
00107
00108 fSigBkgPdfNames.SetOwner();
00109 fBkgPdfNames.SetOwner();
00110 fDatasetsNames.SetOwner();
00111
00112 }
00113
00114
00115 HLFactory::~HLFactory(){
00116
00117
00118 if (fComboSigBkgPdf!=NULL)
00119 delete fComboSigBkgPdf;
00120 if (fComboBkgPdf!=NULL)
00121 delete fComboBkgPdf;
00122 if (fComboDataset!=NULL)
00123 delete fComboDataset;
00124 if (fComboCat!=NULL)
00125 delete fComboCat;
00126
00127 if (fOwnWs)
00128 delete fWs;
00129 }
00130
00131
00132 int HLFactory::AddChannel(const char* label,
00133 const char* SigBkgPdfName,
00134 const char* BkgPdfName,
00135 const char* DatasetName){
00136
00137
00138
00139
00140
00141
00142
00143 if (fCombinationDone){
00144 std::cerr << "Cannot add anymore channels. "
00145 << "Combination already carried out.\n";
00146 return -1;
00147 }
00148
00149 if (SigBkgPdfName!=0){
00150 if (fWs->pdf(SigBkgPdfName)==NULL){
00151 std::cerr << "Pdf " << SigBkgPdfName << " not found in workspace!\n";
00152 return -1;
00153 }
00154 TObjString* name = new TObjString(SigBkgPdfName);
00155 fSigBkgPdfNames.Add(name);
00156 }
00157
00158 if (BkgPdfName!=0){
00159 if (fWs->pdf(BkgPdfName)==NULL){
00160 std::cerr << "Pdf " << BkgPdfName << " not found in workspace!\n";
00161 return -1;
00162 }
00163 TObjString* name = new TObjString(BkgPdfName);
00164 fBkgPdfNames.Add(name);
00165 }
00166
00167 if (DatasetName!=0){
00168 if (fWs->data(DatasetName)==NULL){
00169 std::cerr << "Dataset " << DatasetName << " not found in workspace!\n";
00170 return -1;
00171 }
00172 TObjString* name = new TObjString(DatasetName);
00173 fDatasetsNames.Add(name);
00174 }
00175
00176 if (label!=0){
00177 TObjString* name = new TObjString(label);
00178 fLabelsNames.Add(name);
00179 }
00180 return 0;
00181
00182 }
00183
00184
00185 RooAbsPdf* HLFactory::GetTotSigBkgPdf(){
00186
00187
00188
00189 if (fSigBkgPdfNames.GetSize()==0)
00190 return 0;
00191
00192 if (fComboSigBkgPdf!=NULL)
00193 return fComboSigBkgPdf;
00194
00195 if (!fNamesListsConsistent())
00196 return NULL;
00197
00198 if (fSigBkgPdfNames.GetSize()==1){
00199 TString name(((TObjString*)fSigBkgPdfNames.At(0))->String());
00200 fComboSigBkgPdf=fWs->pdf(name);
00201 return fComboSigBkgPdf;
00202 }
00203
00204 if (!fCombinationDone)
00205 fCreateCategory();
00206
00207 RooArgList pdfs("pdfs");
00208
00209 TIterator* it=fSigBkgPdfNames.MakeIterator();
00210 TObjString* ostring;
00211 TObject* obj;
00212 while ((obj = it->Next())){
00213 ostring=(TObjString*) obj;
00214 pdfs.add( *(fWs->pdf(ostring->String())) );
00215 }
00216 delete it;
00217
00218 TString name(GetName());
00219 name+="_sigbkg";
00220
00221 TString title(GetName());
00222 title+="_sigbkg";
00223
00224 fComboSigBkgPdf=
00225 new RooSimultaneous(name,
00226 title,
00227 pdfs,
00228 *fComboCat);
00229
00230 return fComboSigBkgPdf;
00231
00232 }
00233
00234 RooAbsPdf* HLFactory::GetTotBkgPdf(){
00235
00236
00237
00238
00239 if (fBkgPdfNames.GetSize()==0)
00240 return 0;
00241
00242 if (fComboBkgPdf!=NULL)
00243 return fComboBkgPdf;
00244
00245 if (!fNamesListsConsistent())
00246 return NULL;
00247
00248 if (fBkgPdfNames.GetSize()==1){
00249 fComboBkgPdf=fWs->pdf(((TObjString*)fBkgPdfNames.First())->String());
00250 return fComboBkgPdf;
00251 }
00252
00253 if (!fCombinationDone)
00254 fCreateCategory();
00255
00256 RooArgList pdfs("pdfs");
00257
00258 TIterator* it = fBkgPdfNames.MakeIterator();
00259 TObjString* ostring;
00260 TObject* obj;
00261 while ((obj = it->Next())){
00262 ostring=(TObjString*) obj;
00263 pdfs.add( *(fWs->pdf(ostring->String())) );
00264 }
00265
00266 TString name(GetName());
00267 name+="_bkg";
00268
00269 TString title(GetName());
00270 title+="_bkg";
00271
00272 fComboBkgPdf=
00273 new RooSimultaneous(name,
00274 title,
00275 pdfs,
00276 *fComboCat);
00277
00278 return fComboBkgPdf;
00279
00280 }
00281
00282
00283 RooDataSet* HLFactory::GetTotDataSet(){
00284
00285
00286
00287
00288 if (fDatasetsNames.GetSize()==0)
00289 return 0;
00290
00291 if (fComboDataset!=NULL)
00292 return fComboDataset;
00293
00294 if (!fNamesListsConsistent())
00295 return NULL;
00296
00297 if (fDatasetsNames.GetSize()==1){
00298 fComboDataset=(RooDataSet*)fWs->data(((TObjString*)fDatasetsNames.First())->String());
00299 return fComboDataset;
00300 }
00301
00302 if (!fCombinationDone)
00303 fCreateCategory();
00304
00305 TIterator* it = fDatasetsNames.MakeIterator();
00306 TObjString* ostring;
00307 TObject* obj = it->Next();
00308 ostring = (TObjString*) obj;
00309 fComboDataset = (RooDataSet*) fWs->data(ostring->String()) ;
00310 fComboDataset->Print();
00311 TString dataname(GetName());
00312 fComboDataset = new RooDataSet(*fComboDataset,dataname+"_TotData");
00313 int catindex=0;
00314 fComboCat->setIndex(catindex);
00315 fComboDataset->addColumn(*fComboCat);
00316 while ((obj = it->Next())){
00317 ostring=(TObjString*) obj;
00318 catindex++;
00319 RooDataSet* dummy = new RooDataSet(*(RooDataSet*)fWs->data(ostring->String()),"");
00320 fComboCat->setIndex(catindex);
00321 fComboCat->Print();
00322 dummy->addColumn(*fComboCat);
00323 fComboDataset->append(*dummy);
00324 delete dummy;
00325 }
00326
00327 delete it;
00328 return fComboDataset;
00329
00330 }
00331
00332
00333 RooCategory* HLFactory::GetTotCategory(){
00334
00335
00336
00337 if (fComboCat!=NULL)
00338 return fComboCat;
00339
00340 if (!fNamesListsConsistent())
00341 return NULL;
00342
00343 if (!fCombinationDone)
00344 fCreateCategory();
00345
00346 return fComboCat;
00347
00348 }
00349
00350
00351 int HLFactory::ProcessCard(const char* filename){
00352
00353 return fReadFile(filename,0);
00354 }
00355
00356
00357 int HLFactory::fReadFile(const char*fileName, bool is_included){
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 if (is_included)
00386 fInclusionLevel+=1;
00387 else
00388 fInclusionLevel=0;
00389
00390 const int maxDeepness=50;
00391 if (fInclusionLevel>maxDeepness){
00392 TString warning("The inclusion stack is deeper than ");
00393 warning+=maxDeepness;
00394 warning+=". Is this a recursive inclusion?";
00395 Warning("fReadFile", "%s", warning.Data());
00396 }
00397
00398
00399
00400 std::ifstream ifile(fileName);
00401
00402 if(ifile.fail()){
00403 TString error("File ");
00404 error+=fileName;
00405 error+=" could not be opened.";
00406 Error("fReadFile", "%s", error.Data());
00407 return -1;
00408 }
00409
00410 TString ifileContent("");
00411 ifileContent.ReadFile(ifile);
00412 ifile.close();
00413
00414
00415
00416 TString ifileContentStripped("");
00417
00418 TObjArray* lines_array = ifileContent.Tokenize("\n");
00419 TIterator* lineIt=lines_array->MakeIterator();
00420
00421 bool in_comment=false;
00422 TString line;
00423 TObject* line_o;
00424
00425 while((line_o=(*lineIt)())){
00426 line = (static_cast<TObjString*>(line_o))->GetString();
00427
00428
00429 if (in_comment)
00430 if (line.EndsWith("*/")){
00431 in_comment=false;
00432 if (fVerbose) Info("fReadFile","Out of multiline comment ...");
00433
00434 continue;
00435 }
00436
00437
00438
00439 if ((line.BeginsWith("/*") && line.EndsWith("*/")) ||
00440 line.BeginsWith("//")){
00441 if (fVerbose) Info("fReadFile","In single line comment ...");
00442 continue;
00443 }
00444
00445
00446 if (line.BeginsWith("/*")){
00447 in_comment=true;
00448 if (fVerbose) Info("fReadFile","In multiline comment ...");
00449 continue;
00450 }
00451
00452 ifileContentStripped+=line+"\n";
00453 }
00454
00455 delete lines_array;
00456 delete lineIt;
00457
00458
00459
00460 lines_array = ifileContentStripped.Tokenize(";");
00461 lineIt=lines_array->MakeIterator();
00462 in_comment=false;
00463
00464 const int nNeutrals=2;
00465 TString neutrals[nNeutrals]={"\t"," "};
00466
00467 while((line_o=(*lineIt)())){
00468
00469 line = (static_cast<TObjString*>(line_o))->GetString();
00470
00471
00472 line.Strip(TString::kBoth,' ');
00473
00474
00475 line.ReplaceAll("\n","");
00476
00477
00478 if (line.BeginsWith("echo")){
00479 line = line(5,line.Length()-1);
00480 if (fVerbose)
00481 std::cout << "Echoing line " << line.Data() << std::endl;
00482 std::cout << "[" << GetName() << "] echo: "
00483 << line.Data() << std::endl;
00484 continue;
00485 }
00486
00487
00488 for (int i=0;i<nNeutrals;++i)
00489 line.ReplaceAll(neutrals[i],"");
00490
00491
00492 if (fVerbose) Info("fReadFile","Reading --> %s <--", line.Data());
00493
00494
00495 if (line == ""){
00496 if (fVerbose) Info("fReadFile", "%s", "Empty line: skipping ...");
00497 continue;
00498 }
00499
00500
00501
00502 if (line.BeginsWith("#include")){
00503 line.ReplaceAll("#include","");
00504 if (fVerbose) Info("fReadFile","Reading included file...");
00505 fReadFile(line,true);
00506 continue;
00507 }
00508
00509
00510 if (fVerbose) Info("fReadFile","Parsing the line...");
00511 fParseLine(line);
00512 }
00513
00514 delete lineIt;
00515 delete lines_array;
00516
00517 return 0;
00518 }
00519
00520
00521
00522 void HLFactory::fCreateCategory(){
00523
00524
00525
00526
00527 fCombinationDone=true;
00528
00529 TString name(GetName());
00530 name+="_category";
00531
00532 TString title(GetName());
00533 title+="_category";
00534
00535 fComboCat=new RooCategory(name,title);
00536
00537 TIterator* it=fLabelsNames.MakeIterator();
00538 TObjString* ostring;
00539 TObject* obj;
00540 while ((obj = it->Next())){
00541 ostring=(TObjString*) obj;
00542 fComboCat->defineType(ostring->String());
00543 }
00544
00545 }
00546
00547
00548 bool HLFactory::fNamesListsConsistent(){
00549
00550
00551
00552 if ((fSigBkgPdfNames.GetEntries()==fBkgPdfNames.GetEntries() || fBkgPdfNames.GetEntries()==0) &&
00553 (fSigBkgPdfNames.GetEntries()==fDatasetsNames.GetEntries() || fDatasetsNames.GetEntries()==0) &&
00554 (fSigBkgPdfNames.GetEntries()==fLabelsNames.GetEntries() || fLabelsNames.GetEntries()==0))
00555 return true;
00556 else{
00557 std::cerr << "The number of datasets and models added as channels "
00558 << " is not the same!\n";
00559 return false;
00560 }
00561 }
00562
00563
00564 int HLFactory::fParseLine(TString& line){
00565
00566
00567 if (fVerbose) Info("fParseLine", "Parsing line: %s", line.Data());
00568
00569 TString new_line("");
00570
00571 const int nequals = line.CountChar('=');
00572
00573
00574
00575 if (line.Contains("::") ||
00576 nequals==0 ||
00577 (line.Contains("[") &&
00578 line.Contains("]") &&
00579 nequals>0 &&
00580 ! line.Contains("(") &&
00581 ! line.Contains(")"))) {
00582 fWs->factory(line);
00583 return 0;
00584 }
00585
00586
00587 if (nequals==1 ||
00588 (nequals > 1 && line.Contains("SIMUL"))){
00589
00590
00591
00592 const int equal_index=line.First('=');
00593 const int par_index=line.First('(');
00594 TString o_name(line(0,equal_index));
00595 TString o_class(line(equal_index+1,par_index-equal_index-1));
00596 TString o_descr(line(par_index+1,line.Length()-par_index-2));
00597
00598 if (fVerbose) Info("fParseLine", "o_name=%s o_class=%s o_descr=%s",
00599 o_name.Data(), o_class.Data(), o_descr.Data());
00600
00601
00602
00603 if (o_class =="import"){
00604
00605
00606
00607 TObjArray* descr_array = o_descr.Tokenize(",");
00608
00609 const int n_descr_parts=descr_array->GetEntries();
00610
00611 if (n_descr_parts<2 || n_descr_parts>3)
00612 Error("fParseLine","Import wrong syntax: cannot process %s", o_descr.Data());
00613
00614 TString obj_name (static_cast<TObjString*>(descr_array->At(n_descr_parts-1))->GetString());
00615 TString ws_name("");
00616 TString rootfile_name (static_cast<TObjString*>(descr_array->At(0))->GetString());
00617
00618 TFile* ifile=TFile::Open(rootfile_name);
00619 if (ifile==0)
00620 return 1;
00621
00622 if (n_descr_parts==3){
00623 o_descr.ReplaceAll(",",":");
00624 fWs->import(o_descr);
00625 }
00626 else if(n_descr_parts==2){
00627 if (fVerbose)
00628 Info("fParseLine","Importing %s from %s under the name of %s",
00629 obj_name.Data(), rootfile_name.Data(), o_name.Data());
00630 TObject* the_obj=ifile->Get(obj_name);
00631 fWs->import(*the_obj,o_name);
00632 }
00633 delete ifile;
00634 return 0;
00635 }
00636
00637 new_line=o_class+"::"+o_name+"("+o_descr+")";
00638
00639 if (fVerbose){
00640 std::cout << "DEBUG: line: " << line.Data() << std::endl;
00641 std::cout << "DEBUG: new_line: " << new_line.Data() << std::endl;
00642 }
00643
00644 fWs->factory(new_line);
00645
00646 return 0;
00647 }
00648
00649 else {
00650 fWs->factory(line);
00651 }
00652
00653 return 0;
00654
00655 }
00656
00657
00658
00659