FairRoot/PandaRoot
Functions
combinePlotsPhi.C File Reference

Go to the source code of this file.

Functions

void confgraph (TGraph *g, TString tit, int col=1, int marker=20, double shift=0.)
 
TH1F * createHistoGraph (TGraph *g, TString tit="", double xmin=0, double xmax=0)
 
void combinePlotsPhi (TString fname="anaPhi.root")
 

Function Documentation

void combinePlotsPhi ( TString  fname = "anaPhi.root")

Definition at line 65 of file combinePlotsPhi.C.

References c1, confgraph(), createHistoGraph(), f, h, i, and TString.

66 {
67  gStyle->SetOptStat(0);
68  gStyle->SetOptFit(0);
69 
70  TFile *f=new TFile(fname);
71 
72  TString plotnam[2]={"phi_A","phi_B"};
73  TString lab[4] = { "signal efficiency;EMC missing (polar #theta) [%];signal efficiency [%]", "background efficiency;EMC missing (polar #theta) [%];background efficiency [%]", "signal to noise;EMC missing (polar #theta) [%];S/B", "significance;EMC missing (polar #theta) [%];significance [#sigma]"};
74 
75  TGraphErrors *geffs[2];
76  TGraphErrors *geffb[2];
77  TGraphErrors *gston[2];
78  TGraphErrors *gsign[2];
79 
80  TCanvas *c1=new TCanvas("c1","c1",10,10,1300,1000);
81  c1->Divide(2,2,0.0001,0.0001);
82 
83  int colors[2] = {kBlue,kCyan-2};
84  //int colors[2] = {2,kRed-7};
85 
86  TH1F *h[4]={0};
87 
88  for (int i=0;i<2;++i)
89  {
90  c1->cd(1);
91  geffs[i]=(TGraphErrors*)f->Get(Form("effs_%s",plotnam[i].Data()));
92  confgraph(geffs[i],"",colors[i],20+(i%2)*3,i);
93  if (h[0]==0) h[0]=createHistoGraph(geffs[i],lab[0]);
94  if (i==0) h[0]->Draw();
95  geffs[i]->Draw("P same");
96 
97  c1->cd(2); //gPad->SetLogy();
98  geffb[i]=(TGraphErrors*)f->Get(Form("effb_%s",plotnam[i].Data()));
99  confgraph(geffb[i],"",colors[i],20+(i%2)*3,i);
100  if (h[1]==0) h[1]=createHistoGraph(geffb[i],lab[1]);
101  //h[1]->SetMinimum(1e-8);
102  //h[1]->SetMaximum(0.0002);
103  if (i==0) h[1]->Draw();
104  geffb[i]->Draw("P same");
105 
106  c1->cd(3);
107  gston[i]=(TGraphErrors*)f->Get(Form("sn_%s",plotnam[i].Data()));
108  confgraph(gston[i],"",colors[i],20+(i%2)*3,i);
109  if (h[2]==0) h[2]=createHistoGraph(gston[i],lab[2]);
110  if (i==0) h[2]->Draw();
111  gston[i]->Draw("P same");
112 
113  c1->cd(4);
114  gsign[i]=(TGraphErrors*)f->Get(Form("sign_%s",plotnam[i].Data()));
115  confgraph(gsign[i],"",colors[i],20+(i%2)*3,i);
116  if (h[3]==0) h[3]=createHistoGraph(gsign[i],lab[3]);
117  if (i==0) h[3]->Draw();
118  gsign[i]->Draw("P same");
119  }
120 
121  c1->cd(1);
122  TLegend *leg1=new TLegend(0.16,0.18,0.5,0.34);
123  leg1->AddEntry(geffs[0],"#phi#phi - Setup A","lep");
124  leg1->AddEntry(geffs[1],"#phi#phi - Setup B","lep");
125  leg1->Draw();
126 
127  c1->cd(2);
128  TLegend *leg2=new TLegend(0.16,0.18,0.5,0.34);
129  leg2->AddEntry(geffb[0],"#phi#phi - Setup A","lep");
130  leg2->AddEntry(geffb[1],"#phi#phi - Setup B","lep");
131  leg2->Draw();
132 
133  c1->cd(3);
134  TLegend *leg3=new TLegend(0.16,0.18,0.5,0.34);
135  leg3->AddEntry(gston[0],"#phi#phi - Setup A","lep");
136  leg3->AddEntry(gston[1],"#phi#phi - Setup B","lep");
137  leg3->Draw();
138 
139  c1->cd(4);
140  TLegend *leg4=new TLegend(0.16,0.18,0.5,0.34);
141  leg4->AddEntry(gsign[0],"#phi#phi - Setup A","lep");
142  leg4->AddEntry(gsign[1],"#phi#phi - Setup B","lep");
143  leg4->Draw();
144 
145  c1->SaveAs("fig/comb_phi.gif");
146 
147 }
Int_t i
Definition: run_full.C:25
TH1F * createHistoGraph(TGraph *g, TString tit="", double xmin=0, double xmax=0)
Definition: analyse_etac1.C:60
void confgraph(TGraph *g, TString tit, int col=1, int marker=20)
Definition: analyse_etac1.C:47
TFile * f
Definition: bump_analys.C:12
c1
Definition: plot_dirc.C:35
void confgraph ( TGraph *  g,
TString  tit,
int  col = 1,
int  marker = 20,
double  shift = 0. 
)

Definition at line 4 of file combinePlotsPhi.C.

References col, f1, i, x, and y.

5 {
6  if (tit!="") g->GetHistogram()->SetTitle(tit);
7  g->GetHistogram()->SetMinimum(0);
8  g->SetLineColor(col);
9  g->SetMarkerColor(col);
10  g->SetMarkerStyle(marker);
11  g->SetMarkerSize(1.5);
12  g->SetLineWidth(2);
13 
14  TF1 *f1 = g->GetFunction("f1");
15 
16  if (f1)
17  {
18  f1->SetLineColor(col);
19  f1->SetLineStyle(2);
20  }
21 
22  for (int i=0;i<g->GetN();++i)
23  {
24  double x,y;
25  g->GetPoint(i,x,y);
26  g->SetPoint(i,x+shift,y);
27 
28  // cout <<x<<" "<<y<<endl;
29  }
30 }
Int_t i
Definition: run_full.C:25
TF1 * f1
Definition: reco_analys2.C:50
int col
Definition: anaLmdDigi.C:67
TFile * g
Double_t x
Double_t y
TH1F* createHistoGraph ( TGraph *  g,
TString  tit = "",
double  xmin = 0,
double  xmax = 0 
)

Definition at line 34 of file combinePlotsPhi.C.

References cnt, h, xmax, and xmin.

35 {
36  static int cnt=0;
37  if (tit=="") tit=g->GetTitle();
38 
39  double ymax = TMath::MaxElement(g->GetN(), g->GetY());
40  double dymax = g->GetErrorY(TMath::LocMax(g->GetN(), g->GetEY()));
41 
42  if (xmin>=xmax)
43  {
44  xmin = TMath::MinElement(g->GetN(), g->GetX());
45  xmax = TMath::MaxElement(g->GetN(), g->GetX());
46  }
47  TH1F *h=new TH1F(Form("h%03d",cnt++),tit,g->GetN(),xmin-0.5*(xmax-xmin)/g->GetN(), xmax+0.5*(xmax-xmin)/g->GetN());
48 
49  h->GetXaxis()->SetLabelSize(0.05);
50  h->GetXaxis()->SetTitleSize(0.05);
51  h->GetYaxis()->SetLabelSize(0.05);
52  h->GetYaxis()->SetTitleSize(0.05);
53  h->GetYaxis()->SetTitleOffset(1.4);
54 
55  h->SetMaximum((ymax+dymax)*1.05);
56  h->SetMinimum(0);
57 
58  h->SetStats(0);
59 
60  return h;
61 }
Double_t xmax
TFile * g
Int_t cnt
Definition: hist-t7.C:106
Double_t xmin