FairRoot/PandaRoot
emc_correction_hist.C
Go to the documentation of this file.
1 // This macro produce histograms for energy-position correction of EMC clusters
2 // It reads input file with simulated-reconstructed data,
3 // create array of histograms with ratio of reconstructed
4 // and generated energy and difference of reconstructed and simulated polar angle
5 // for given energy-theta interval, then in each histogram mean value is calculated,
6 // which is stored in 2D histogram and used for correction
7 // Output file name have pattern:
8 // emc_correction_<particle>-_<version>.root, e.g. "emc_correction_gamma_1.root"
9 // Input parametr InputFile1 can be root file or text file containing list of root files
10 
11 // Define interval in energy and theta, in which output histogram is provided
12 Double_t energyIntervals[]= {0,0.1,0.2,0.3,0.4,0.5,0.6,0.7,0.8,0.9,1.0,1.1,1.2,1.5,2.0,2.5,3.0,3.5,4.0,5.0,6.0,7.0,8.0,9.0,10.0};
13 Double_t thetaIntervalsBarrel[]= {22.,23.,24.,25.,30.,35.,40.,45.,50.,55.,60.,65.,70.,75.,80.,85.,90.,95.,100.,105.,110.,115.,120.,125.,130.,135.,137., 139., 140., 141.};
14 Double_t thetaIntervalsFwd[]= {5.,6.,7.,8.,9.,10.,11.,12.,13.,14.,15.,20.,21.,22.};
15 Double_t thetaIntervalsBwd[]= {147.,148.,149.,150.,151.,152.,154.,156.,158.,160.,161.,162.,163.,164.,165.,166.,168.,170.,172.};
16 Double_t thetaIntervalsShashlyk[]= {0.,0.5,1.,1.5,2.,2.5,3.,3.5,4.,4.5,5.,5.5,6.,6.5,7.,7.5,8.,8.5,9.,9.5,10.};
17 
19 Int_t nrThetaIntervals[4]; //Barrel, Fwd, Bwd, Shashlyk
20 
21 // version is used for output file name
22 int emc_correction_hist(Int_t version, TString InputFile1, Bool_t debug=false)
23 {
24  TStopwatch timer;
25  timer.Start();
26 
27  TChain *c=new TChain("pndsim");
28  TString particle="";
29  if (InputFile1.Contains(".root"))
30  {
31  TObjArray* tmp=InputFile1.Tokenize("_");
32  TObjString *os=(TObjString *)tmp->At(2);
33  particle=os->GetString();
34  TFile *f1=new TFile(InputFile1);
35  if (f1->IsZombie())
36  {
37  std::cout<<"Input file "<<InputFile1<<" does not exists"<<std::endl;
38  abort();
39  }
40  c->Add(InputFile1);
41  }
42  else if (InputFile1.Contains(".txt"))
43  {
44  string file_name;
45  ifstream infile(InputFile1.Data(), ios::in);
46 
47  while (getline(infile,file_name, '\n'))
48  {
49  if (std::string::npos != file_name.find(".root"))
50  {
51  TFile *f2=new TFile(file_name.c_str());
52  if (f2->IsZombie())
53  {
54  std::cout<<"Input file "<<file_name<<" does not exists"<<std::endl;
55  abort();
56  }
57  c->Add(file_name.c_str());
58  if (particle=="")
59  {
60  TObjArray* tmp=TString(file_name).Tokenize("_");
61  TObjString *os=(TObjString *)tmp->At(2);
62  particle=os->GetString();
63  }
64  }
65  }
66  }
67  else
68  {
69  std::cout<<"Wrong input file: "<<InputFile1<<std::endl;
70  abort();
71  }
72 
73  TString OutputFile="emc_correction_hist_";
74  OutputFile = OutputFile +particle+"_";
75  OutputFile +=version;
76  OutputFile += ".root";
77 
80  nrThetaIntervals[1]=sizeof(thetaIntervalsFwd)/sizeof(Double_t);
81  nrThetaIntervals[2]=sizeof(thetaIntervalsBwd)/sizeof(Double_t);
83 
84  TH2F hisEnergyRatioBarrel, hisEnergyRatioFwd, hisEnergyRatioBwd, hisEnergyRatioShashlyk;
85  TH2F hisThetaDiffBarrel, hisThetaDiffFwd, hisThetaDiffBwd, hisThetaDiffShashlyk;
86 
87  TH2F hisEnergyRatioBarrel("hisEnergyRatioBarrel","Ene_MC/Ene_reco: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[0]-1),thetaIntervalsBarrel);
88  hisEnergyRatioBarrel.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
89  hisEnergyRatioBarrel.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
90 
91  TH2F hisThetaDiffBarrel("hisThetaDiffBarrel","Theta_MC - Theta_reco: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[0]-1),thetaIntervalsBarrel);
92  hisThetaDiffBarrel.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
93  hisThetaDiffBarrel.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
94 
95  TH2F hisEnergyRatioFwd("hisEnergyRatioFwd","Ene_MC/Ene_reco: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[1]-1),thetaIntervalsFwd);
96  hisEnergyRatioFwd.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
97  hisEnergyRatioFwd.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
98 
99  TH2F hisThetaDiffFwd("hisThetaDiffFwd","Theta_MC - Theta_reco: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[1]-1),thetaIntervalsFwd);
100  hisThetaDiffFwd.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
101  hisThetaDiffFwd.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
102 
103  TH2F hisEnergyRatioBwd("hisEnergyRatioBwd","Ene_MC/Ene_reco: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[2]-1),thetaIntervalsBwd);
104  hisEnergyRatioBwd.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
105  hisEnergyRatioBwd.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
106 
107  TH2F hisThetaDiffBwd("hisThetaDiffBwd","Theta_MC - Theta_reco: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[2]-1),thetaIntervalsBwd);
108  hisThetaDiffBwd.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
109  hisThetaDiffBwd.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
110 
111  TH2F hisEnergyRatioShashlyk("hisEnergyRatioShashlyk","Ene_MC/Ene_reco _Shashlyk_: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[3]-1),thetaIntervalsShashlyk);
112  hisEnergyRatioShashlyk.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
113  hisEnergyRatioShashlyk.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
114 
115  TH2F hisThetaDiffShashlyk("hisThetaDiffShashlyk","Theta_MC-Theta_reco _Shashlyk_: GetMean()",(nrEnergyIntervals-1),energyIntervals,(nrThetaIntervals[3]-1),thetaIntervalsShashlyk);
116  hisThetaDiffShashlyk.GetXaxis()->SetTitle("Cluster Reconstructed Photon Energy (GeV)");
117  hisThetaDiffShashlyk.GetYaxis()->SetTitle("Cluster Reconstructed #theta Photon Angle (#circ)");
118 
119 
120  // energy: "E_reco / E_MC"
121  TH1F hisEnergyRatio1[100][100]; // barrel
122  TH1F hisEnergyRatio2[100][100]; // forward endcap
123  TH1F hisEnergyRatio3[100][100]; // backward endcap
124  TH1F hisEnergyRatio4[100][100]; // shashlyk
125 
126  // Theta: "Theta_MC - Theta_reco"
127  TH1F hisThetaDiff1[100][100]; // barrel
128  TH1F hisThetaDiff2[100][100]; // forward endcap
129  TH1F hisThetaDiff3[100][100]; // backward endcap
130  TH1F hisThetaDiff4[100][100]; // shashlyk
131 
132  char cmd[64];
133 
134  for (Int_t i=0; i<(nrEnergyIntervals-1); i++)
135  {
136  for (Int_t j=0; j<(nrThetaIntervals[0]-1); j++)
137  {
138  sprintf(cmd,"hisEnergyRatio1_%i_%i",i,j);
139  hisEnergyRatio1[i][j].SetName(cmd);
140  hisEnergyRatio1[i][j].SetTitle(cmd);
141  hisEnergyRatio1[i][j].SetBins(100,0.,2.);
142  sprintf(cmd,"hisThetaDiff1_%i_%i",i,j);
143  hisThetaDiff1[i][j].SetName(cmd);
144  hisThetaDiff1[i][j].SetTitle(cmd);
145  hisThetaDiff1[i][j].SetBins(100,-5.,5);
146  hisEnergyRatioBarrel.SetBinContent(i+1,j+1,1.);
147  hisThetaDiffBarrel.SetBinContent(i+1,j+1,0.);
148  }
149  for (Int_t j=0; j<(nrThetaIntervals[1]-1); j++)
150  {
151  sprintf(cmd,"hisEnergyRatio2_%i_%i",i,j);
152  hisEnergyRatio2[i][j].SetName(cmd);
153  hisEnergyRatio2[i][j].SetTitle(cmd);
154  hisEnergyRatio2[i][j].SetBins(100,0.,2.);
155  sprintf(cmd,"hisThetaDiff2_%i_%i",i,j);
156  hisThetaDiff2[i][j].SetName(cmd);
157  hisThetaDiff2[i][j].SetTitle(cmd);
158  hisThetaDiff2[i][j].SetBins(100,-5.,5);
159  hisEnergyRatioFwd.SetBinContent(i+1,j+1,1.);
160  hisThetaDiffFwd.SetBinContent(i+1,j+1,0.);
161  }
162  for (Int_t j=0; j<(nrThetaIntervals[2]-1); j++)
163  {
164  sprintf(cmd,"hisEnergyRatio3_%i_%i",i,j);
165  hisEnergyRatio3[i][j].SetName(cmd);
166  hisEnergyRatio3[i][j].SetTitle(cmd);
167  hisEnergyRatio3[i][j].SetBins(100,0.,2.);
168  sprintf(cmd,"hisThetaDiff3_%i_%i",i,j);
169  hisThetaDiff3[i][j].SetName(cmd);
170  hisThetaDiff3[i][j].SetTitle(cmd);
171  hisThetaDiff3[i][j].SetBins(100,-5.,5);
172  hisEnergyRatioBwd.SetBinContent(i+1,j+1,1.);
173  hisThetaDiffBwd.SetBinContent(i+1,j+1,0.);
174  }
175  for (Int_t j=0; j<(nrThetaIntervals[3]-1); j++)
176  {
177  sprintf(cmd,"hisEnergyRatio4_%i_%i",i,j);
178  hisEnergyRatio4[i][j].SetName(cmd);
179  hisEnergyRatio4[i][j].SetTitle(cmd);
180  hisEnergyRatio4[i][j].SetBins(100,0.,2.);
181  sprintf(cmd,"hisThetaDiff4_%i_%i",i,j);
182  hisThetaDiff4[i][j].SetName(cmd);
183  hisThetaDiff4[i][j].SetTitle(cmd);
184  hisThetaDiff4[i][j].SetBins(100,-5.,5);
185  hisEnergyRatioShashlyk.SetBinContent(i+1,j+1,1.);
186  hisThetaDiffShashlyk.SetBinContent(i+1,j+1,0.);
187  }
188  }
189 
192 
193  TClonesArray* track_array=new TClonesArray("PndMCTrack");
194  c->SetBranchAddress("MCTrack",&track_array);
195  TClonesArray* cluster_array=new TClonesArray("PndEmcCluster");
196  c->SetBranchAddress("EmcCluster",&cluster_array);
197 
198 
199  Double_t cluster_energy, cluster_theta, cluster_phi; //position of the cluster
200  Double_t thMC, phiMC, enMC;
201  for (Int_t j=0; j< c->GetEntries(); j++)
202  {
203  if (j%1000==0)
204  cout<<"Event "<<j<<endl;
205 
206  c->GetEntry(j);
207 
208  PndMCTrack *track=(PndMCTrack*)track_array->At(0);
209  TLorentzVector p4mom=track->Get4Momentum();
210  TVector3 mc_momentum=track->GetMomentum();
211 
212  thMC = mc_momentum.Theta()*(180./TMath::Pi());
213  phiMC = mc_momentum.Phi()*(180./TMath::Pi());
214 
215  enMC = p4mom.E();
216 
217  // Select cluster of highest energy
218  if (cluster_array->GetEntriesFast()>0)
219  {
220  Int_t idWithHighestEnergy = 0;
221  Double_t highestEnergy = -1.;
222 
223  for (Int_t i=0; i<cluster_array->GetEntriesFast(); i++)
224  {
225  PndEmcCluster *cluster;
226  cluster=(PndEmcCluster*)cluster_array->At(i);
227  cluster_energy=cluster->energy();
228 
229  if (cluster_energy>highestEnergy)
230  {
231  idWithHighestEnergy = i;
232  highestEnergy = cluster_energy;
233  }
234  }
235 
236  // Lets analyze that cluster!
237  PndEmcCluster *cluster=(PndEmcCluster*)cluster_array->At(idWithHighestEnergy);
238  TVector3 cluster_pos=cluster->where();
239  cluster_theta=cluster_pos.Theta()*180./TMath::Pi();
240  cluster_phi=cluster_pos.Phi()*180./TMath::Pi();
241  cluster_energy=cluster->energy();
242  Int_t module=cluster->GetModule();
243 
244  Int_t thetaBin=GetThetaBin(cluster_theta, module); // bin number for cluster_theta
245  Int_t energyBin=GetEnergyBin(cluster_energy, module);// bin number for cluster_energy
246 
247  if (cluster_theta >= 141. && cluster_theta < 147.) continue; // Avoid the egges between barrel and bwendcap...
248 
249  Double_t en_div, phi_diff, theta_diff;
250  if (thetaBin<0 || energyBin<0)
251  {
252  cout<<"!!!!!!!!!!!!!!!!!!!!!!"<<endl;
253  cout<<"Module="<<module<<endl;
254  cout << "theta_reco =" << cluster_theta << "/ energy_reco =" << cluster_energy << endl;
255  cout << "thetaBin= "<<thetaBin<<" & energyBin= "<<energyBin<<endl;
256  }
257  else
258  {
259  if (cluster_theta > (thMC-5) && cluster_theta <=(thMC+5) ){
260 
261  en_div = cluster_energy/enMC;
262  phi_diff = phiMC-cluster_phi;
263  theta_diff= thMC-cluster_theta;
264 
265  if (fabs(phi_diff)<2.5 || fabs(phi_diff+360)<2.5 || fabs(phi_diff-360)<2.5){
266  switch (module)
267  {
268  case 1:
269  case 2:
270  {
271  // avoid beam pipe region
272  if ((fabs(phiMC) <5.)||fabs(phiMC-180.)<5.||fabs(phiMC-360.)<5.) break;
273  if (en_div > 0.7 && en_div <1.3){
274  hisEnergyRatio1[energyBin][thetaBin].Fill(en_div);
275  hisThetaDiff1[energyBin][thetaBin].Fill(theta_diff);
276  }
277  break;
278  }
279  case 3:
280  {
281  if (cluster_theta <6.) break; // Avoid the egges for FwEndCap...
282  if (en_div > 0.7 && en_div <1.3){
283  hisEnergyRatio2[energyBin][thetaBin].Fill(en_div);
284  hisThetaDiff2[energyBin][thetaBin].Fill(theta_diff);
285  }
286  break;
287  }
288  case 4:
289  {
290  if (en_div > 0.7 && en_div <1.3){
291  hisEnergyRatio3[energyBin][thetaBin].Fill(en_div);
292  hisThetaDiff3[energyBin][thetaBin].Fill(theta_diff);
293  }
294  break;
295  }
296  case 5:
297  {
298  if (en_div > 0.7 && en_div <1.3){
299  hisEnergyRatio4[energyBin][thetaBin].Fill(en_div);
300  hisThetaDiff4[energyBin][thetaBin].Fill(theta_diff);
301  }
302  break;
303  }
304  default:
305  {
306  abort();
307  }
308  }
309 
310  }
311  }
312  }
313  }
314  }
315 
316 
317  // Loop over energy and theta bins to GetMean() values from energy and theta histograms
318  // and write GetMean() value to the 2D histograms
319 
320  TFile fout(OutputFile,"recreate"); // output file with corrections
321  Double_t mean_th_dif, mean_en_ratio;
322 
323  for (Int_t i=0; i<(nrEnergyIntervals-1); i++)
324  {
325  for (Int_t j=0; j<(nrThetaIntervals[0]-1); j++)
326  {
327  //
328  // Energy & Theta for Barrel
329  //
330  if (hisThetaDiff1[i][j].GetEntries() == 0) {
331  mean_th_dif = 0;
332  }else{
333  mean_th_dif = hisThetaDiff1[i][j].GetMean(); // (Theta_MC - Theta_reco)
334  }
335 
336  hisThetaDiffBarrel.SetBinContent(i+1,j+1,mean_th_dif); // Theta_MC - Theta_Reco: GetMean()
337 
338  if ( hisEnergyRatio1[i][j].GetEntries() == 0) {
339  mean_en_ratio = 1;
340  }else{
341  mean_en_ratio = hisEnergyRatio1[i][j].GetMean(); // (E_MC / E_reco)
342  }
343 
344  hisEnergyRatioBarrel.SetBinContent(i+1,j+1,mean_en_ratio); // Mean
345  }
346  for (Int_t j=0; j<(nrThetaIntervals[1]-1); j++)
347  {
348  //
349  // Energy & Theta for Forward Edncap
350  //
351  if (hisThetaDiff2[i][j].GetEntries() == 0) {
352  mean_th_dif = 0;
353  }else{
354  mean_th_dif = hisThetaDiff2[i][j].GetMean(); // (Theta_MC - Theta_reco)
355  }
356 
357  hisThetaDiffFwd.SetBinContent(i+1,j+1,mean_th_dif); // Theta_MC - Theta_Reco: GetMean()
358 
359  if ( hisEnergyRatio2[i][j].GetEntries() == 0) {
360  mean_en_ratio = 1;
361  }else{
362  mean_en_ratio = hisEnergyRatio2[i][j].GetMean(); // (E_MC / E_reco)
363  }
364 
365  hisEnergyRatioFwd.SetBinContent(i+1,j+1,mean_en_ratio); // Mean
366  }
367  for (Int_t j=0; j<(nrThetaIntervals[2]-1); j++)
368  {
369  //
370  // Energy & Theta for Backward Endcap
371  //
372  if (hisThetaDiff3[i][j].GetEntries() == 0) {
373  mean_th_dif = 0;
374  }else{
375  mean_th_dif = hisThetaDiff3[i][j].GetMean(); // (Theta_MC - Theta_reco)
376  }
377 
378  hisThetaDiffBwd.SetBinContent(i+1,j+1,mean_th_dif); // Theta_MC - Theta_Reco: GetMean()
379 
380  if ( hisEnergyRatio3[i][j].GetEntries() == 0) {
381  mean_en_ratio = 1;
382  }else{
383  mean_en_ratio = hisEnergyRatio3[i][j].GetMean(); // (E_MC / E_reco)
384  }
385 
386  hisEnergyRatioBwd.SetBinContent(i+1,j+1,mean_en_ratio); // Mean
387  }
388  for (Int_t j=0; j<(nrThetaIntervals[3]-1); j++)
389  {
390  //
391  // Energy & theta for Shashlyk
392  //
393  if (hisThetaDiff4[i][j].GetEntries() == 0){
394  mean_th_dif = 0;
395  }else{
396  mean_th_dif = hisThetaDiff4[i][j].GetMean();
397  }
398 
399  hisThetaDiffShashlyk.SetBinContent(i+1,j+1,mean_th_dif); // Theta_MC - Theta_Reco: GetMean()
400 
401 
402  if (hisEnergyRatio4[i][j].GetEntries() == 0) {
403  mean_en_ratio = 1;
404  }else{
405  mean_en_ratio = hisEnergyRatio4[i][j].GetMean();
406  }
407 
408  hisEnergyRatioShashlyk.SetBinContent(i+1,j+1,mean_en_ratio);
409  }
410  }
411 
412  hisEnergyRatioBarrel.Write();
413  hisThetaDiffBarrel.Write();
414  hisEnergyRatioFwd.Write();
415  hisThetaDiffFwd.Write();
416  hisEnergyRatioBwd.Write();
417  hisThetaDiffBwd.Write();
418  hisEnergyRatioShashlyk.Write();
419  hisThetaDiffShashlyk.Write();
420 
421  if (debug)
422  {
423  for (Int_t i=0; i<(nrEnergyIntervals-1); i++)
424  {
425  for (Int_t j=0; j<(nrThetaIntervals[0]-1); j++)
426  {
427  TH1F *h1=(TH1F *)hisEnergyRatio1[i][j].Clone();
428  h1->Write();
429  TH1F *h1=(TH1F *)hisThetaDiff1[i][j].Clone();
430  h1->Write();
431  }
432  for (Int_t j=0; j<(nrThetaIntervals[1]-1); j++)
433  {
434  TH1F *h1=(TH1F *)hisEnergyRatio2[i][j].Clone();
435  h1->Write();
436  TH1F *h1=(TH1F *)hisThetaDiff2[i][j].Clone();
437  h1->Write();
438  }
439  for (Int_t j=0; j<(nrThetaIntervals[2]-1); j++)
440  {
441  TH1F *h1=(TH1F *)hisEnergyRatio3[i][j].Clone();
442  h1->Write();
443  TH1F *h1=(TH1F *)hisThetaDiff3[i][j].Clone();
444  h1->Write();
445  }
446  for (Int_t j=0; j<(nrThetaIntervals[3]-1); j++)
447  {
448  TH1F *h1=(TH1F *)hisEnergyRatio4[i][j].Clone();
449  h1->Write();
450  TH1F *h1=(TH1F *)hisThetaDiff4[i][j].Clone();
451  h1->Write();
452  }
453  }
454 
455  }
456  fout.Close();
457 
458  timer.Stop();
459  Double_t rtime = timer.RealTime();
460  Double_t ctime = timer.CpuTime();
461  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
462 
463  return 0;
464 }
465 
466 Int_t GetThetaBin(Double_t val, Int_t module)
467 {
468  Double_t *thetaIntervals;
469  Int_t nrThetaInt;
470  switch (module)
471  {
472  case 1:
473  case 2:
474  nrThetaInt=nrThetaIntervals[0];
475  thetaIntervals=thetaIntervalsBarrel;
476  break;
477  case 3:
478  nrThetaInt=nrThetaIntervals[1];
479  thetaIntervals=thetaIntervalsFwd;
480  break;
481  case 4:
482  nrThetaInt=nrThetaIntervals[2];
483  thetaIntervals=thetaIntervalsBwd;
484  break;
485  case 5:
486  nrThetaInt=nrThetaIntervals[3];
487  thetaIntervals=thetaIntervalsShashlyk;
488  break;
489  default:
490  std::cout<<"Wrong EMC module: "<<module<<std::endl;
491  abort();
492  }
493  for (Int_t i=0; i<(nrThetaInt-1); i++)
494  {
495  if (val>=thetaIntervals[i] && val<thetaIntervals[i+1])
496  {
497  return (i);
498  }
499  }
500  return -1;
501 }
502 
503 Int_t GetEnergyBin(Double_t val, Int_t module)
504 {
505  for (Int_t i=0; i<(nrEnergyIntervals-1); i++)
506  {
507  if (val>=energyIntervals[i] && val<energyIntervals[i+1])
508  return (i);
509  }
510  return -1;
511 }
Double_t thetaIntervalsBwd[]
printf("RealTime=%f seconds, CpuTime=%f seconds\n", rtime, ctime)
TVector3 where() const
Int_t i
Definition: run_full.C:25
Short_t GetModule() const
TF1 * f1
Definition: reco_analys2.C:50
PndRiemannTrack track
Definition: RiemannTest.C:33
Double_t val[nBoxes][nFEBox]
Definition: createCalib.C:11
Double_t thetaIntervalsBarrel[]
TVector3 GetMomentum() const
Definition: PndMCTrack.h:78
TLorentzVector Get4Momentum() const
Definition: PndMCTrack.cxx:102
const int particle
Double_t cluster_theta
TClonesArray * cluster_array
Double_t
TStopwatch timer
Definition: hit_dirc.C:51
Double_t cluster_energy
Int_t nrEnergyIntervals
Double_t cluster_phi
Int_t nrThetaIntervals[4]
friend F32vec4 fabs(const F32vec4 &a)
Definition: P4_F32vec4.h:47
a cluster (group of neighboring crystals) of hit emc crystals
Definition: PndEmcCluster.h:29
Double_t energyIntervals[]
Double_t ctime
Definition: hit_dirc.C:114
TString cmd
Definition: runMvdTpcDigi.C:31
TFile * f2
Double_t thetaIntervalsShashlyk[]
Int_t GetEnergyBin(Double_t val, Int_t module)
cout<<"will loop over "<< t-> GetEntries()
Definition: root2ascii.C:17
virtual Double_t energy() const
TClonesArray * track_array
Definition: plot_rk.C:34
TLorentzVector p4mom
Double_t rtime
Definition: hit_dirc.C:113
Double_t Pi
TFile infile("dedx_out.root","READ")
Double_t thetaIntervalsFwd[]
int emc_correction_hist(Int_t version, TString InputFile1, Bool_t debug=false)
string file_name
Int_t GetThetaBin(Double_t val, Int_t module)