/* below is the content of the file testUnfold7binning.xml, which is required as input to run this example. */ #include #include #include #include #include #include #include #include #include "TUnfoldBinningXML.h" using std::cout; void testUnfold7b() { // switch on histogram errors TH1::SetDefaultSumw2(); //======================================================= // Step 1: open file to save histograms and binning schemes TFile *outputFile=new TFile("testUnfold7_histograms.root","recreate"); //======================================================= // Step 2: read binning from XML // and save them to output file TUnfoldBinning *fineBinningRoot,*coarseBinningRoot; // read binning schemes in XML format TDOMParser parser; TString dir = gSystem->UnixPathName(gSystem->GetDirName(__FILE__)); Int_t error = parser.ParseFile(dir+"/testUnfold7binning.xml"); if(error) { cout<<"error="<Write(); coarseBinningRoot->Write(); if(fineBinningRoot) { fineBinningRoot->PrintStream(cout); } else { cout<<"could not read 'detector' binning\n"; } if(coarseBinningRoot) { coarseBinningRoot->PrintStream(cout); } else { cout<<"could not read 'generator' binning\n"; } TUnfoldBinning const *fineBinning=fineBinningRoot;//->FindNode("ptfine"); TUnfoldBinning const *coarseBinning=coarseBinningRoot;//->FindNode("ptcoarse"); //======================================================= // Step 3: book and fill data histograms Float_t ptRec[3],ptGen[3],weight; Int_t isTriggered,isSignal; outputFile->cd(); TH1 *histDataRecF=fineBinning->CreateHistogram("histDataRecF"); TH1 *histDataRecC=coarseBinning->CreateHistogram("histDataRecC"); TH1 *histDataBgrF=fineBinning->CreateHistogram("histDataBgrF"); TH1 *histDataBgrC=coarseBinning->CreateHistogram("histDataBgrC"); TH1 *histDataGen=coarseBinning->CreateHistogram("histDataGen"); for (auto histo : {histDataRecF, histDataRecC, histDataBgrF, histDataBgrC, histDataGen}) { histo->SetDirectory(outputFile); } TFile *dataFile=new TFile("testUnfold7_data.root"); TTree *dataTree=(TTree *) dataFile->Get("data"); if(!dataTree) { cout<<"could not read 'data' tree\n"; } dataTree->ResetBranchAddresses(); dataTree->SetBranchAddress("ptrec",ptRec); //dataTree->SetBranchAddress("discr",&discr); // for real data, only the triggered events are available dataTree->SetBranchAddress("istriggered",&isTriggered); // data truth parameters dataTree->SetBranchAddress("ptgen",ptGen); dataTree->SetBranchAddress("issignal",&isSignal); dataTree->SetBranchStatus("*",1); cout<<"loop over data events\n"; #define VAR_REC (ptRec[2]) #define VAR_GEN (ptGen[2]) for(Int_t ievent=0;ieventGetEntriesFast();ievent++) { if(dataTree->GetEntry(ievent)<=0) break; // fill histogram with reconstructed quantities if(isTriggered) { int binF=fineBinning->GetGlobalBinNumber(VAR_REC); int binC=coarseBinning->GetGlobalBinNumber(VAR_REC); histDataRecF->Fill(binF); histDataRecC->Fill(binC); if(!isSignal) { histDataBgrF->Fill(binF); histDataBgrC->Fill(binC); } } // fill histogram with data truth parameters if(isSignal) { int binGen=coarseBinning->GetGlobalBinNumber(VAR_GEN); histDataGen->Fill(binGen); } } delete dataTree; delete dataFile; //======================================================= // Step 4: book and fill histogram of migrations // it receives events from both signal MC and background MC TH2 *histMcsigGenRecF=TUnfoldBinning::CreateHistogramOfMigrations (coarseBinning,fineBinning,"histMcsigGenRecF"); TH2 *histMcsigGenRecC=TUnfoldBinning::CreateHistogramOfMigrations (coarseBinning,coarseBinning,"histMcsigGenRecC"); TH1 *histMcsigRecF=fineBinning->CreateHistogram("histMcsigRecF"); TH1 *histMcsigRecC=coarseBinning->CreateHistogram("histMcsigRecC"); TH1 *histMcsigGen=coarseBinning->CreateHistogram("histMcsigGen"); for (auto histo : std::initializer_list{histMcsigGenRecF, histMcsigGenRecC, histMcsigRecF, histMcsigRecC, histMcsigGen}) { histo->SetDirectory(outputFile); } TFile *signalFile=new TFile("testUnfold7_signal.root"); TTree *signalTree=(TTree *) signalFile->Get("signal"); if(!signalTree) { cout<<"could not read 'signal' tree\n"; } signalTree->ResetBranchAddresses(); signalTree->SetBranchAddress("ptrec",&ptRec); //signalTree->SetBranchAddress("discr",&discr); signalTree->SetBranchAddress("istriggered",&isTriggered); signalTree->SetBranchAddress("ptgen",&ptGen); signalTree->SetBranchAddress("weight",&weight); signalTree->SetBranchStatus("*",1); cout<<"loop over MC signal events\n"; for(Int_t ievent=0;ieventGetEntriesFast();ievent++) { if(signalTree->GetEntry(ievent)<=0) break; int binC=0,binF=0; if(isTriggered) { binF=fineBinning->GetGlobalBinNumber(VAR_REC); binC=coarseBinning->GetGlobalBinNumber(VAR_REC); } int binGen=coarseBinning->GetGlobalBinNumber(VAR_GEN); histMcsigGenRecF->Fill(binGen,binF,weight); histMcsigGenRecC->Fill(binGen,binC,weight); histMcsigRecF->Fill(binF,weight); histMcsigRecC->Fill(binC,weight); histMcsigGen->Fill(binGen,weight); } delete signalTree; delete signalFile; TH1 *histMcbgrRecF=fineBinning->CreateHistogram("histMcbgrRecF"); TH1 *histMcbgrRecC=coarseBinning->CreateHistogram("histMcbgrRecC"); histMcbgrRecF->SetDirectory(outputFile); histMcbgrRecC->SetDirectory(outputFile); TFile *bgrFile=new TFile("testUnfold7_background.root"); TTree *bgrTree=(TTree *) bgrFile->Get("background"); if(!bgrTree) { cout<<"could not read 'background' tree\n"; } bgrTree->ResetBranchAddresses(); bgrTree->SetBranchAddress("ptrec",&ptRec); //bgrTree->SetBranchAddress("discr",&discr); bgrTree->SetBranchAddress("istriggered",&isTriggered); bgrTree->SetBranchAddress("weight",&weight); bgrTree->SetBranchStatus("*",1); cout<<"loop over MC background events\n"; for(Int_t ievent=0;ieventGetEntriesFast();ievent++) { if(bgrTree->GetEntry(ievent)<=0) break; // here, for background only reconstructed quantities are known // and only the reconstructed events are relevant if(isTriggered) { int binF=fineBinning->GetGlobalBinNumber(VAR_REC); int binC=coarseBinning->GetGlobalBinNumber(VAR_REC); histMcbgrRecF->Fill(binF,weight); histMcbgrRecC->Fill(binC,weight); } } delete bgrTree; delete bgrFile; outputFile->Write(); delete outputFile; }