FairRoot/PandaRoot
PndDiscTaskReconstruction.cxx
Go to the documentation of this file.
1 //-------------------------------------------------------------------------
2 // Author: Mustafa Schmidt (Mustafa.A.Schmidt@physik.uni-giessen.de)
3 // Changes:
4 // Date: 30.11.2015
5 // Description: Track Reconstruction
6 //-------------------------------------------------------------------------
7 
8 
10 #include "PndDiscDigitizedHit.h"
11 #include "PndDiscParticleMCPoint.h"
12 #include "PndDiscSensorMCPoint.h"
13 #include "PndDiscReconResult.h"
14 
15 #include "FairRunAna.h"
16 #include "FairRootManager.h"
17 #include "FairEventHeader.h"
18 
19 #include "TClonesArray.h"
20 #include "TString.h"
21 #include "TRandom.h"
22 #include "TMath.h"
23 
24 #include <cstdlib>
25 
26 
27 double pos_fel_x[108] =
28 {
29  3.200016,
30  8.995571,
31  14.791126,
32  20.586681,
33  26.382236,
34  32.17779,
35  37.973345,
36  43.7689,
37  49.564455,
38  55.109197,
39  59.355735,
40  63.602272,
41  67.84881,
42  72.095348,
43  76.341886,
44  80.588423,
45  84.834961,
46  89.081499,
47  92.318168,
48  93.869134,
49  95.420099,
50  96.971064,
51  98.52203,
52  100.072995,
53  101.623961,
54  103.174926,
55  104.725891,
56  104.742558,
57  103.189644,
58  101.63673,
59  100.083816,
60  98.530901,
61  96.977987,
62  95.425073,
63  93.872158,
64  92.319244,
65  89.121382,
66  84.882642,
67  80.643902,
68  76.405162,
69  72.166422,
70  67.927682,
71  63.688942,
72  59.450202,
73  55.211462,
74  49.635225,
75  43.839148,
76  38.043072,
77  32.246995,
78  26.450918,
79  20.654841,
80  14.858764,
81  9.062688,
82  3.266611,
83  -3.200016,
84  -8.995571,
85  -14.791126,
86  -20.586681,
87  -26.382236,
88  -32.17779,
89  -37.973345,
90  -43.7689,
91  -49.564455,
92  -55.109197,
93  -59.355735,
94  -63.602272,
95  -67.84881,
96  -72.095348,
97  -76.341886,
98  -80.588423,
99  -84.834961,
100  -89.081499,
101  -92.318168,
102  -93.869134,
103  -95.420099,
104  -96.971064,
105  -98.52203,
106  -100.072995,
107  -101.623961,
108  -103.174926,
109  -104.725891,
110  -104.742558,
111  -103.189644,
112  -101.63673,
113  -100.083816,
114  -98.530901,
115  -96.977987,
116  -95.425073,
117  -93.872158,
118  -92.319244,
119  -89.121382,
120  -84.882642,
121  -80.643902,
122  -76.405162,
123  -72.166422,
124  -67.927682,
125  -63.688942,
126  -59.450202,
127  -55.211462,
128  -49.635225,
129  -43.839148,
130  -38.043072,
131  -32.246995,
132  -26.450918,
133  -20.654841,
134  -14.858764,
135  -9.062688,
136  -3.266611
137 };
138 
139 double pos_fel_y[108] =
140 {
141  104.742558,
142  103.189644,
143  101.63673,
144  100.083816,
145  98.530901,
146  96.977987,
147  95.425073,
148  93.872158,
149  92.319244,
150  89.121382,
151  84.882642,
152  80.643902,
153  76.405162,
154  72.166422,
155  67.927682,
156  63.688942,
157  59.450202,
158  55.211462,
159  49.635225,
160  43.839148,
161  38.043072,
162  32.246995,
163  26.450918,
164  20.654841,
165  14.858764,
166  9.062688,
167  3.266611,
168  -3.200016,
169  -8.995571,
170  -14.791126,
171  -20.586681,
172  -26.382236,
173  -32.17779,
174  -37.973345,
175  -43.7689,
176  -49.564455,
177  -55.109197,
178  -59.355735,
179  -63.602272,
180  -67.84881,
181  -72.095348,
182  -76.341886,
183  -80.588423,
184  -84.834961,
185  -89.081499,
186  -92.318168,
187  -93.869134,
188  -95.420099,
189  -96.971064,
190  -98.52203,
191  -100.072995,
192  -101.623961,
193  -103.174926,
194  -104.725891,
195  -104.742558,
196  -103.189644,
197  -101.63673,
198  -100.083816,
199  -98.530901,
200  -96.977987,
201  -95.425073,
202  -93.872158,
203  -92.319244,
204  -89.121382,
205  -84.882642,
206  -80.643902,
207  -76.405162,
208  -72.166422,
209  -67.927682,
210  -63.688942,
211  -59.450202,
212  -55.211462,
213  -49.635225,
214  -43.839148,
215  -38.043072,
216  -32.246995,
217  -26.450918,
218  -20.654841,
219  -14.858764,
220  -9.062688,
221  -3.266611,
222  3.200016,
223  8.995571,
224  14.791126,
225  20.586681,
226  26.382236,
227  32.17779,
228  37.973345,
229  43.7689,
230  49.564455,
231  55.109197,
232  59.355735,
233  63.602272,
234  67.84881,
235  72.095348,
236  76.341886,
237  80.588423,
238  84.834961,
239  89.081499,
240  92.318168,
241  93.869134,
242  95.420099,
243  96.971064,
244  98.52203,
245  100.072995,
246  101.623961,
247  103.174926,
248  104.725891
249 };
250 
251 
252 //Angles of tilted Sensors
253 double angle_sensor[108] =
254 {
255  75,
256  75,
257  75,
258  75,
259  75,
260  75,
261  75,
262  75,
263  75,
264  45,
265  45,
266  45,
267  45,
268  45,
269  45,
270  45,
271  45,
272  45,
273  14.98,
274  14.98,
275  14.98,
276  14.98,
277  14.98,
278  14.98,
279  14.98,
280  14.98,
281  14.98,
282  345,
283  345,
284  345,
285  345,
286  345,
287  345,
288  345,
289  345,
290  345,
291  315,
292  315,
293  315,
294  315,
295  315,
296  315,
297  315,
298  315,
299  315,
300  284.98,
301  284.98,
302  284.98,
303  284.98,
304  284.98,
305  284.98,
306  284.98,
307  284.98,
308  284.98,
309  255,
310  255,
311  255,
312  255,
313  255,
314  255,
315  255,
316  255,
317  255,
318  225.05,
319  225.05,
320  225.05,
321  225.05,
322  225.05,
323  225.05,
324  225.05,
325  225.05,
326  225.05,
327  194.98,
328  194.98,
329  194.98,
330  194.98,
331  194.98,
332  194.98,
333  194.98,
334  194.98,
335  194.98,
336  165,
337  165,
338  165,
339  165,
340  165,
341  165,
342  165,
343  165,
344  165,
345  135.05,
346  135.05,
347  135.05,
348  135.05,
349  135.05,
350  135.05,
351  135.05,
352  135.05,
353  135.05,
354  104.98,
355  104.98,
356  104.98,
357  104.98,
358  104.98,
359  104.98,
360  104.98,
361  104.98,
362  104.98,
363 };
364 
365 
366 //Rest masses of pions, kaons and protons in MeV
367 double E0[3] = {134.976, 493.677, 938.272};
368 
369 //---------------------------------------------------------------------------
370 
371 double PndDiscTaskReconstruction::gauss(int x, double mean_value, double rms)
372 {
373  return 1/sqrt(2*TMath::Pi()*rms*rms)*exp(-0.5*((x-mean_value)/rms)*((x-mean_value)/rms));
374 }
375 
376 //---------------------------------------------------------------------------
377 
378 double PndDiscTaskReconstruction::mean(std::vector<double> values)
379 {
380  double average = 0;
381 
382  for(size_t i = 0; i < values.size(); i++)
383  {
384  average += values[i];
385  }
386 
387  return average/double(values.size());
388 }
389 
390 
391 //---------------------------------------------------------------------------
392 
393 double PndDiscTaskReconstruction::deviation(double mean_value, std::vector<double> values)
394 {
395  double average = 0;
396 
397  for(size_t i = 0; i < values.size(); i++)
398  {
399  average += (values[i]-mean_value)*(values[i]-mean_value);
400  }
401 
402  average = sqrt(1/(double(values.size())-1)*average);
403 
404  return average;
405 }
406 
407 
408 //---------------------------------------------------------------------------
409 
410 PndDiscTaskReconstruction::PndDiscTaskReconstruction() : PndPersistencyTask("DiscDircTaskReconstruction"), tclarr_digits(NULL), tclarr_particles(NULL), tclarr_tracks(NULL), is_time_based(kTRUE), flag_export_patterns(kFALSE), average_wl(420.), minimum_wl(385.)
411 {
412  branch_name_digits = "DiscDigit";
413  folder_name_digits = "DiscDIRC";
414 }
415 
416 //---------------------------------------------------------------------------
417 
419 {
420 }
421 
422 
423 //--------------------------------------------------------------------------
424 
425 
427 {
428  LOG(INFO) << "DiscDircTaskReconstruction::ReInit()";
429  return kSUCCESS;
430 }
431 
432 //--------------------------------------------------------------------------
433 
434 
436 {
437  LOG(INFO) << "PndDiscTaskReconstruction::Init()";
438 
439  // Get IO manager instance:
440  FairRootManager* io_manager = FairRootManager::Instance();
441  if(!io_manager)
442  {
443  LOG(FATAL) << "FairRootManager instance is NULL !!!";
444  return kFATAL;
445  }
446 
447  // Get branch with digitized hits:
448  tclarr_digits = (TClonesArray*) io_manager->GetObject(branch_name_digits);
449  if(!tclarr_digits)
450  {
451  LOG(ERROR) << "Branch " << branch_name_digits.Data() << " is not accessible through FairRootManager.";
452  return kERROR;
453  }
454 
455  tclarr_particles = (TClonesArray*) io_manager->GetObject("DiscMCTruthTracks");
456  if(!tclarr_particles)
457  {
458  LOG(ERROR) << "No DiscDIRC_ParticleMCPoint collection registered with FairRootManager";
459  return kFATAL;
460  }
461 
462 
463  tclarr_recon_results = new TClonesArray("PndDiscReconResult");
464  tclarr_digits_out = new TClonesArray("PndDiscDigitizedHit");
465  tclarr_particles_out = new TClonesArray("PndDiscParticleMCPoint");
466 
467  io_manager->Register("DiscPatternPrediction", "DiscDircDetector", tclarr_recon_results, GetPersistency());
468  io_manager->Register("DiscDigitizedHit","DiscDircDetector", tclarr_digits_out, GetPersistency());
469  io_manager->Register("DiscRealTracks","DiscDircDetector", tclarr_particles_out, GetPersistency());
470 
471 
472 
473  static bool recon_initialized = false;
474  if(recon_initialized) return kSUCCESS;
475 
476  return kSUCCESS;
477 }
478 
479 
480 //----------------------------------------------------------------------
481 
482 
483 
487 {
488  //Reading out entries
489  int particles = tclarr_particles->GetEntries();
490  int entries = tclarr_digits->GetEntries();
491 
492  // time window around track time:
493  //const double time_window_t1 = -5.0; //[R.K. 01/2017] unused variable?
494  //const double time_window_t2 = 40.0; //[R.K. 01/2017] unused variable?
495 
496 
497  //Information of primary particle
498  std::vector<double> x0,y0,px,py,pz,p;
499  std::vector<double> theta,phi;
500  std::vector<double> t;
501 
502  std::vector<std::vector<double > > mean_time;
503  std::vector<std::vector<double > > deviation_time;
504 
505 
506 
507  //double average; //[R.K. 01/2017] unused variable?
508 
509  //-----------------------------------------------------------------------
510  // Reading out primary particle information
511  //-----------------------------------------------------------------------
512 
513  for(int i = 0; i < particles; i++)
514  {
516 
517  x0.push_back(particle_mc_point->GetX());
518  y0.push_back(particle_mc_point->GetY());
519 
520  px.push_back(particle_mc_point->GetPx());
521  py.push_back(particle_mc_point->GetPy());
522  pz.push_back(particle_mc_point->GetPz());
523 
524  p.push_back(sqrt(px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i]));
525 
526  t.push_back(particle_mc_point->GetTime());
527 
528  theta.push_back(acos(pz[i]/sqrt(px[i]*px[i]+py[i]*py[i]+pz[i]*pz[i])));
529 
530  std::cout << "Primary particles:" << std::endl;
531  std::cout << "Particle pos: x = "<< x0[i] << std::endl;
532  std::cout << "Particle pos: y = "<< y0[i] << std::endl;
533  std::cout << "Particle mom: px = "<< px[i] << std::endl;
534  std::cout << "Particle mom: py = "<< py[i] << std::endl;
535  std::cout << "Particle mom: pz = "<< pz[i] << std::endl;
536  std::cout << "Particle abs mom: p = "<< p[i] << std::endl;
537  std::cout << "Particle time: t = "<< t[i] << std::endl;
538  std::cout << "Particle angle: theta[deg] = " << theta[i]*180/TMath::Pi() << std::endl;
539  std::cout << std::endl;
540  }
541 
542  //-------------------------------------------------------------------------
543  // Caculating theoretical hit pattern for pion, kaon, proton
544  //-------------------------------------------------------------------------
545 
546  PndDiscReconResult *recon_result = new PndDiscReconResult();
547 
548  for(int i = 0; i < particles; i++)
549  {
550  for(int j = 0; j < 3; j++)
551  {
552  std::cout << "Particle: " << j << std::endl;
553  for(int k = 0; k < 108; k++)
554  {
555  double E = sqrt(p[i]*1000*p[i]*1000 + E0[j]*E0[j]);
556  double gamma = E/E0[j];
557  double beta = sqrt(1-1/(gamma*gamma));
558  double n = 1.47;
559  double cherenkov = acos(1/(n*beta));
560 
561  double dx = pos_fel_x[k] - x0[i];
562  double dy = pos_fel_y[k] - y0[i];
563 
564  double phirel = acos((dx*px[i]+dy*py[i])/(sqrt(dx*dx+dy*dy)*sqrt(px[i]*px[i]+py[i]*py[i])));
565 
566  double angle;
567 
568  if(div(k,27).quot == 0 || div(k,27).quot == 3)
569  {
570  angle = angle_sensor[k]*TMath::Pi()/180;
571  }
572  if(div(k,27).quot == 1 || div(k,27).quot == 2)
573  {
574  angle = 2*TMath::Pi() - angle_sensor[k]*TMath::Pi()/180;
575  }
576 
577  double alpha = acos(dx/(sqrt(dx*dx+dy*dy))) - angle;
578 
579  double delta = sin(theta[i])*cos(phirel);
580  double lambda = delta*delta+cos(theta[i])*cos(theta[i]);
581  double varphi = acos(delta*cos(cherenkov)/lambda+sqrt((cos(theta[i])*cos(theta[i])-cos(cherenkov)*cos(cherenkov))/lambda+(delta*cos(cherenkov)/lambda)*(delta*cos(cherenkov)/lambda)));
582  double varphidash = atan(tan(varphi)/cos(alpha));
583 
584  int pixel = (TMath::Pi()/2 - varphidash-0.8530931655)/0.0035321867;
585 
586  double time = (sqrt(dx*dx+dy*dy))/100/3e8*n/cos(varphidash)*1e9;
587 
588  //Results are only valid for pixels between 0 and 100
589  if(pixel >= 0 && pixel <= 100)
590  {
591 
592  std::cout << "Sensor: " << k << " Pixel: " << pixel << " Time: " << time << std::endl;
593 
594  recon_result->hypothesis = j;
595  recon_result->sensor = k;
596  recon_result->pixel = pixel;
597  recon_result->time = time;
598  recon_result->particle = i;
599 
600  new ((*tclarr_recon_results)[tclarr_recon_results->GetEntriesFast()]) PndDiscReconResult(*recon_result);
601  }
602  else
603  {
604 
605 
606  recon_result->hypothesis = j;
607  recon_result->sensor = k;
608  recon_result->pixel = -1;
609  recon_result->time = time;
610  recon_result->particle = i;
611 
612  new ((*tclarr_recon_results)[tclarr_recon_results->GetEntriesFast()]) PndDiscReconResult(*recon_result);
613  }
614  }
615  std::cout << std::endl;
616  }
617  }
618 
619 
620 
621  //-------------------------------------------------------------------------
622  // Copying digitized hits and particle information
623  //-------------------------------------------------------------------------
624 
625  for(int i = 0; i < entries; i++)
626  {
628  digits = new((*tclarr_digits_out)[tclarr_digits_out->GetEntriesFast()]) PndDiscDigitizedHit(*digits);
629  }
630 
631  for(int i = 0; i < particles; i++)
632  {
634  particles_out = new((*tclarr_particles_out)[tclarr_particles_out->GetEntriesFast()]) PndDiscParticleMCPoint(*particles_out);
635  }
636 
637 
638 
639  //Reading out hit information and calculating Cherenkov angle
640  for(int i = 0; i < entries; i++)
641  {
642  PndDiscDigitizedHit* sensor_mc_point = (PndDiscDigitizedHit*)tclarr_digits->At(i);
643 
644  int detector_id = sensor_mc_point->GetDetectorID();
645  int sensor_id = 27*detector_id + sensor_mc_point->GetReadoutID(); //Sensor ID of hit
646  int pixel = sensor_mc_point->GetPixelNumber(); // Pixel number of hit
647 
648  double time = sensor_mc_point->GetTdcTime();
649 
650  //Difference vector between particle and sensor element
651  //double dx = pos_fel_x[sensor_id] - x0[0];
652  //double dy = pos_fel_y[sensor_id] - y0[0];
653 
654  //double angle = 0; //[R.K. 01/2017] unused variable
655 
656  //if(detector_id == 0 || detector_id == 3) //[R.K. 01/2017] unused variable
657  //{
658  //angle = angle_sensor[sensor_id]*TMath::Pi()/180; //[R.K. 01/2017] unused variable
659  //}
660  //if(detector_id == 1 || detector_id == 2) //[R.K. 01/2017] unused variable
661  //{
662  //angle = 2*TMath::Pi() - angle_sensor[sensor_id]*TMath::Pi()/180; //[R.K. 01/2017] unused variable
663  //}
664 
665  //double phirel = acos((dx*px[0]+dy*py[0])/(sqrt(dx*dx+dy*dy)*sqrt(px[0]*px[0]+py[0]*py[0])));
666  //double alpha = acos(dx/(sqrt(dx*dx+dy*dy))) - angle;
667  //double varphidash = TMath::Pi()/2-(0.00347196*pixel + 0.852934);
668  //double varphi = atan(tan(varphidash)*cos(alpha));
669 
670  //double cherenkov = acos(sin(theta[0])*cos(phirel)*cos(varphi)+cos(theta[0])*sin(varphi)); //Reconstructed Cherenkov angle
671 
672  std::cout << "Hit informations:" << std::endl;
673  std::cout << "TDC Time [ns]: " << time << std::endl;
674  std::cout << "Detector ID: " << detector_id << std::endl;
675  std::cout << "Sensor ID: " << sensor_id << std::endl;
676  std::cout << "Sensor coor: x = " << pos_fel_x[sensor_id] << ", y = " << pos_fel_y[sensor_id] << std::endl;
677  std::cout << "Pixel no: " << pixel << std::endl;
678  std::cout << std::endl;
679  }
680 
681 
682  std::cout << "DiscDIRC_TaskReconstruction::Exec\n";
683 }
684 
685 
687 {
688  // called after all Tasks did their Exec() and the data is copied to the file
689  if(tclarr_recon_results != NULL) tclarr_recon_results->Clear();
690  if(tclarr_digits_out != NULL) tclarr_digits_out->Clear();
691  if(tclarr_particles_out != NULL) tclarr_particles_out->Clear();
692 
693  FinishEvents();
694 }
695 
696 
698 {
699 }
700 
701 
703 
friend F32vec4 acos(const F32vec4 &a)
Definition: P4_F32vec4.h:113
Double_t x0
Definition: checkhelixhit.C:70
Double_t p
Definition: anasim.C:58
friend F32vec4 cos(const F32vec4 &a)
Definition: P4_F32vec4.h:112
double gauss(int x, double mean, double rms)
double dy
friend F32vec4 exp(const F32vec4 &a)
Definition: P4_F32vec4.h:109
Int_t i
Definition: run_full.C:25
const Int_t & GetReadoutID() const
Double_t lambda(Double_t x, Double_t y, Double_t z)
Definition: drawdal.C:48
double angle_sensor[108]
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
double rms
Definition: plot_eta_c.C:98
friend F32vec4 sin(const F32vec4 &a)
Definition: P4_F32vec4.h:111
TString branch_name_digits
Branch name where digitized hits are stored.
int n
virtual void Exec(Option_t *opt)
TString folder_name_digits
Folder name for digits.
const Double_t & GetTdcTime() const
double mean(std::vector< double > values)
Double_t y0
Definition: checkhelixhit.C:71
TClonesArray * tclarr_particles
to cache the pointer to particle MC TClonesArray returned by IO manager.
const Int_t & GetDetectorID() const
double pos_fel_y[108]
const Int_t & GetPixelNumber() const
double dx
double E0[3]
double deviation(double mean, std::vector< double > values)
Double_t x
TClonesArray * tclarr_recon_results
results of reconstruction (pattern matching)
ClassImp(PndAnaContFact)
TTree * t
Definition: bump_analys.C:13
double alpha
Definition: f_Init.h:9
Double_t angle
double pos_fel_x[108]
Double_t Pi
double pz[39]
Definition: pipisigmas.h:14