31 #include <RooFitResult.h> 32 #include <RooRealVar.h> 33 #include <RooAbsPdf.h> 34 #include <RooDataSet.h> 54 if(PDF==
NULL || SeparationVariable==
NULL)
59 for(
unsigned int i=0;
i < Regions.size();
i++)
62 cout<<
"Integrating over Range: "<<Regions[
i].RegionName<<
" from "<<Regions[
i].min<<
" to "<<Regions[
i].max<<endl;
63 intPDF = PDF->createIntegral(*SeparationVariable,
Range(Regions[
i].RegionName.c_str()));
64 yield += intPDF->getVal();
66 cout<<
"Current yield: "<<yield<<endl;
73 histo->SetName(name.c_str());
74 histo->SetTitle(title.c_str());
75 if(axis_label ==
"GeV/c^2")
76 axis_label =
"Mass (" + axis_label +
")";
77 if(axis_label ==
"GeV/c")
78 axis_label =
"Momentum (" + axis_label +
")";
79 histo->GetXaxis()->SetTitle(axis_label.c_str());
83 if(Data==
NULL || SeparationVariable==
NULL)
85 cerr <<
"ERROR: Data or SeparationVariable is NULL returning now!\n";
88 TH1F* SideBandHist = (TH1F*)BaseHistos[index]->Clone();
89 setHistOptions(SideBandHist,(
string)variable->GetName()+
"Sideband",(
string)SideBandHist->GetTitle() +
" Sideband",(
string)variable->getUnit());
91 TH1F* SignalHist = (TH1F*)BaseHistos[index]->Clone();
92 setHistOptions(SignalHist,(
string)variable->GetName()+
"SignalHist",(
string)SignalHist->GetTitle() +
" Raw Signal",(
string)variable->getUnit());
98 TIterator* iter = (TIterator*) Data->get()->createIterator();
100 RooRealVar *sep_var=
NULL;
101 while((var = (RooAbsArg*)iter->Next()))
103 if((
string)var->GetName()==(
string)SeparationVariable->GetName())
105 sep_var = (RooRealVar*)var;
109 for(
int i=0;
i < Data->numEntries();
i++)
112 Double_t
value = variable->getVal();
113 Double_t cutval = sep_var->getVal();
115 for(
unsigned int j=0; j < SideBandRegions.size(); j++)
117 if(cutval > SideBandRegions[j].
min && cutval < SideBandRegions[j].
max)
118 SideBandHist->Fill(value);
120 for(
unsigned int j=0; j < SignalRegions.size(); j++)
122 if(cutval > SignalRegions[j].
min && cutval < SignalRegions[j].
max)
123 SignalHist->Fill(value);
127 SignalHist->Sumw2(); SideBandHist->Sumw2();
129 RawHistos.push_back(*SignalHist);
131 SignalHist->Add(SideBandHist, -stsratio);
133 SignalHist->SetTitle(((
string)BaseHistos[index]->GetTitle() +
" SBS Signal").c_str());
135 SBSHistos.push_back(*SignalHist);
137 SideBandHistos.push_back(*SideBandHist);
139 if(SideBandHist)
delete SideBandHist;
144 TCanvas genericCanvas;
146 outname = outname +
".eps";
147 genericCanvas.SaveAs(outname.c_str());
148 outname.replace(outname.size()-3,3,
"gif");
149 genericCanvas.SaveAs(outname.c_str());
155 if(SeparationVariable==
NULL)
157 cerr <<
"ERROR: printResults, SeparationVariable is NULL!\n";
161 for(
unsigned int i=0;
i < RawHistos.size(); ++
i)
163 filename=prefix +
"Raw_" + (
string)RawHistos[
i].GetName();
166 for(
unsigned int i=0;
i < SBSHistos.size(); ++
i)
168 filename=prefix +
"SBS_" + (
string)RawHistos[
i].GetName();
173 for(
unsigned int i=0;
i < SideBandHistos.size(); ++
i)
175 filename=prefix +
"SideBand_" + (
string)RawHistos[
i].GetName();
180 string outname = prefix + (
string)SeparationVariable->GetName() +
"_fitted.eps";
183 RooPlot *SepVarFrame = SeparationVariable->frame();
184 Data->plotOn(SepVarFrame);
185 ModelPDF->plotOn(SepVarFrame);
188 Canvas.SaveAs(outname.c_str());
189 outname.replace(outname.size()-3,3,
"gif");
190 Canvas.SaveAs(outname.c_str());
193 cerr <<
"ERROR: printResults, Data or ModelPDF is NULL!\n";
195 string result_outname = prefix +
"_fit_results.txt";
199 cout <<
"ERROR: Could not open file for writing!\n";
205 #if ROOT_VERSION_CODE <= ROOT_VERSION(5,19,0) 208 fit_result->printMultiline(
output,kTRUE);
219 TFile
output(outname.c_str(),
"UPDATE");
227 TIter nextkey(
output.GetListOfKeys());
229 TDirectory* curDir=
NULL;
230 while((key=(TKey*)nextkey.Next()))
235 TObject *
obj = key->ReadObj();
236 if(obj->IsA()->InheritsFrom(
"TDirectory"))
237 dirname=obj->GetName();
244 curDir =
output.mkdir(
"run0",
"Run 0");
251 Int_t run_num = dirname.Atoi();
254 curDir =
output.mkdir(dirname.Data(),(
"Run "+
stringify(run_num)).c_str());
255 output.cd(dirname.Data());
258 curDir =
output.GetDirectory(
"",kTRUE,
"");
266 for(
unsigned int i=0;
i < SideBandHistos.size(); ++
i)
268 SideBandHistos[
i].SetDirectory(curDir);
269 SideBandHistos[
i].Write();
271 for(
unsigned int i=0;
i < RawHistos.size(); ++
i)
273 RawHistos[
i].SetDirectory(curDir);
274 RawHistos[
i].Write();
276 for(
unsigned int i=0;
i < SBSHistos.size(); ++
i)
278 SBSHistos[
i].SetDirectory(curDir);
279 SBSHistos[
i].Write();
281 if(Data!=
NULL && ModelPDF!=
NULL && BackgroundPDF!=
NULL && SeparationVariable!=
NULL)
283 RooPlot *sep_varFrame = SeparationVariable->frame();
284 Data->plotOn(sep_varFrame);
285 ModelPDF->plotOn(sep_varFrame);
286 BackgroundPDF->plotOn(sep_varFrame);
287 sep_varFrame->Write();
290 cerr <<
"ERROR: saveResults, did not save RooPlot of data and fit\n";
303 cerr <<
"ERROR: print_plot, Data or ModelPDF are NULL\n";
307 RooPlot *genericFrame = printVar->frame();
308 Data->plotOn(genericFrame);
309 ModelPDF->plotOn(genericFrame);
310 TCanvas genericCanvas;
311 genericFrame->Draw();
312 outname = outname +
".eps";
313 genericCanvas.SaveAs(outname.c_str());
314 outname.replace(outname.size()-3,3,
"gif");
315 genericCanvas.SaveAs(outname.c_str());
321 SeparationVariable(0),
330 SignalSidebandRatio(0)
335 RooAbsPdf *bkg_shape,
338 const vector<TH1F*>&
base,
446 cout <<
"Beginning SideBand Subtraction\n";
454 cerr <<
"ERROR: doGobalFit, no ModelPDF, SeparationVariable or Data specified\n";
464 cout <<
"Finished fitting background!\n";
473 TIterator* iter = (TIterator*)
Data->get()->createIterator();
475 while((variable = (RooAbsArg*)iter->Next()))
487 if(iter)
delete iter;
492 Int_t binMin = Total.FindBin(leftRegion.
max,0.0,0.0);
493 Int_t binMax = Total.FindBin(leftRegion.
min,0.0,0.0);
494 double numLeft = Total.Integral( binMin, binMax );
496 binMin = Total.FindBin(rightRegion.
max,0.0,0.0);
497 binMax = Total.FindBin(rightRegion.
min,0.0,0.0);
498 double numRight = Total.Integral( binMin, binMax );
500 const unsigned int nbinsx = Total.GetNbinsX();
501 const double x1 = (leftRegion.
max + leftRegion.
min)/2.0;
502 const double x2 = (rightRegion.
max + rightRegion.
min)/2.0;
504 const double y1 = numLeft/(leftRegion.
max - leftRegion.
min);
505 const double y2 = numRight/(rightRegion.
max - rightRegion.
min);
507 const double Slope = (y2 - y1)/(x2 - x1);
511 TF1
function(
"sbs_function_line",
"[0]*x + [1]",Total.GetMinimum(), Total.GetMaximum());
512 for (
unsigned int binx = 1; binx <= nbinsx; ++binx )
514 double binWidth = Total.GetBinWidth(binx);
515 function.SetParameter(0,Slope*binWidth);
516 function.SetParameter(1,Intercept*binWidth);
518 double xx = Total.GetBinCenter(binx);
519 double cu = Total.GetBinContent(binx) -
function.Eval(xx);
521 double error1 = Total.GetBinError(binx);
523 Result.SetBinContent(binx, cu);
524 Result.SetBinError(binx, error1);
526 Result.SetEntries(Result.Integral() );
RooAbsPdf * BackgroundPDF
Double_t getYield(const std::vector< SbsRegion > &Regions, RooAbsPdf *PDF)
RooFitResult * getFitResult()
std::vector< TH1F > SideBandHistos
std::vector< TH1F > getSBSHistos()
void saveResults(std::string outname)
RooFitResult * fit_result
Double_t SignalSidebandRatio
std::vector< TH1F > getRawHistos()
static void print_histo(TH1F *plot, string outname)
std::vector< TH1F * > getBaseHistos()
int doSubtraction(RooRealVar *variable, Double_t stsratio, Int_t index)
std::string stringify(const T &t)
static void setHistOptions(TH1F *histo, string name, string title, string axis_label)
std::vector< TH1F > SBSHistos
PixelRecoRange< float > Range
RooRealVar * SeparationVariable
std::vector< SbsRegion > SideBandRegions
base
Make Sure CMSSW is Setup ##.
void printResults(std::string prefix="")
void addSignalRegion(Double_t min, Double_t max)
std::vector< SbsRegion > SignalRegions
char data[epos_bytes_allocation]
void addSideBandRegion(Double_t min, Double_t max)
void print_plot(RooRealVar *printVar, std::string outname)
void doFastSubtraction(TH1F &Total, TH1F &Result, SbsRegion &leftRegion, SbsRegion &rightRegion)
void setDataSet(RooDataSet *newData)
std::vector< TH1F * > BaseHistos
std::vector< TH1F > RawHistos