]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFClusterFinder.cxx
Including TFile.h
[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   for (ii=0; ii<nDigits; ii++) {
182     AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
183     dig[0]=d->GetSector();
184     dig[1]=d->GetPlate();
185     dig[2]=d->GetStrip();
186     dig[3]=d->GetPadz();
187     dig[4]=d->GetPadx();
188
189     for (jj=0; jj<3; jj++) g[jj] = 0.;
190     fTOFGeometry->GetPos(dig,g);
191
192     h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
193     h[1] = TMath::ATan2(g[1],g[0]);
194     h[2] = g[2];
195     h[3] = d->GetTdc();
196     h[4] = d->GetAdc();
197
198     AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii);
199     InsertCluster(tofCluster);
200
201   }
202
203   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
204
205   FillRecPoint();
206
207   fTreeR->Fill();
208   ResetRecpoint();
209
210   fTOFLoader = fRunLoader->GetLoader("TOFLoader");  
211   fTOFLoader->WriteRecPoints("OVERWRITE");
212
213 }
214 //______________________________________________________________________________
215
216 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
217                                            TTree *clustersTree)
218 {
219   //
220   // Converts RAW data to recpoints for TOF
221   //
222
223   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
224
225   ResetRecpoint();
226
227   Int_t bufsize = 32000;
228   clustersTree->Branch("TOF", &fRecPoints, bufsize);
229
230   Int_t ii = 0;
231   Int_t indexDDL = 0;
232
233   Int_t detectorIndex[5];
234   Float_t position[3];
235   Double_t cylindricalPosition[5];
236
237   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
238
239     rawReader->Reset();
240     AliTOFRawStream tofInput(rawReader);
241     rawReader->Select(5, indexDDL, indexDDL);
242
243     while(tofInput.Next()) {
244
245       detectorIndex[0] = tofInput.GetSector();
246       detectorIndex[1] = tofInput.GetPlate();
247       detectorIndex[2] = tofInput.GetStrip();
248       detectorIndex[3] = tofInput.GetPadZ();
249       detectorIndex[4] = tofInput.GetPadX();
250       
251       for (ii=0; ii<3; ii++) position[ii] =  0.;
252
253       fTOFGeometry->GetPos(detectorIndex, position);
254
255       cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
256       cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
257       cylindricalPosition[2] = position[2];
258       cylindricalPosition[3] = tofInput.GetTofBin();
259       cylindricalPosition[4] = tofInput.GetADCbin();
260
261       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
262       InsertCluster(tofCluster);
263
264     } // while loop
265
266   } // loop on DDL files
267
268   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
269
270   FillRecPoint();
271
272   clustersTree->Fill();
273   ResetRecpoint();
274
275 }
276 //______________________________________________________________________________
277
278 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
279 {
280   //
281   // Converts RAW data to recpoints for TOF
282   //
283
284   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
285
286   fRunLoader->GetEvent(iEvent);
287
288   AliDebug(2,Form(" Event number %2i ", iEvent));
289
290   fTreeR = fTOFLoader->TreeR();
291
292   if (fTreeR == 0x0){
293     fTOFLoader->MakeTree("R");
294     fTreeR = fTOFLoader->TreeR();
295   }
296
297   Int_t bufsize = 32000;
298   fTreeR->Branch("TOF", &fRecPoints, bufsize);
299
300   Int_t ii = 0;
301   Int_t indexDDL = 0;
302
303   Int_t detectorIndex[5];
304   Float_t position[3];
305   Double_t cylindricalPosition[5];
306
307   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
308
309     rawReader->Reset();
310     AliTOFRawStream tofInput(rawReader);
311     rawReader->Select(5, indexDDL, indexDDL);
312
313     while(tofInput.Next()) {
314
315       detectorIndex[0] = (Int_t)tofInput.GetSector();
316       detectorIndex[1] = (Int_t)tofInput.GetPlate();
317       detectorIndex[2] = (Int_t)tofInput.GetStrip();
318       detectorIndex[3] = (Int_t)tofInput.GetPadZ();
319       detectorIndex[4] = (Int_t)tofInput.GetPadX();
320
321       for (ii=0; ii<3; ii++) position[ii] =  0.;
322
323       fTOFGeometry->GetPos(detectorIndex, position);
324       
325       cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
326       cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
327       cylindricalPosition[2] = (Double_t)position[2];
328       cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
329       cylindricalPosition[4] = (Double_t)tofInput.GetADCbin();
330
331       AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
332       InsertCluster(tofCluster);
333
334     } // while loop
335
336   } // DDL Loop
337
338   AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
339
340   FillRecPoint();
341
342   fTreeR->Fill();
343   ResetRecpoint();
344
345   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
346   fTOFLoader->WriteRecPoints("OVERWRITE");
347   
348 }
349 //______________________________________________________________________________
350
351 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
352 {
353   //
354   // Converts RAW data to MC digits for TOF
355   //
356   //             (temporary solution)
357   //
358
359   fRunLoader->GetEvent(iEvent);
360
361   fTreeD = fTOFLoader->TreeD();
362   if (fTreeD)
363     {
364     AliInfo("AliTOFClusterFinder: TreeD re-creation");
365     fTreeD = 0x0;
366     fTOFLoader->MakeTree("D");
367     fTreeD = fTOFLoader->TreeD();
368     }
369
370
371   fTreeR = fTOFLoader->TreeD();
372   if (fTreeD == 0x0)
373     {
374       fTOFLoader->MakeTree("D");
375       fTreeD = fTOFLoader->TreeD();
376     }
377
378   TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
379   Int_t bufsize = 32000;
380   fTreeD->Branch("TOF", &tofDigits, bufsize);
381
382
383   const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
384
385   fRunLoader->GetEvent(iEvent);
386
387   AliDebug(2,Form(" Event number %2i ", iEvent));
388
389   Int_t indexDDL = 0;
390
391   Int_t detectorIndex[5];
392   Float_t digit[2];
393
394   for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
395
396     rawReader->Reset();
397     AliTOFRawStream tofInput(rawReader);
398     rawReader->Select(5, indexDDL, indexDDL);
399
400     while(tofInput.Next()) {
401
402       detectorIndex[0] = tofInput.GetSector();
403       detectorIndex[1] = tofInput.GetPlate();
404       detectorIndex[2] = tofInput.GetStrip();
405       detectorIndex[3] = tofInput.GetPadX();
406       detectorIndex[4] = tofInput.GetPadZ();
407       
408       digit[0] = (Float_t)tofInput.GetTofBin();
409       digit[1] = (Float_t)tofInput.GetADCbin();
410
411       Int_t tracknum[3]={-1,-1,-1};
412
413       TClonesArray &aDigits = *tofDigits;
414       Int_t last=tofDigits->GetEntriesFast();
415       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
416
417     } // while loop
418
419   } // DDL Loop
420
421   fTreeD->Fill();
422
423   fTOFLoader = fRunLoader->GetLoader("TOFLoader");
424   fTOFLoader->WriteDigits("OVERWRITE");
425   
426 }
427 //______________________________________________________________________________
428
429 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
430   //---------------------------------------------------------------------------//
431   // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
432   //---------------------------------------------------------------------------//
433   if (fNumberOfTofClusters==kTofMaxCluster) {
434     AliError("Too many clusters !");
435     return 1;
436   }
437
438   if (fNumberOfTofClusters==0) {
439     fTofClusters[fNumberOfTofClusters++] = tofCluster;
440     return 0;
441   }
442
443   Int_t ii = FindClusterIndex(tofCluster->GetZ());
444   memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
445   fTofClusters[ii] = tofCluster;
446   fNumberOfTofClusters++;
447   
448   return 0;
449
450 }
451 //_________________________________________________________________________
452
453 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
454   //--------------------------------------------------------------------
455   // This function returns the index of the nearest cluster 
456   //--------------------------------------------------------------------
457   if (fNumberOfTofClusters==0) return 0;
458   if (z <= fTofClusters[0]->GetZ()) return 0;
459   if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
460   Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
461   for (; b<e; m=(b+e)/2) {
462     if (z > fTofClusters[m]->GetZ()) b=m+1;
463     else e=m;
464   }
465
466   return m;
467
468 }
469 //_________________________________________________________________________
470
471 void AliTOFClusterFinder::FillRecPoint()
472 {
473   //
474   // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
475   // in Z) in the global TClonesArray of AliTOFcluster,
476   // i.e. fRecPoints.
477   //
478
479   Int_t ii, jj;
480
481   Int_t detectorIndex[5];
482   Double_t cylindricalPosition[5];
483   Int_t trackLabels[3];
484   Int_t digitIndex = -1;
485
486   TClonesArray &lRecPoints = *fRecPoints;
487   
488   for (ii=0; ii<fNumberOfTofClusters; ii++) {
489
490     digitIndex = fTofClusters[ii]->GetIndex();
491     for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
492     for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
493     cylindricalPosition[0] = fTofClusters[ii]->GetR();
494     cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
495     cylindricalPosition[2] = fTofClusters[ii]->GetZ();
496     cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
497     cylindricalPosition[4] = fTofClusters[ii]->GetADC();
498
499     new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex);
500
501     //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]));
502
503   } // loop on clusters
504
505 }
506 //______________________________________________________________________________
507
508 void AliTOFClusterFinder::ResetRecpoint()
509 {
510   //
511   // Clear the list of reconstructed points
512   //
513
514   fNumberOfTofClusters = 0;
515   if (fRecPoints) fRecPoints->Clear();
516
517 }
518 //______________________________________________________________________________
519
520 void AliTOFClusterFinder::Load()
521 {
522   //
523   // Load TOF.Digits.root and TOF.RecPoints.root files
524   //
525
526   fTOFLoader->LoadDigits("READ");
527   fTOFLoader->LoadRecPoints("recreate");
528
529 }
530 //______________________________________________________________________________
531
532 void AliTOFClusterFinder::LoadClusters()
533 {
534   //
535   // Load TOF.RecPoints.root file
536   //
537
538   fTOFLoader->LoadRecPoints("recreate");
539
540 }
541 //______________________________________________________________________________
542
543 void AliTOFClusterFinder::UnLoad()
544 {
545   //
546   // Unload TOF.Digits.root and TOF.RecPoints.root files
547   //
548
549   fTOFLoader->UnloadDigits();
550   fTOFLoader->UnloadRecPoints();
551
552 }
553 //______________________________________________________________________________
554
555 void AliTOFClusterFinder::UnLoadClusters()
556 {
557   //
558   // Unload TOF.RecPoints.root file
559   //
560
561   fTOFLoader->UnloadRecPoints();
562
563 }
564 //______________________________________________________________________________