Hall-D Software  alpha
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
bcal_inv_mass.C
Go to the documentation of this file.
1 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_300
2 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_500
3 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_700
4 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_900
5 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_v_E
6 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_v_z_lowE
7 // hnamepath: /bcal_inv_mass/bcal_diphoton_mass_v_z_highE
8 
9 {
10 
11  TDirectory *dir = (TDirectory*)gDirectory->FindObjectAny("bcal_inv_mass");
12  if(dir) dir->cd();
13 
14  TH1I* bcal_diphoton_mass_300 = (TH1I*)gDirectory->FindObjectAny("bcal_diphoton_mass_300");
15  TH1I* bcal_diphoton_mass_500 = (TH1I*)gDirectory->FindObjectAny("bcal_diphoton_mass_500");
16  TH1I* bcal_diphoton_mass_700 = (TH1I*)gDirectory->FindObjectAny("bcal_diphoton_mass_700");
17  TH1I* bcal_diphoton_mass_900 = (TH1I*)gDirectory->FindObjectAny("bcal_diphoton_mass_900");
18  TH2I* bcal_diphoton_mass_v_E = (TH2I*)gDirectory->FindObjectAny("bcal_diphoton_mass_v_E");
19  TH2I* bcal_diphoton_mass_v_z_lowE = (TH2I*)gDirectory->FindObjectAny("bcal_diphoton_mass_v_z_lowE");
20  TH2I* bcal_diphoton_mass_v_z_highE = (TH2I*)gDirectory->FindObjectAny("bcal_diphoton_mass_v_z_highE");
21 
22  int polnumber = 3;
23  float fit_low = 0.07;
24  float fit_high = 0.20;
25  float par_300[15];
26  float par_500[15];
27  float par_700[15];
28  float par_900[15];
29 
30  TCanvas *c1 = NULL;
31  if(gPad == NULL){
32  c1 = new TCanvas( "c1", "BCAL_inv_mass_plot", 800, 800 );
33  c1->cd(0);
34  c1->Draw();
35  c1->Update();
36  } else {
37  c1 = gPad->GetCanvas();
38  }
39 
40  if( !gPad ) return;
41  c1->Divide(2,2);
42 
43  if( bcal_diphoton_mass_300 ){
44 
45  bcal_diphoton_mass_300->SetStats(0);
46  double max_300 = bcal_diphoton_mass_300->GetMaximum();
47  c1->cd(1);
48  bcal_diphoton_mass_300->Draw();
49  bcal_diphoton_mass_300->SetLineWidth(2);
50  TF1 *fitfunc_300 = new TF1("fitfunc_300",Form("gaus(0)+pol%i(3)",polnumber),fit_low,fit_high);
51  fitfunc_300->SetParameters(max_300/2, 0.1, 0.009);
52  for (int i=0; i<=polnumber; i++) {
53  fitfunc_300->SetParameter(3+i,par_300[i]);
54  }
55  fitfunc_300->SetParNames("height","mean","sigma");
56  fitfunc_300->SetParLimits(0,max_300/5.,max_300+max_300/10.);
57  fitfunc_300->SetParLimits(1,0.06,0.18);
58  fitfunc_300->SetParLimits(2,0.006,0.02);
59  fitfunc_300->SetLineWidth(2);
60  bcal_diphoton_mass_300->Fit(fitfunc_300,"RQ");
61  for (int i=0; i<3+polnumber; i++) {
62  par_300[i] = fitfunc_300->GetParameter(i);
63  }
64  fitfunc_300->Draw("same");
65  TPaveText *pt_300 = new TPaveText(0.6, 0.65, 0.99, 0.89, "NDC");
66  pt_300->SetFillColor(0);
67  pt_300->AddText(Form("M_{#pi^{0}} = %.3f MeV",par_300[1]*1000));
68  pt_300->AddText(Form("#sigma_{#pi^{0}} = %.3f MeV",par_300[2]*1000));
69  pt_300->AddText(Form("#sigma/M = %.3f %%",(par_300[2]/par_300[1])*100));
70  pt_300->Draw();
71  }
72  if( bcal_diphoton_mass_500 ){
73  bcal_diphoton_mass_500->SetStats(0);
74  double max_500 = bcal_diphoton_mass_500->GetMaximum();
75  c1->cd(2);
76  bcal_diphoton_mass_500->Draw();
77  bcal_diphoton_mass_500->SetLineWidth(2);
78  TF1 *fitfunc_500 = new TF1("fitfunc_500",Form("gaus(0)+pol%i(3)",polnumber),fit_low,fit_high);
79  fitfunc_500->SetParameters(max_500/5, 0.135, 0.01);
80  for (int i=0; i<=polnumber; i++) {
81  fitfunc_500->SetParameter(3+i,par_500[i]);
82  }
83  fitfunc_500->SetParNames("height","mean","sigma");
84  fitfunc_500->SetParLimits(0,0.,max_500+max_500/10.);
85  fitfunc_500->SetParLimits(1,0.02,0.30);
86  fitfunc_500->SetParLimits(2,0.001,0.050);
87  fitfunc_500->SetLineWidth(2);
88  bcal_diphoton_mass_500->Fit(fitfunc_500,"RQ");
89  for (int i=0; i<3+polnumber; i++) {
90  par_500[i] = fitfunc_500->GetParameter(i);
91  }
92  fitfunc_500->Draw("same");
93 
94  TPaveText *pt_500 = new TPaveText(0.6, 0.65, 0.99, 0.89, "NDC");
95  pt_500->SetFillColor(0);
96  pt_500->AddText(Form("M_{#pi^{0}} = %.3f MeV",par_500[1]*1000));
97  pt_500->AddText(Form("#sigma_{#pi^{0}} = %.3f MeV",par_500[2]*1000));
98  pt_500->AddText(Form("#sigma/M = %.3f %%",(par_500[2]/par_500[1])*100));
99  pt_500->Draw();
100  }
101  if( bcal_diphoton_mass_700 ){
102 
103  bcal_diphoton_mass_700->SetStats(0);
104  double max_700 = bcal_diphoton_mass_700->GetMaximum();
105  c1->cd(3);
106  bcal_diphoton_mass_700->Draw();
107  bcal_diphoton_mass_700->SetLineWidth(2);
108  TF1 *fitfunc_700 = new TF1("fitfunc_700",Form("gaus(0)+pol%i(3)",polnumber),fit_low,fit_high);
109  fitfunc_700->SetParameters(max_700/5, 0.135, 0.01);
110  for (int i=0; i<=polnumber; i++) {
111  fitfunc_700->SetParameter(3+i,par_700[i]);
112  }
113  fitfunc_700->SetParNames("height","mean","sigma");
114  fitfunc_700->SetParLimits(0,0.,max_700+max_700/10.);
115  fitfunc_700->SetParLimits(1,0.02,0.30);
116  fitfunc_700->SetParLimits(2,0.001,0.050);
117  fitfunc_700->SetLineWidth(2);
118  bcal_diphoton_mass_700->Fit(fitfunc_700,"RQ");
119  for (int i=0; i<3+polnumber; i++) {
120  par_700[i] = fitfunc_700->GetParameter(i);
121  }
122  fitfunc_700->Draw("same");
123 
124  TPaveText *pt_700 = new TPaveText(0.6, 0.65, 0.99, 0.89, "NDC");
125  pt_700->SetFillColor(0);
126  pt_700->AddText(Form("M_{#pi^{0}} = %.3f MeV",par_700[1]*1000));
127  pt_700->AddText(Form("#sigma_{#pi^{0}} = %.3f MeV",par_700[2]*1000));
128  pt_700->AddText(Form("#sigma/M = %.3f %%",(par_700[2]/par_700[1])*100));
129  pt_700->Draw();
130  }
131  if( bcal_diphoton_mass_900 ){
132 
133  bcal_diphoton_mass_900->SetStats(0);
134  double max_900 = bcal_diphoton_mass_900->GetMaximum();
135  c1->cd(4);
136  bcal_diphoton_mass_900->Draw();
137  bcal_diphoton_mass_900->SetLineWidth(2);
138  TF1 *fitfunc_900 = new TF1("fitfunc_900",Form("gaus(0)+pol%i(3)",polnumber),fit_low,fit_high);
139  fitfunc_900->SetParameters(max_900/5, 0.135, 0.01);
140  for (int i=0; i<=polnumber; i++) {
141  fitfunc_900->SetParameter(3+i,par_900[i]);
142  }
143  fitfunc_900->SetParNames("height","mean","sigma");
144  fitfunc_900->SetParLimits(0,0.,max_900+max_900/10.);
145  fitfunc_900->SetParLimits(1,0.02,0.30);
146  fitfunc_900->SetParLimits(2,0.001,0.050);
147  fitfunc_900->SetLineWidth(2);
148  bcal_diphoton_mass_900->Fit(fitfunc_900,"RQ");
149  for (int i=0; i<3+polnumber; i++) {
150  par_900[i] = fitfunc_900->GetParameter(i);
151  }
152  fitfunc_900->Draw("same");
153 
154  TPaveText *pt_900 = new TPaveText(0.6, 0.65, 0.99, 0.89, "NDC");
155  pt_900->SetFillColor(0);
156  pt_900->AddText(Form("M_{#pi^{0}} = %.3f MeV",par_900[1]*1000));
157  pt_900->AddText(Form("#sigma_{#pi^{0}} = %.3f MeV",par_900[2]*1000));
158  pt_900->AddText(Form("#sigma/M = %.3f %%",(par_900[2]/par_900[1])*100));
159  pt_900->Draw();
160  }
161 
162 }
static TH1I * bcal_diphoton_mass_900
static TH1I * bcal_diphoton_mass_500
Double_t c1[2][NMODULES]
Definition: tw_corr.C:68
static TH2I * bcal_diphoton_mass_v_z_highE
static TH2I * bcal_diphoton_mass_v_E
static TH1I * bcal_diphoton_mass_300
static TH2I * bcal_diphoton_mass_v_z_lowE
TDirectory * dir
Definition: bcal_hist_eff.C:31
static TH1I * bcal_diphoton_mass_700