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