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