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>
60 #include "AliRunLoader.h"
61 #include "AliLoader.h"
63 #include "AliTOFdigit.h"
64 #include "AliTOFcluster.h"
65 #include "AliTOFGeometry.h"
66 #include "AliTOFGeometryV4.h"
67 #include "AliTOFGeometryV5.h"
68 #include "AliTOFRawStream.h"
70 #include "AliTOFClusterFinder.h"
72 ClassImp(AliTOFClusterFinder)
74 AliTOFClusterFinder::AliTOFClusterFinder():
79 fDigits(new TClonesArray("AliTOFdigit", 4000)),
80 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
81 fNumberOfTofClusters(0)
87 fTOFGeometry = new AliTOFGeometry();
90 //______________________________________________________________________________
92 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
93 fRunLoader(runLoader),
94 fTOFLoader(runLoader->GetLoader("TOFLoader")),
97 fDigits(new TClonesArray("AliTOFdigit", 4000)),
98 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
99 fNumberOfTofClusters(0)
105 runLoader->CdGAFile();
106 TFile *in=(TFile*)gFile;
108 fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
111 //______________________________________________________________________________
113 AliTOFClusterFinder::~AliTOFClusterFinder()
128 fRecPoints->Delete();
136 //______________________________________________________________________________
138 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
141 // Converts digits to recpoints for TOF
144 fRunLoader->GetEvent(iEvent);
146 fTreeD = fTOFLoader->TreeD();
149 AliFatal("AliTOFClusterFinder: Can not get TreeD");
152 TBranch *branch = fTreeD->GetBranch("TOF");
154 AliError("can't get the branch with the TOF digits !");
158 TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
159 branch->SetAddress(&digits);
163 fTreeR = fTOFLoader->TreeR();
166 fTOFLoader->MakeTree("R");
167 fTreeR = fTOFLoader->TreeR();
170 Int_t bufsize = 32000;
171 fTreeR->Branch("TOF", &fRecPoints, bufsize);
174 Int_t nDigits = digits->GetEntriesFast();
175 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
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();
189 for (jj=0; jj<3; jj++) g[jj] = 0.;
190 fTOFGeometry->GetPos(dig,g);
192 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
193 h[1] = TMath::ATan2(g[1],g[0]);
198 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii);
199 InsertCluster(tofCluster);
203 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
210 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
211 fTOFLoader->WriteRecPoints("OVERWRITE");
214 //______________________________________________________________________________
216 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
220 // Converts RAW data to recpoints for TOF
223 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
227 Int_t bufsize = 32000;
228 clustersTree->Branch("TOF", &fRecPoints, bufsize);
233 Int_t detectorIndex[5];
235 Double_t cylindricalPosition[5];
237 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
240 AliTOFRawStream tofInput(rawReader);
241 rawReader->Select(5, indexDDL, indexDDL);
243 while(tofInput.Next()) {
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();
251 for (ii=0; ii<3; ii++) position[ii] = 0.;
253 fTOFGeometry->GetPos(detectorIndex, position);
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();
261 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
262 InsertCluster(tofCluster);
266 } // loop on DDL files
268 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
272 clustersTree->Fill();
276 //______________________________________________________________________________
278 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
281 // Converts RAW data to recpoints for TOF
284 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
286 fRunLoader->GetEvent(iEvent);
288 AliDebug(2,Form(" Event number %2i ", iEvent));
290 fTreeR = fTOFLoader->TreeR();
293 fTOFLoader->MakeTree("R");
294 fTreeR = fTOFLoader->TreeR();
297 Int_t bufsize = 32000;
298 fTreeR->Branch("TOF", &fRecPoints, bufsize);
303 Int_t detectorIndex[5];
305 Double_t cylindricalPosition[5];
307 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
310 AliTOFRawStream tofInput(rawReader);
311 rawReader->Select(5, indexDDL, indexDDL);
313 while(tofInput.Next()) {
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();
321 for (ii=0; ii<3; ii++) position[ii] = 0.;
323 fTOFGeometry->GetPos(detectorIndex, position);
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();
331 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
332 InsertCluster(tofCluster);
338 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
345 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
346 fTOFLoader->WriteRecPoints("OVERWRITE");
349 //______________________________________________________________________________
351 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
354 // Converts RAW data to MC digits for TOF
356 // (temporary solution)
359 fRunLoader->GetEvent(iEvent);
361 fTreeD = fTOFLoader->TreeD();
364 AliInfo("AliTOFClusterFinder: TreeD re-creation");
366 fTOFLoader->MakeTree("D");
367 fTreeD = fTOFLoader->TreeD();
371 fTreeR = fTOFLoader->TreeD();
374 fTOFLoader->MakeTree("D");
375 fTreeD = fTOFLoader->TreeD();
378 TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
379 Int_t bufsize = 32000;
380 fTreeD->Branch("TOF", &tofDigits, bufsize);
383 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
385 fRunLoader->GetEvent(iEvent);
387 AliDebug(2,Form(" Event number %2i ", iEvent));
391 Int_t detectorIndex[5];
394 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
397 AliTOFRawStream tofInput(rawReader);
398 rawReader->Select(5, indexDDL, indexDDL);
400 while(tofInput.Next()) {
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();
408 digit[0] = (Float_t)tofInput.GetTofBin();
409 digit[1] = (Float_t)tofInput.GetADCbin();
411 Int_t tracknum[3]={-1,-1,-1};
413 TClonesArray &aDigits = *tofDigits;
414 Int_t last=tofDigits->GetEntriesFast();
415 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
423 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
424 fTOFLoader->WriteDigits("OVERWRITE");
427 //______________________________________________________________________________
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 !");
438 if (fNumberOfTofClusters==0) {
439 fTofClusters[fNumberOfTofClusters++] = tofCluster;
443 Int_t ii = FindClusterIndex(tofCluster->GetZ());
444 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
445 fTofClusters[ii] = tofCluster;
446 fNumberOfTofClusters++;
451 //_________________________________________________________________________
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;
469 //_________________________________________________________________________
471 void AliTOFClusterFinder::FillRecPoint()
474 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
475 // in Z) in the global TClonesArray of AliTOFcluster,
481 Int_t detectorIndex[5];
482 Double_t cylindricalPosition[5];
483 Int_t trackLabels[3];
484 Int_t digitIndex = -1;
486 TClonesArray &lRecPoints = *fRecPoints;
488 for (ii=0; ii<fNumberOfTofClusters; ii++) {
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();
499 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex);
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]));
503 } // loop on clusters
506 //______________________________________________________________________________
508 void AliTOFClusterFinder::ResetRecpoint()
511 // Clear the list of reconstructed points
514 fNumberOfTofClusters = 0;
515 if (fRecPoints) fRecPoints->Clear();
518 //______________________________________________________________________________
520 void AliTOFClusterFinder::Load()
523 // Load TOF.Digits.root and TOF.RecPoints.root files
526 fTOFLoader->LoadDigits("READ");
527 fTOFLoader->LoadRecPoints("recreate");
530 //______________________________________________________________________________
532 void AliTOFClusterFinder::LoadClusters()
535 // Load TOF.RecPoints.root file
538 fTOFLoader->LoadRecPoints("recreate");
541 //______________________________________________________________________________
543 void AliTOFClusterFinder::UnLoad()
546 // Unload TOF.Digits.root and TOF.RecPoints.root files
549 fTOFLoader->UnloadDigits();
550 fTOFLoader->UnloadRecPoints();
553 //______________________________________________________________________________
555 void AliTOFClusterFinder::UnLoadClusters()
558 // Unload TOF.RecPoints.root file
561 fTOFLoader->UnloadRecPoints();
564 //______________________________________________________________________________