1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
24 Revision 0.02 2005/07/27 A. De Caro
25 Implement public method
26 Digits2RecPoint(Int_t)
27 to convert digits in clusters
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
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
43 ////////////////////////////////////////////////////////////////
45 // Class for TOF cluster finder //
47 // Starting from Raw Data, create rec points, //
48 // fill TreeR for TOF, //
49 // write TOF.RecPoints.root file //
51 ////////////////////////////////////////////////////////////////
53 #include <Riostream.h>
55 #include <TObjArray.h>
56 #include <TClonesArray.h>
59 #include "AliRunLoader.h"
60 #include "AliLoader.h"
62 #include "AliTOFdigit.h"
63 #include "AliTOFcluster.h"
64 #include "AliTOFGeometry.h"
65 #include "AliTOFGeometryV4.h"
66 #include "AliTOFGeometryV5.h"
67 #include "AliTOFRawStream.h"
69 #include "AliTOFClusterFinder.h"
71 ClassImp(AliTOFClusterFinder)
73 AliTOFClusterFinder::AliTOFClusterFinder():
78 fDigits(new TClonesArray("AliTOFdigit", 4000)),
79 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
80 fNumberOfTofClusters(0)
86 fTOFGeometry = new AliTOFGeometry();
89 //______________________________________________________________________________
91 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
92 fRunLoader(runLoader),
93 fTOFLoader(runLoader->GetLoader("TOFLoader")),
96 fDigits(new TClonesArray("AliTOFdigit", 4000)),
97 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
98 fNumberOfTofClusters(0)
104 runLoader->CdGAFile();
105 TFile *in=(TFile*)gFile;
107 fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
110 //______________________________________________________________________________
112 AliTOFClusterFinder::~AliTOFClusterFinder()
127 fRecPoints->Delete();
135 //______________________________________________________________________________
137 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
140 // Converts digits to recpoints for TOF
143 fRunLoader->GetEvent(iEvent);
145 fTreeD = fTOFLoader->TreeD();
148 AliFatal("AliTOFClusterFinder: Can not get TreeD");
151 TBranch *branch = fTreeD->GetBranch("TOF");
153 AliError("can't get the branch with the TOF digits !");
157 TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
158 branch->SetAddress(&digits);
162 fTreeR = fTOFLoader->TreeR();
165 fTOFLoader->MakeTree("R");
166 fTreeR = fTOFLoader->TreeR();
169 Int_t bufsize = 32000;
170 fTreeR->Branch("TOF", &fRecPoints, bufsize);
173 Int_t nDigits = digits->GetEntriesFast();
174 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
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();
188 for (jj=0; jj<3; jj++) g[jj] = 0.;
189 fTOFGeometry->GetPos(dig,g);
191 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
192 h[1] = TMath::ATan2(g[1],g[0]);
197 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii);
198 InsertCluster(tofCluster);
202 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
209 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
210 fTOFLoader->WriteRecPoints("OVERWRITE");
213 //______________________________________________________________________________
215 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
219 // Converts RAW data to recpoints for TOF
222 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
226 Int_t bufsize = 32000;
227 clustersTree->Branch("TOF", &fRecPoints, bufsize);
232 Int_t detectorIndex[5];
234 Double_t cylindricalPosition[5];
236 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
239 AliTOFRawStream tofInput(rawReader);
240 rawReader->Select(5, indexDDL, indexDDL);
242 while(tofInput.Next()) {
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();
250 for (ii=0; ii<3; ii++) position[ii] = 0.;
252 fTOFGeometry->GetPos(detectorIndex, position);
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();
260 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
261 InsertCluster(tofCluster);
265 } // loop on DDL files
267 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
271 clustersTree->Fill();
275 //______________________________________________________________________________
277 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
280 // Converts RAW data to recpoints for TOF
283 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
285 fRunLoader->GetEvent(iEvent);
287 AliDebug(2,Form(" Event number %2i ", iEvent));
289 fTreeR = fTOFLoader->TreeR();
292 fTOFLoader->MakeTree("R");
293 fTreeR = fTOFLoader->TreeR();
296 Int_t bufsize = 32000;
297 fTreeR->Branch("TOF", &fRecPoints, bufsize);
302 Int_t detectorIndex[5];
304 Double_t cylindricalPosition[5];
306 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
309 AliTOFRawStream tofInput(rawReader);
310 rawReader->Select(5, indexDDL, indexDDL);
312 while(tofInput.Next()) {
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();
320 for (ii=0; ii<3; ii++) position[ii] = 0.;
322 fTOFGeometry->GetPos(detectorIndex, position);
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();
330 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
331 InsertCluster(tofCluster);
337 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
344 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
345 fTOFLoader->WriteRecPoints("OVERWRITE");
348 //______________________________________________________________________________
350 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
353 // Converts RAW data to MC digits for TOF
355 // (temporary solution)
358 fRunLoader->GetEvent(iEvent);
360 fTreeD = fTOFLoader->TreeD();
363 AliInfo("AliTOFClusterFinder: TreeD re-creation");
365 fTOFLoader->MakeTree("D");
366 fTreeD = fTOFLoader->TreeD();
370 fTreeR = fTOFLoader->TreeD();
373 fTOFLoader->MakeTree("D");
374 fTreeD = fTOFLoader->TreeD();
377 TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
378 Int_t bufsize = 32000;
379 fTreeD->Branch("TOF", &tofDigits, bufsize);
382 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
384 fRunLoader->GetEvent(iEvent);
386 AliDebug(2,Form(" Event number %2i ", iEvent));
390 Int_t detectorIndex[5];
393 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
396 AliTOFRawStream tofInput(rawReader);
397 rawReader->Select(5, indexDDL, indexDDL);
399 while(tofInput.Next()) {
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();
407 digit[0] = (Float_t)tofInput.GetTofBin();
408 digit[1] = (Float_t)tofInput.GetADCbin();
410 Int_t tracknum[3]={-1,-1,-1};
412 TClonesArray &aDigits = *tofDigits;
413 Int_t last=tofDigits->GetEntriesFast();
414 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
422 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
423 fTOFLoader->WriteDigits("OVERWRITE");
426 //______________________________________________________________________________
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 !");
437 if (fNumberOfTofClusters==0) {
438 fTofClusters[fNumberOfTofClusters++] = tofCluster;
442 Int_t ii = FindClusterIndex(tofCluster->GetZ());
443 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
444 fTofClusters[ii] = tofCluster;
445 fNumberOfTofClusters++;
450 //_________________________________________________________________________
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;
468 //_________________________________________________________________________
470 void AliTOFClusterFinder::FillRecPoint()
473 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
474 // in Z) in the global TClonesArray of AliTOFcluster,
480 Int_t detectorIndex[5];
481 Double_t cylindricalPosition[5];
482 Int_t trackLabels[3];
483 Int_t digitIndex = -1;
485 TClonesArray &lRecPoints = *fRecPoints;
487 for (ii=0; ii<fNumberOfTofClusters; ii++) {
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();
498 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex);
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]));
502 } // loop on clusters
505 //______________________________________________________________________________
507 void AliTOFClusterFinder::ResetRecpoint()
510 // Clear the list of reconstructed points
513 fNumberOfTofClusters = 0;
514 if (fRecPoints) fRecPoints->Clear();
517 //______________________________________________________________________________
519 void AliTOFClusterFinder::Load()
522 // Load TOF.Digits.root and TOF.RecPoints.root files
525 fTOFLoader->LoadDigits("READ");
526 fTOFLoader->LoadRecPoints("recreate");
529 //______________________________________________________________________________
531 void AliTOFClusterFinder::LoadClusters()
534 // Load TOF.RecPoints.root file
537 fTOFLoader->LoadRecPoints("recreate");
540 //______________________________________________________________________________
542 void AliTOFClusterFinder::UnLoad()
545 // Unload TOF.Digits.root and TOF.RecPoints.root files
548 fTOFLoader->UnloadDigits();
549 fTOFLoader->UnloadRecPoints();
552 //______________________________________________________________________________
554 void AliTOFClusterFinder::UnLoadClusters()
557 // Unload TOF.RecPoints.root file
560 fTOFLoader->UnloadRecPoints();
563 //______________________________________________________________________________