]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFClusterFinder.cxx
Coding conventions
[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 "TClonesArray.h"
54 #include "TFile.h"
55 #include "TTree.h"
56
57 #include "AliLoader.h"
58 #include "AliLog.h"
59 #include "AliRawReader.h"
60 #include "AliRunLoader.h"
61
62 #include "AliTOFcalib.h"
63 #include "AliTOFCal.h"
64 #include "AliTOFChannel.h"
65 #include "AliTOFClusterFinder.h"
66 #include "AliTOFcluster.h"
67 #include "AliTOFdigit.h"
68 #include "AliTOFGeometryV5.h"
69 #include "AliTOFGeometry.h"
70 #include "AliTOFRawStream.h"
71
72 extern TFile *gFile;
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(const AliTOFClusterFinder &source)
117   :TObject()
118 {
119   // copy constructor
120   this->fDigits=source.fDigits;
121   this->fRecPoints=source.fRecPoints;
122   this->fTOFGeometry=source.fTOFGeometry;
123
124 }
125
126 //------------------------------------------------------------------------
127   AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
128 {
129   // ass. op.
130   this->fDigits=source.fDigits;
131   this->fRecPoints=source.fRecPoints;
132   this->fTOFGeometry=source.fTOFGeometry;
133   return *this;
134
135 }
136 //______________________________________________________________________________
137
138 AliTOFClusterFinder::~AliTOFClusterFinder()
139 {
140
141   //
142   // Destructor
143   //
144
145   if (fDigits)
146     {
147       fDigits->Delete();
148       delete fDigits;
149       fDigits=0;
150     }
151   if (fRecPoints)
152     {
153       fRecPoints->Delete();
154       delete fRecPoints;
155       fRecPoints=0;
156     }
157
158   delete fTOFGeometry;
159
160 }
161 //______________________________________________________________________________
162
163 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
164 {
165   //
166   // Converts digits to recpoints for TOF
167   //
168
169   fRunLoader->GetEvent(iEvent);
170
171   fTreeD = fTOFLoader->TreeD();
172   if (fTreeD == 0x0)
173     {
174       AliFatal("AliTOFClusterFinder: Can not get TreeD");
175     }
176
177   TBranch *branch = fTreeD->GetBranch("TOF");
178   if (!branch) { 
179     AliError("can't get the branch with the TOF digits !");
180     return;
181   }
182
183   TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
184   branch->SetAddress(&digits);
185
186   ResetRecpoint();
187
188   fTreeR = fTOFLoader->TreeR();
189   if (fTreeR == 0x0)
190     {
191       fTOFLoader->MakeTree("R");
192       fTreeR = fTOFLoader->TreeR();
193     }
194
195   Int_t bufsize = 32000;
196   fTreeR->Branch("TOF", &fRecPoints, bufsize);
197
198   fTreeD->GetEvent(0);
199   Int_t nDigits = digits->GetEntriesFast();
200   AliDebug(2,Form("Number of TOF digits: %d",nDigits));
201
202   Int_t ii, jj;
203   Int_t dig[5];
204   Float_t g[3];
205   Double_t h[5];
206   Float_t tToT;
207   Double_t tTdcND;
208   for (ii=0; ii<nDigits; ii++) {
209     AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
210     dig[0]=d->GetSector();
211     dig[1]=d->GetPlate();
212     dig[2]=d->GetStrip();
213     dig[3]=d->GetPadz();
214     dig[4]=d->GetPadx();
215
216     for (jj=0; jj<3; jj++) g[jj] = 0.;
217     fTOFGeometry->GetPos(dig,g);
218
219     h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
220     h[1] = TMath::ATan2(g[1],g[0]);
221     h[2] = g[2];
222     h[3] = d->GetTdc();
223     h[4] = d->GetAdc();
224     tToT = d->GetToT();
225     tTdcND = d->GetTdcND();
226
227     AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
228     InsertCluster(tofCluster);
229
230   }
231
232   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
233
234   CalibrateRecPoint();
235   FillRecPoint();
236
237   fTreeR->Fill();
238   ResetRecpoint();
239
240   fTOFLoader = fRunLoader->GetLoader("TOFLoader");  
241   fTOFLoader->WriteRecPoints("OVERWRITE");
242
243 }
244 //______________________________________________________________________________
245
246 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
247                                            TTree *clustersTree)
248 {
249   //
250   // Converts RAW data to recpoints for TOF
251   //
252
253   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
254
255   ResetRecpoint();
256
257   Int_t bufsize = 32000;
258   clustersTree->Branch("TOF", &fRecPoints, bufsize);
259
260   Int_t ii = 0;
261   Int_t indexDDL = 0;
262
263   Int_t detectorIndex[5];
264   Float_t position[3];
265   Double_t cylindricalPosition[5];
266   Float_t tToT;
267   Double_t tTdcND;
268
269   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
270
271     rawReader->Reset();
272     AliTOFRawStream tofInput(rawReader);
273     rawReader->Select(5, indexDDL, indexDDL);
274
275     while(tofInput.Next()) {
276
277       detectorIndex[0] = tofInput.GetSector();
278       detectorIndex[1] = tofInput.GetPlate();
279       detectorIndex[2] = tofInput.GetStrip();
280       detectorIndex[3] = tofInput.GetPadZ();
281       detectorIndex[4] = tofInput.GetPadX();
282       
283       for (ii=0; ii<3; ii++) position[ii] =  0.;
284
285       fTOFGeometry->GetPos(detectorIndex, position);
286
287       cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
288       cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
289       cylindricalPosition[2] = position[2];
290       cylindricalPosition[3] = tofInput.GetTofBin();
291       cylindricalPosition[4] = tofInput.GetADCbin();
292       tToT = tofInput.GetADCbin();
293       tTdcND = -1.;
294       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
295       tofCluster->SetToT(tToT);
296       tofCluster->SetTDCND(tTdcND);
297       InsertCluster(tofCluster);
298
299     } // while loop
300
301   } // loop on DDL files
302
303   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
304
305   CalibrateRecPoint();
306   FillRecPoint();
307
308   clustersTree->Fill();
309
310   ResetRecpoint();
311
312 }
313 //______________________________________________________________________________
314
315 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
316 {
317   //
318   // Converts RAW data to recpoints for TOF
319   //
320
321   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
322
323   fRunLoader->GetEvent(iEvent);
324
325   AliDebug(2,Form(" Event number %2i ", iEvent));
326
327   fTreeR = fTOFLoader->TreeR();
328
329   if (fTreeR == 0x0){
330     fTOFLoader->MakeTree("R");
331     fTreeR = fTOFLoader->TreeR();
332   }
333
334   Int_t bufsize = 32000;
335   fTreeR->Branch("TOF", &fRecPoints, bufsize);
336
337   Int_t ii = 0;
338   Int_t indexDDL = 0;
339
340   Int_t detectorIndex[5];
341   Float_t position[3];
342   Double_t cylindricalPosition[5];
343   Float_t tToT;
344   Double_t tTdcND;
345
346   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
347
348     rawReader->Reset();
349     AliTOFRawStream tofInput(rawReader);
350     rawReader->Select(5, indexDDL, indexDDL);
351
352     while(tofInput.Next()) {
353
354       detectorIndex[0] = (Int_t)tofInput.GetSector();
355       detectorIndex[1] = (Int_t)tofInput.GetPlate();
356       detectorIndex[2] = (Int_t)tofInput.GetStrip();
357       detectorIndex[3] = (Int_t)tofInput.GetPadZ();
358       detectorIndex[4] = (Int_t)tofInput.GetPadX();
359
360       for (ii=0; ii<3; ii++) position[ii] =  0.;
361
362       fTOFGeometry->GetPos(detectorIndex, position);
363       
364       cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
365       cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
366       cylindricalPosition[2] = (Double_t)position[2];
367       cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
368       cylindricalPosition[4] = (Double_t)tofInput.GetADCbin();
369       tToT = tofInput.GetADCbin();
370       tTdcND = -1.;
371     
372       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
373       tofCluster->SetToT(tToT);
374       tofCluster->SetTDCND(tTdcND);
375       InsertCluster(tofCluster);
376
377     } // while loop
378
379   } // DDL Loop
380
381   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
382
383   CalibrateRecPoint();
384   FillRecPoint();
385
386   fTreeR->Fill();
387   ResetRecpoint();
388
389   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
390   fTOFLoader->WriteRecPoints("OVERWRITE");
391   
392 }
393 //______________________________________________________________________________
394
395 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
396 {
397   //
398   // Converts RAW data to MC digits for TOF
399   //
400   //             (temporary solution)
401   //
402
403   fRunLoader->GetEvent(iEvent);
404
405   fTreeD = fTOFLoader->TreeD();
406   if (fTreeD)
407     {
408     AliInfo("AliTOFClusterFinder: TreeD re-creation");
409     fTreeD = 0x0;
410     fTOFLoader->MakeTree("D");
411     fTreeD = fTOFLoader->TreeD();
412     }
413
414
415   fTreeR = fTOFLoader->TreeD();
416   if (fTreeD == 0x0)
417     {
418       fTOFLoader->MakeTree("D");
419       fTreeD = fTOFLoader->TreeD();
420     }
421
422   TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
423   Int_t bufsize = 32000;
424   fTreeD->Branch("TOF", &tofDigits, bufsize);
425
426
427   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
428
429   fRunLoader->GetEvent(iEvent);
430
431   AliDebug(2,Form(" Event number %2i ", iEvent));
432
433   Int_t indexDDL = 0;
434
435   Int_t detectorIndex[5];
436   Float_t digit[2];
437
438   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
439
440     rawReader->Reset();
441     AliTOFRawStream tofInput(rawReader);
442     rawReader->Select(5, indexDDL, indexDDL);
443
444     while(tofInput.Next()) {
445
446       detectorIndex[0] = tofInput.GetSector();
447       detectorIndex[1] = tofInput.GetPlate();
448       detectorIndex[2] = tofInput.GetStrip();
449       detectorIndex[3] = tofInput.GetPadX();
450       detectorIndex[4] = tofInput.GetPadZ();
451       
452       digit[0] = (Float_t)tofInput.GetTofBin();
453       digit[1] = (Float_t)tofInput.GetADCbin();
454
455       Int_t tracknum[3]={-1,-1,-1};
456
457       TClonesArray &aDigits = *tofDigits;
458       Int_t last=tofDigits->GetEntriesFast();
459       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
460
461     } // while loop
462
463   } // DDL Loop
464
465   fTreeD->Fill();
466
467   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
468   fTOFLoader->WriteDigits("OVERWRITE");
469   
470 }
471 //______________________________________________________________________________
472
473 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
474   //---------------------------------------------------------------------------//
475   // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
476   //---------------------------------------------------------------------------//
477   if (fNumberOfTofClusters==kTofMaxCluster) {
478     AliError("Too many clusters !");
479     return 1;
480   }
481
482   if (fNumberOfTofClusters==0) {
483     fTofClusters[fNumberOfTofClusters++] = tofCluster;
484     return 0;
485   }
486
487   Int_t ii = FindClusterIndex(tofCluster->GetZ());
488   memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
489   fTofClusters[ii] = tofCluster;
490   fNumberOfTofClusters++;
491   
492   return 0;
493
494 }
495 //_________________________________________________________________________
496
497 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
498   //--------------------------------------------------------------------
499   // This function returns the index of the nearest cluster 
500   //--------------------------------------------------------------------
501   if (fNumberOfTofClusters==0) return 0;
502   if (z <= fTofClusters[0]->GetZ()) return 0;
503   if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
504   Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
505   for (; b<e; m=(b+e)/2) {
506     if (z > fTofClusters[m]->GetZ()) b=m+1;
507     else e=m;
508   }
509
510   return m;
511
512 }
513 //_________________________________________________________________________
514
515 void AliTOFClusterFinder::FillRecPoint()
516 {
517   //
518   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
519   // in Z) in the global TClonesArray of AliTOFcluster,
520   // i.e. fRecPoints.
521   //
522
523   Int_t ii, jj;
524
525   Int_t detectorIndex[5];
526   Double_t cylindricalPosition[5];
527   Int_t trackLabels[3];
528   Int_t digitIndex = -1;
529   Float_t tToT;
530   Double_t tTdcND;
531
532   TClonesArray &lRecPoints = *fRecPoints;
533   
534   for (ii=0; ii<fNumberOfTofClusters; ii++) {
535
536     digitIndex = fTofClusters[ii]->GetIndex();
537     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
538     for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
539     cylindricalPosition[0] = fTofClusters[ii]->GetR();
540     cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
541     cylindricalPosition[2] = fTofClusters[ii]->GetZ();
542     cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
543     cylindricalPosition[4] = fTofClusters[ii]->GetADC();
544     tToT = fTofClusters[ii]->GetToT();
545     tTdcND = fTofClusters[ii]->GetTDCND();
546
547     new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND);
548
549     //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]));
550
551   } // loop on clusters
552
553 }
554
555 //_________________________________________________________________________
556 void AliTOFClusterFinder::CalibrateRecPoint()
557 {
558   //
559   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
560   // in Z) in the global TClonesArray of AliTOFcluster,
561   // i.e. fRecPoints.
562   //
563
564   Int_t ii, jj;
565
566   Int_t detectorIndex[5];
567   Int_t digitIndex = -1;
568   Float_t tToT;
569   Float_t tdcCorr;
570   AliInfo(" Calibrating TOF Clusters: ")
571   AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
572   // calib->ReadParFromCDB("TOF/Calib",0); // original
573   calib->ReadParFromCDB("TOF/Calib",-1); // Use AliCDBManager's run number
574   AliTOFCal *calTOFArray = calib->GetTOFCalArray();  
575
576   for (ii=0; ii<fNumberOfTofClusters; ii++) {
577     digitIndex = fTofClusters[ii]->GetIndex();
578     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
579
580     Int_t index = calib->GetIndex(detectorIndex);
581      
582     AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
583     Float_t par[6];
584     for (Int_t j = 0; j<6; j++){
585       par[j]=calChannel->GetSlewPar(j);
586     }
587     tToT = fTofClusters[ii]->GetToT();
588     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;
589     tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
590     tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
591     fTofClusters[ii]->SetTDC(tdcCorr);
592
593   } // loop on clusters
594
595   delete calib;
596 }
597 //______________________________________________________________________________
598
599 void AliTOFClusterFinder::ResetRecpoint()
600 {
601   //
602   // Clear the list of reconstructed points
603   //
604
605   fNumberOfTofClusters = 0;
606   if (fRecPoints) fRecPoints->Clear();
607
608 }
609 //______________________________________________________________________________
610
611 void AliTOFClusterFinder::Load()
612 {
613   //
614   // Load TOF.Digits.root and TOF.RecPoints.root files
615   //
616
617   fTOFLoader->LoadDigits("READ");
618   fTOFLoader->LoadRecPoints("recreate");
619
620 }
621 //______________________________________________________________________________
622
623 void AliTOFClusterFinder::LoadClusters()
624 {
625   //
626   // Load TOF.RecPoints.root file
627   //
628
629   fTOFLoader->LoadRecPoints("recreate");
630
631 }
632 //______________________________________________________________________________
633
634 void AliTOFClusterFinder::UnLoad()
635 {
636   //
637   // Unload TOF.Digits.root and TOF.RecPoints.root files
638   //
639
640   fTOFLoader->UnloadDigits();
641   fTOFLoader->UnloadRecPoints();
642
643 }
644 //______________________________________________________________________________
645
646 void AliTOFClusterFinder::UnLoadClusters()
647 {
648   //
649   // Unload TOF.RecPoints.root file
650   //
651
652   fTOFLoader->UnloadRecPoints();
653
654 }
655 //______________________________________________________________________________