]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFClusterFinder.cxx
New TOF geometry description (V5) -G. Cara Romeo and A. De Caro
[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 <Riostream.h>
54 #include <TTree.h>
55 #include <TObjArray.h>
56 #include <TClonesArray.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
69 #include "AliTOFClusterFinder.h"
70
71 ClassImp(AliTOFClusterFinder)
72
73 AliTOFClusterFinder::AliTOFClusterFinder():
74   fRunLoader(0),
75   fTOFLoader(0),
76   fTreeD(0),
77   fTreeR(0),
78   fDigits(new TClonesArray("AliTOFdigit", 4000)),
79   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
80   fNumberOfTofClusters(0)
81 {
82 //
83 // Constructor
84 //
85
86   fTOFGeometry = new AliTOFGeometry();
87
88 }
89 //______________________________________________________________________________
90
91 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
92   fRunLoader(runLoader),
93   fTOFLoader(runLoader->GetLoader("TOFLoader")),
94   fTreeD(0),
95   fTreeR(0),
96   fDigits(new TClonesArray("AliTOFdigit", 4000)),
97   fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
98   fNumberOfTofClusters(0)
99 {
100 //
101 // Constructor
102 //
103
104   runLoader->CdGAFile();
105   TFile *in=(TFile*)gFile;
106   in->cd();
107   fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
108
109 }
110 //______________________________________________________________________________
111
112 AliTOFClusterFinder::~AliTOFClusterFinder()
113 {
114
115   //
116   // Destructor
117   //
118
119   if (fDigits)
120     {
121       fDigits->Delete();
122       delete fDigits;
123       fDigits=0;
124     }
125   if (fRecPoints)
126     {
127       fRecPoints->Delete();
128       delete fRecPoints;
129       fRecPoints=0;
130     }
131
132   delete fTOFGeometry;
133
134 }
135 //______________________________________________________________________________
136
137 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
138 {
139   //
140   // Converts digits to recpoints for TOF
141   //
142
143   fRunLoader->GetEvent(iEvent);
144
145   fTreeD = fTOFLoader->TreeD();
146   if (fTreeD == 0x0)
147     {
148       AliFatal("AliTOFClusterFinder: Can not get TreeD");
149     }
150
151   TBranch *branch = fTreeD->GetBranch("TOF");
152   if (!branch) { 
153     AliError("can't get the branch with the TOF digits !");
154     return;
155   }
156
157   TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
158   branch->SetAddress(&digits);
159
160   ResetRecpoint();
161
162   fTreeR = fTOFLoader->TreeR();
163   if (fTreeR == 0x0)
164     {
165       fTOFLoader->MakeTree("R");
166       fTreeR = fTOFLoader->TreeR();
167     }
168
169   Int_t bufsize = 32000;
170   fTreeR->Branch("TOF", &fRecPoints, bufsize);
171
172   fTreeD->GetEvent(0);
173   Int_t nDigits = digits->GetEntriesFast();
174   AliDebug(2,Form("Number of TOF digits: %d",nDigits));
175
176   Int_t ii, jj;
177   Int_t dig[5];
178   Float_t g[3];
179   Double_t h[5];
180   for (ii=0; ii<nDigits; ii++) {
181     AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
182     dig[0]=d->GetSector();
183     dig[1]=d->GetPlate();
184     dig[2]=d->GetStrip();
185     dig[3]=d->GetPadz();
186     dig[4]=d->GetPadx();
187
188     for (jj=0; jj<3; jj++) g[jj] = 0.;
189     fTOFGeometry->GetPos(dig,g);
190
191     h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
192     h[1] = TMath::ATan2(g[1],g[0]);
193     h[2] = g[2];
194     h[3] = d->GetTdc();
195     h[4] = d->GetAdc();
196
197     AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii);
198     InsertCluster(tofCluster);
199
200   }
201
202   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
203
204   FillRecPoint();
205
206   fTreeR->Fill();
207   ResetRecpoint();
208
209   fTOFLoader = fRunLoader->GetLoader("TOFLoader");  
210   fTOFLoader->WriteRecPoints("OVERWRITE");
211
212 }
213 //______________________________________________________________________________
214
215 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
216                                            TTree *clustersTree)
217 {
218   //
219   // Converts RAW data to recpoints for TOF
220   //
221
222   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
223
224   ResetRecpoint();
225
226   Int_t bufsize = 32000;
227   clustersTree->Branch("TOF", &fRecPoints, bufsize);
228
229   Int_t ii = 0;
230   Int_t indexDDL = 0;
231
232   Int_t detectorIndex[5];
233   Float_t position[3];
234   Double_t cylindricalPosition[5];
235
236   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
237
238     rawReader->Reset();
239     AliTOFRawStream tofInput(rawReader);
240     rawReader->Select(5, indexDDL, indexDDL);
241
242     while(tofInput.Next()) {
243
244       detectorIndex[0] = tofInput.GetSector();
245       detectorIndex[1] = tofInput.GetPlate();
246       detectorIndex[2] = tofInput.GetStrip();
247       detectorIndex[3] = tofInput.GetPadZ();
248       detectorIndex[4] = tofInput.GetPadX();
249       
250       for (ii=0; ii<3; ii++) position[ii] =  0.;
251
252       fTOFGeometry->GetPos(detectorIndex, position);
253
254       cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
255       cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
256       cylindricalPosition[2] = position[2];
257       cylindricalPosition[3] = tofInput.GetTofBin();
258       cylindricalPosition[4] = tofInput.GetADCbin();
259
260       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
261       InsertCluster(tofCluster);
262
263     } // while loop
264
265   } // loop on DDL files
266
267   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
268
269   FillRecPoint();
270
271   clustersTree->Fill();
272   ResetRecpoint();
273
274 }
275 //______________________________________________________________________________
276
277 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
278 {
279   //
280   // Converts RAW data to recpoints for TOF
281   //
282
283   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
284
285   fRunLoader->GetEvent(iEvent);
286
287   AliDebug(2,Form(" Event number %2i ", iEvent));
288
289   fTreeR = fTOFLoader->TreeR();
290
291   if (fTreeR == 0x0){
292     fTOFLoader->MakeTree("R");
293     fTreeR = fTOFLoader->TreeR();
294   }
295
296   Int_t bufsize = 32000;
297   fTreeR->Branch("TOF", &fRecPoints, bufsize);
298
299   Int_t ii = 0;
300   Int_t indexDDL = 0;
301
302   Int_t detectorIndex[5];
303   Float_t position[3];
304   Double_t cylindricalPosition[5];
305
306   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
307
308     rawReader->Reset();
309     AliTOFRawStream tofInput(rawReader);
310     rawReader->Select(5, indexDDL, indexDDL);
311
312     while(tofInput.Next()) {
313
314       detectorIndex[0] = (Int_t)tofInput.GetSector();
315       detectorIndex[1] = (Int_t)tofInput.GetPlate();
316       detectorIndex[2] = (Int_t)tofInput.GetStrip();
317       detectorIndex[3] = (Int_t)tofInput.GetPadZ();
318       detectorIndex[4] = (Int_t)tofInput.GetPadX();
319
320       for (ii=0; ii<3; ii++) position[ii] =  0.;
321
322       fTOFGeometry->GetPos(detectorIndex, position);
323       
324       cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
325       cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
326       cylindricalPosition[2] = (Double_t)position[2];
327       cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
328       cylindricalPosition[4] = (Double_t)tofInput.GetADCbin();
329
330       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
331       InsertCluster(tofCluster);
332
333     } // while loop
334
335   } // DDL Loop
336
337   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
338
339   FillRecPoint();
340
341   fTreeR->Fill();
342   ResetRecpoint();
343
344   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
345   fTOFLoader->WriteRecPoints("OVERWRITE");
346   
347 }
348 //______________________________________________________________________________
349
350 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
351 {
352   //
353   // Converts RAW data to MC digits for TOF
354   //
355   //             (temporary solution)
356   //
357
358   fRunLoader->GetEvent(iEvent);
359
360   fTreeD = fTOFLoader->TreeD();
361   if (fTreeD)
362     {
363     AliInfo("AliTOFClusterFinder: TreeD re-creation");
364     fTreeD = 0x0;
365     fTOFLoader->MakeTree("D");
366     fTreeD = fTOFLoader->TreeD();
367     }
368
369
370   fTreeR = fTOFLoader->TreeD();
371   if (fTreeD == 0x0)
372     {
373       fTOFLoader->MakeTree("D");
374       fTreeD = fTOFLoader->TreeD();
375     }
376
377   TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
378   Int_t bufsize = 32000;
379   fTreeD->Branch("TOF", &tofDigits, bufsize);
380
381
382   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
383
384   fRunLoader->GetEvent(iEvent);
385
386   AliDebug(2,Form(" Event number %2i ", iEvent));
387
388   Int_t indexDDL = 0;
389
390   Int_t detectorIndex[5];
391   Float_t digit[2];
392
393   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
394
395     rawReader->Reset();
396     AliTOFRawStream tofInput(rawReader);
397     rawReader->Select(5, indexDDL, indexDDL);
398
399     while(tofInput.Next()) {
400
401       detectorIndex[0] = tofInput.GetSector();
402       detectorIndex[1] = tofInput.GetPlate();
403       detectorIndex[2] = tofInput.GetStrip();
404       detectorIndex[3] = tofInput.GetPadX();
405       detectorIndex[4] = tofInput.GetPadZ();
406       
407       digit[0] = (Float_t)tofInput.GetTofBin();
408       digit[1] = (Float_t)tofInput.GetADCbin();
409
410       Int_t tracknum[3]={-1,-1,-1};
411
412       TClonesArray &aDigits = *tofDigits;
413       Int_t last=tofDigits->GetEntriesFast();
414       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
415
416     } // while loop
417
418   } // DDL Loop
419
420   fTreeD->Fill();
421
422   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
423   fTOFLoader->WriteDigits("OVERWRITE");
424   
425 }
426 //______________________________________________________________________________
427
428 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
429   //---------------------------------------------------------------------------//
430   // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
431   //---------------------------------------------------------------------------//
432   if (fNumberOfTofClusters==kTofMaxCluster) {
433     AliError("Too many clusters !");
434     return 1;
435   }
436
437   if (fNumberOfTofClusters==0) {
438     fTofClusters[fNumberOfTofClusters++] = tofCluster;
439     return 0;
440   }
441
442   Int_t ii = FindClusterIndex(tofCluster->GetZ());
443   memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
444   fTofClusters[ii] = tofCluster;
445   fNumberOfTofClusters++;
446   
447   return 0;
448
449 }
450 //_________________________________________________________________________
451
452 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
453   //--------------------------------------------------------------------
454   // This function returns the index of the nearest cluster 
455   //--------------------------------------------------------------------
456   if (fNumberOfTofClusters==0) return 0;
457   if (z <= fTofClusters[0]->GetZ()) return 0;
458   if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
459   Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
460   for (; b<e; m=(b+e)/2) {
461     if (z > fTofClusters[m]->GetZ()) b=m+1;
462     else e=m;
463   }
464
465   return m;
466
467 }
468 //_________________________________________________________________________
469
470 void AliTOFClusterFinder::FillRecPoint()
471 {
472   //
473   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
474   // in Z) in the global TClonesArray of AliTOFcluster,
475   // i.e. fRecPoints.
476   //
477
478   Int_t ii, jj;
479
480   Int_t detectorIndex[5];
481   Double_t cylindricalPosition[5];
482   Int_t trackLabels[3];
483   Int_t digitIndex = -1;
484
485   TClonesArray &lRecPoints = *fRecPoints;
486   
487   for (ii=0; ii<fNumberOfTofClusters; ii++) {
488
489     digitIndex = fTofClusters[ii]->GetIndex();
490     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
491     for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
492     cylindricalPosition[0] = fTofClusters[ii]->GetR();
493     cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
494     cylindricalPosition[2] = fTofClusters[ii]->GetZ();
495     cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
496     cylindricalPosition[4] = fTofClusters[ii]->GetADC();
497
498     new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex);
499
500     //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]));
501
502   } // loop on clusters
503
504 }
505 //______________________________________________________________________________
506
507 void AliTOFClusterFinder::ResetRecpoint()
508 {
509   //
510   // Clear the list of reconstructed points
511   //
512
513   fNumberOfTofClusters = 0;
514   if (fRecPoints) fRecPoints->Clear();
515
516 }
517 //______________________________________________________________________________
518
519 void AliTOFClusterFinder::Load()
520 {
521   //
522   // Load TOF.Digits.root and TOF.RecPoints.root files
523   //
524
525   fTOFLoader->LoadDigits("READ");
526   fTOFLoader->LoadRecPoints("recreate");
527
528 }
529 //______________________________________________________________________________
530
531 void AliTOFClusterFinder::LoadClusters()
532 {
533   //
534   // Load TOF.RecPoints.root file
535   //
536
537   fTOFLoader->LoadRecPoints("recreate");
538
539 }
540 //______________________________________________________________________________
541
542 void AliTOFClusterFinder::UnLoad()
543 {
544   //
545   // Unload TOF.Digits.root and TOF.RecPoints.root files
546   //
547
548   fTOFLoader->UnloadDigits();
549   fTOFLoader->UnloadRecPoints();
550
551 }
552 //______________________________________________________________________________
553
554 void AliTOFClusterFinder::UnLoadClusters()
555 {
556   //
557   // Unload TOF.RecPoints.root file
558   //
559
560   fTOFLoader->UnloadRecPoints();
561
562 }
563 //______________________________________________________________________________