]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFClusterFinder.cxx
f40af02b585ab26c4aeec1541b781334fa53c28c
[u/mrichter/AliRoot.git] / TOF / AliTOFClusterFinder.cxx
1 /***************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* 
17
18 Revision 0.03  2005/07/28 A. De Caro
19          Implement public method
20          Raw2Digits(Int_t, AliRawReader *)
21          to convert digits from raw data in MC digits
22          (temporary solution)
23
24 Revision 0.02  2005/07/27 A. De Caro
25          Implement public method
26          Digits2RecPoint(Int_t)
27          to convert digits in clusters
28
29 Revision 0.02  2005/07/26 A. De Caro
30          Implement private methods
31          InsertCluster(AliTOFcluster *)
32          FindClusterIndex(Double_t)
33          originally implemented in AliTOFtracker
34          by S. Arcelli and C. Zampolli
35
36 Revision 0.01  2005/07/25 A. De Caro
37          Implement public methods
38          Digits2RecPoint(AliRawReader *, TTree *)
39          Digits2RecPoint(Int_t, AliRawReader *)
40          to convert raw data in clusters
41  */
42
43 ////////////////////////////////////////////////////////////////
44 //                                                            //
45 //         Class for TOF cluster finder                       //
46 //                                                            //
47 // Starting from Raw Data, create rec points,                 //
48 //                         fill TreeR for TOF,                //
49 //                         write TOF.RecPoints.root file      //
50 //                                                            //
51 ////////////////////////////////////////////////////////////////
52
53 #include <TTree.h>
54 #include <TObjArray.h>
55 #include <TClonesArray.h>
56 #include <TFile.h>
57
58 #include "AliLog.h"
59 #include "AliRunLoader.h"
60 #include "AliLoader.h"
61
62 #include "AliTOFdigit.h"
63 #include "AliTOFcluster.h"
64 #include "AliTOFGeometry.h"
65 #include "AliTOFGeometryV4.h"
66 #include "AliTOFGeometryV5.h"
67 #include "AliTOFRawStream.h"
68 #include "AliTOFcalib.h"
69 #include "AliTOFCal.h"
70 #include "AliTOFChannel.h"
71
72 #include "AliTOFClusterFinder.h"
73
74 ClassImp(AliTOFClusterFinder)
75
76 AliTOFClusterFinder::AliTOFClusterFinder():
77   fRunLoader(0),
78   fTOFLoader(0),
79   fTreeD(0),
80   fTreeR(0),
81   fDigits(new TClonesArray("AliTOFdigit", 4000)),
82   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
83   fNumberOfTofClusters(0)
84 {
85 //
86 // Constructor
87 //
88
89   fTOFGeometry = new AliTOFGeometryV5();
90   AliInfo("V5 TOF Geometry is taken as the default");
91
92 }
93 //______________________________________________________________________________
94
95 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
96   fRunLoader(runLoader),
97   fTOFLoader(runLoader->GetLoader("TOFLoader")),
98   fTreeD(0),
99   fTreeR(0),
100   fDigits(new TClonesArray("AliTOFdigit", 4000)),
101   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
102   fNumberOfTofClusters(0)
103 {
104 //
105 // Constructor
106 //
107
108   runLoader->CdGAFile();
109   TFile *in=(TFile*)gFile;
110   in->cd();
111   fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
112
113 }
114 //______________________________________________________________________________
115
116 AliTOFClusterFinder::~AliTOFClusterFinder()
117 {
118
119   //
120   // Destructor
121   //
122
123   if (fDigits)
124     {
125       fDigits->Delete();
126       delete fDigits;
127       fDigits=0;
128     }
129   if (fRecPoints)
130     {
131       fRecPoints->Delete();
132       delete fRecPoints;
133       fRecPoints=0;
134     }
135
136   delete fTOFGeometry;
137
138 }
139 //______________________________________________________________________________
140
141 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
142 {
143   //
144   // Converts digits to recpoints for TOF
145   //
146
147   fRunLoader->GetEvent(iEvent);
148
149   fTreeD = fTOFLoader->TreeD();
150   if (fTreeD == 0x0)
151     {
152       AliFatal("AliTOFClusterFinder: Can not get TreeD");
153     }
154
155   TBranch *branch = fTreeD->GetBranch("TOF");
156   if (!branch) { 
157     AliError("can't get the branch with the TOF digits !");
158     return;
159   }
160
161   TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
162   branch->SetAddress(&digits);
163
164   ResetRecpoint();
165
166   fTreeR = fTOFLoader->TreeR();
167   if (fTreeR == 0x0)
168     {
169       fTOFLoader->MakeTree("R");
170       fTreeR = fTOFLoader->TreeR();
171     }
172
173   Int_t bufsize = 32000;
174   fTreeR->Branch("TOF", &fRecPoints, bufsize);
175
176   fTreeD->GetEvent(0);
177   Int_t nDigits = digits->GetEntriesFast();
178   AliDebug(2,Form("Number of TOF digits: %d",nDigits));
179
180   Int_t ii, jj;
181   Int_t dig[5];
182   Float_t g[3];
183   Double_t h[5];
184   Float_t tToT;
185   Double_t tTdcND;
186   for (ii=0; ii<nDigits; ii++) {
187     AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
188     dig[0]=d->GetSector();
189     dig[1]=d->GetPlate();
190     dig[2]=d->GetStrip();
191     dig[3]=d->GetPadz();
192     dig[4]=d->GetPadx();
193
194     for (jj=0; jj<3; jj++) g[jj] = 0.;
195     fTOFGeometry->GetPos(dig,g);
196
197     h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
198     h[1] = TMath::ATan2(g[1],g[0]);
199     h[2] = g[2];
200     h[3] = d->GetTdc();
201     h[4] = d->GetAdc();
202     tToT = d->GetToT();
203     tTdcND = d->GetTdcND();
204
205     AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
206     InsertCluster(tofCluster);
207
208   }
209
210   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
211
212   CalibrateRecPoint();
213   FillRecPoint();
214
215   fTreeR->Fill();
216   ResetRecpoint();
217
218   fTOFLoader = fRunLoader->GetLoader("TOFLoader");  
219   fTOFLoader->WriteRecPoints("OVERWRITE");
220
221 }
222 //______________________________________________________________________________
223
224 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
225                                            TTree *clustersTree)
226 {
227   //
228   // Converts RAW data to recpoints for TOF
229   //
230
231   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
232
233   ResetRecpoint();
234
235   Int_t bufsize = 32000;
236   clustersTree->Branch("TOF", &fRecPoints, bufsize);
237
238   Int_t ii = 0;
239   Int_t indexDDL = 0;
240
241   Int_t detectorIndex[5];
242   Float_t position[3];
243   Double_t cylindricalPosition[5];
244   Float_t tToT;
245   Double_t tTdcND;
246
247   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
248
249     rawReader->Reset();
250     AliTOFRawStream tofInput(rawReader);
251     rawReader->Select(5, indexDDL, indexDDL);
252
253     while(tofInput.Next()) {
254
255       detectorIndex[0] = tofInput.GetSector();
256       detectorIndex[1] = tofInput.GetPlate();
257       detectorIndex[2] = tofInput.GetStrip();
258       detectorIndex[3] = tofInput.GetPadZ();
259       detectorIndex[4] = tofInput.GetPadX();
260       
261       for (ii=0; ii<3; ii++) position[ii] =  0.;
262
263       fTOFGeometry->GetPos(detectorIndex, position);
264
265       cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
266       cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
267       cylindricalPosition[2] = position[2];
268       cylindricalPosition[3] = tofInput.GetTofBin();
269       cylindricalPosition[4] = tofInput.GetADCbin();
270       tToT = tofInput.GetADCbin();
271       tTdcND = -1.;
272       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
273       tofCluster->SetToT(tToT);
274       tofCluster->SetTDCND(tTdcND);
275       InsertCluster(tofCluster);
276
277     } // while loop
278
279   } // loop on DDL files
280
281   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
282
283   CalibrateRecPoint();
284   FillRecPoint();
285
286   clustersTree->Fill();
287
288   ResetRecpoint();
289
290 }
291 //______________________________________________________________________________
292
293 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
294 {
295   //
296   // Converts RAW data to recpoints for TOF
297   //
298
299   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
300
301   fRunLoader->GetEvent(iEvent);
302
303   AliDebug(2,Form(" Event number %2i ", iEvent));
304
305   fTreeR = fTOFLoader->TreeR();
306
307   if (fTreeR == 0x0){
308     fTOFLoader->MakeTree("R");
309     fTreeR = fTOFLoader->TreeR();
310   }
311
312   Int_t bufsize = 32000;
313   fTreeR->Branch("TOF", &fRecPoints, bufsize);
314
315   Int_t ii = 0;
316   Int_t indexDDL = 0;
317
318   Int_t detectorIndex[5];
319   Float_t position[3];
320   Double_t cylindricalPosition[5];
321   Float_t tToT;
322   Double_t tTdcND;
323
324   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
325
326     rawReader->Reset();
327     AliTOFRawStream tofInput(rawReader);
328     rawReader->Select(5, indexDDL, indexDDL);
329
330     while(tofInput.Next()) {
331
332       detectorIndex[0] = (Int_t)tofInput.GetSector();
333       detectorIndex[1] = (Int_t)tofInput.GetPlate();
334       detectorIndex[2] = (Int_t)tofInput.GetStrip();
335       detectorIndex[3] = (Int_t)tofInput.GetPadZ();
336       detectorIndex[4] = (Int_t)tofInput.GetPadX();
337
338       for (ii=0; ii<3; ii++) position[ii] =  0.;
339
340       fTOFGeometry->GetPos(detectorIndex, position);
341       
342       cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
343       cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
344       cylindricalPosition[2] = (Double_t)position[2];
345       cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
346       cylindricalPosition[4] = (Double_t)tofInput.GetADCbin();
347       tToT = tofInput.GetADCbin();
348       tTdcND = -1.;
349     
350       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
351       tofCluster->SetToT(tToT);
352       tofCluster->SetTDCND(tTdcND);
353       InsertCluster(tofCluster);
354
355     } // while loop
356
357   } // DDL Loop
358
359   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
360
361   CalibrateRecPoint();
362   FillRecPoint();
363
364   fTreeR->Fill();
365   ResetRecpoint();
366
367   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
368   fTOFLoader->WriteRecPoints("OVERWRITE");
369   
370 }
371 //______________________________________________________________________________
372
373 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
374 {
375   //
376   // Converts RAW data to MC digits for TOF
377   //
378   //             (temporary solution)
379   //
380
381   fRunLoader->GetEvent(iEvent);
382
383   fTreeD = fTOFLoader->TreeD();
384   if (fTreeD)
385     {
386     AliInfo("AliTOFClusterFinder: TreeD re-creation");
387     fTreeD = 0x0;
388     fTOFLoader->MakeTree("D");
389     fTreeD = fTOFLoader->TreeD();
390     }
391
392
393   fTreeR = fTOFLoader->TreeD();
394   if (fTreeD == 0x0)
395     {
396       fTOFLoader->MakeTree("D");
397       fTreeD = fTOFLoader->TreeD();
398     }
399
400   TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
401   Int_t bufsize = 32000;
402   fTreeD->Branch("TOF", &tofDigits, bufsize);
403
404
405   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
406
407   fRunLoader->GetEvent(iEvent);
408
409   AliDebug(2,Form(" Event number %2i ", iEvent));
410
411   Int_t indexDDL = 0;
412
413   Int_t detectorIndex[5];
414   Float_t digit[2];
415
416   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
417
418     rawReader->Reset();
419     AliTOFRawStream tofInput(rawReader);
420     rawReader->Select(5, indexDDL, indexDDL);
421
422     while(tofInput.Next()) {
423
424       detectorIndex[0] = tofInput.GetSector();
425       detectorIndex[1] = tofInput.GetPlate();
426       detectorIndex[2] = tofInput.GetStrip();
427       detectorIndex[3] = tofInput.GetPadX();
428       detectorIndex[4] = tofInput.GetPadZ();
429       
430       digit[0] = (Float_t)tofInput.GetTofBin();
431       digit[1] = (Float_t)tofInput.GetADCbin();
432
433       Int_t tracknum[3]={-1,-1,-1};
434
435       TClonesArray &aDigits = *tofDigits;
436       Int_t last=tofDigits->GetEntriesFast();
437       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
438
439     } // while loop
440
441   } // DDL Loop
442
443   fTreeD->Fill();
444
445   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
446   fTOFLoader->WriteDigits("OVERWRITE");
447   
448 }
449 //______________________________________________________________________________
450
451 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
452   //---------------------------------------------------------------------------//
453   // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
454   //---------------------------------------------------------------------------//
455   if (fNumberOfTofClusters==kTofMaxCluster) {
456     AliError("Too many clusters !");
457     return 1;
458   }
459
460   if (fNumberOfTofClusters==0) {
461     fTofClusters[fNumberOfTofClusters++] = tofCluster;
462     return 0;
463   }
464
465   Int_t ii = FindClusterIndex(tofCluster->GetZ());
466   memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
467   fTofClusters[ii] = tofCluster;
468   fNumberOfTofClusters++;
469   
470   return 0;
471
472 }
473 //_________________________________________________________________________
474
475 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
476   //--------------------------------------------------------------------
477   // This function returns the index of the nearest cluster 
478   //--------------------------------------------------------------------
479   if (fNumberOfTofClusters==0) return 0;
480   if (z <= fTofClusters[0]->GetZ()) return 0;
481   if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
482   Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
483   for (; b<e; m=(b+e)/2) {
484     if (z > fTofClusters[m]->GetZ()) b=m+1;
485     else e=m;
486   }
487
488   return m;
489
490 }
491 //_________________________________________________________________________
492
493 void AliTOFClusterFinder::FillRecPoint()
494 {
495   //
496   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
497   // in Z) in the global TClonesArray of AliTOFcluster,
498   // i.e. fRecPoints.
499   //
500
501   Int_t ii, jj;
502
503   Int_t detectorIndex[5];
504   Double_t cylindricalPosition[5];
505   Int_t trackLabels[3];
506   Int_t digitIndex = -1;
507   Float_t tToT;
508   Double_t tTdcND;
509
510   TClonesArray &lRecPoints = *fRecPoints;
511   
512   for (ii=0; ii<fNumberOfTofClusters; ii++) {
513
514     digitIndex = fTofClusters[ii]->GetIndex();
515     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
516     for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
517     cylindricalPosition[0] = fTofClusters[ii]->GetR();
518     cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
519     cylindricalPosition[2] = fTofClusters[ii]->GetZ();
520     cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
521     cylindricalPosition[4] = fTofClusters[ii]->GetADC();
522     tToT = fTofClusters[ii]->GetToT();
523     tTdcND = fTofClusters[ii]->GetTDCND();
524
525     new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND);
526
527     //AliInfo(Form("%3i  %3i  %f %f %f %f %f  %2i %2i %2i %1i %2i",ii,digitIndex, cylindricalPosition[2],cylindricalPosition[0],cylindricalPosition[1],cylindricalPosition[3],cylindricalPosition[4],detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
528
529   } // loop on clusters
530
531 }
532
533 //_________________________________________________________________________
534 void AliTOFClusterFinder::CalibrateRecPoint()
535 {
536   //
537   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
538   // in Z) in the global TClonesArray of AliTOFcluster,
539   // i.e. fRecPoints.
540   //
541
542   Int_t ii, jj;
543
544   Int_t detectorIndex[5];
545   Int_t digitIndex = -1;
546   Float_t tToT;
547   Float_t tdcCorr;
548   AliInfo(" Calibrating TOF Clusters: ")
549   AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
550   calib->ReadParFromCDB("TOF/Calib",0);
551   AliTOFCal *calTOFArray = calib->GetTOFCalArray();  
552
553   for (ii=0; ii<fNumberOfTofClusters; ii++) {
554     digitIndex = fTofClusters[ii]->GetIndex();
555     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
556
557     Int_t index = calib->GetIndex(detectorIndex);
558      
559     AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
560     Float_t par[6];
561     for (Int_t j = 0; j<6; j++){
562       par[j]=calChannel->GetSlewPar(j);
563     }
564     tToT = fTofClusters[ii]->GetToT();
565     Float_t timeCorr=par[0]+par[1]*tToT+par[2]*tToT*tToT+par[3]*tToT*tToT*tToT+par[4]*tToT*tToT*tToT*tToT+par[5]*tToT*tToT*tToT*tToT*tToT;
566     tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
567     tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
568     fTofClusters[ii]->SetTDC(tdcCorr);
569
570   } // loop on clusters
571
572   delete calib;
573 }
574 //______________________________________________________________________________
575
576 void AliTOFClusterFinder::ResetRecpoint()
577 {
578   //
579   // Clear the list of reconstructed points
580   //
581
582   fNumberOfTofClusters = 0;
583   if (fRecPoints) fRecPoints->Clear();
584
585 }
586 //______________________________________________________________________________
587
588 void AliTOFClusterFinder::Load()
589 {
590   //
591   // Load TOF.Digits.root and TOF.RecPoints.root files
592   //
593
594   fTOFLoader->LoadDigits("READ");
595   fTOFLoader->LoadRecPoints("recreate");
596
597 }
598 //______________________________________________________________________________
599
600 void AliTOFClusterFinder::LoadClusters()
601 {
602   //
603   // Load TOF.RecPoints.root file
604   //
605
606   fTOFLoader->LoadRecPoints("recreate");
607
608 }
609 //______________________________________________________________________________
610
611 void AliTOFClusterFinder::UnLoad()
612 {
613   //
614   // Unload TOF.Digits.root and TOF.RecPoints.root files
615   //
616
617   fTOFLoader->UnloadDigits();
618   fTOFLoader->UnloadRecPoints();
619
620 }
621 //______________________________________________________________________________
622
623 void AliTOFClusterFinder::UnLoadClusters()
624 {
625   //
626   // Unload TOF.RecPoints.root file
627   //
628
629   fTOFLoader->UnloadRecPoints();
630
631 }
632 //______________________________________________________________________________