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 "TClonesArray.h"
57 #include "AliLoader.h"
59 #include "AliRawReader.h"
60 #include "AliRunLoader.h"
62 #include "AliTOFcalib.h"
63 #include "AliTOFCal.h"
64 #include "AliTOFChannel.h"
65 #include "AliTOFClusterFinder.h"
66 #include "AliTOFcluster.h"
67 #include "AliTOFdigit.h"
68 #include "AliTOFGeometryV5.h"
69 #include "AliTOFGeometry.h"
70 #include "AliTOFRawStream.h"
74 ClassImp(AliTOFClusterFinder)
76 AliTOFClusterFinder::AliTOFClusterFinder():
81 fDigits(new TClonesArray("AliTOFdigit", 4000)),
82 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
83 fNumberOfTofClusters(0)
89 fTOFGeometry = new AliTOFGeometryV5();
90 AliInfo("V5 TOF Geometry is taken as the default");
93 //______________________________________________________________________________
95 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
96 fRunLoader(runLoader),
97 fTOFLoader(runLoader->GetLoader("TOFLoader")),
100 fDigits(new TClonesArray("AliTOFdigit", 4000)),
101 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
102 fNumberOfTofClusters(0)
108 runLoader->CdGAFile();
109 TFile *in=(TFile*)gFile;
111 fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
115 //------------------------------------------------------------------------
116 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
120 this->fDigits=source.fDigits;
121 this->fRecPoints=source.fRecPoints;
122 this->fTOFGeometry=source.fTOFGeometry;
126 //------------------------------------------------------------------------
127 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
130 this->fDigits=source.fDigits;
131 this->fRecPoints=source.fRecPoints;
132 this->fTOFGeometry=source.fTOFGeometry;
136 //______________________________________________________________________________
138 AliTOFClusterFinder::~AliTOFClusterFinder()
153 fRecPoints->Delete();
161 //______________________________________________________________________________
163 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
166 // Converts digits to recpoints for TOF
169 fRunLoader->GetEvent(iEvent);
171 fTreeD = fTOFLoader->TreeD();
174 AliFatal("AliTOFClusterFinder: Can not get TreeD");
177 TBranch *branch = fTreeD->GetBranch("TOF");
179 AliError("can't get the branch with the TOF digits !");
183 TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
184 branch->SetAddress(&digits);
188 fTreeR = fTOFLoader->TreeR();
191 fTOFLoader->MakeTree("R");
192 fTreeR = fTOFLoader->TreeR();
195 Int_t bufsize = 32000;
196 fTreeR->Branch("TOF", &fRecPoints, bufsize);
199 Int_t nDigits = digits->GetEntriesFast();
200 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
208 for (ii=0; ii<nDigits; ii++) {
209 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
210 dig[0]=d->GetSector();
211 dig[1]=d->GetPlate();
212 dig[2]=d->GetStrip();
216 for (jj=0; jj<3; jj++) g[jj] = 0.;
217 fTOFGeometry->GetPos(dig,g);
219 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
220 h[1] = TMath::ATan2(g[1],g[0]);
225 tTdcND = d->GetTdcND();
227 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
228 InsertCluster(tofCluster);
232 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
240 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
241 fTOFLoader->WriteRecPoints("OVERWRITE");
244 //______________________________________________________________________________
246 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
250 // Converts RAW data to recpoints for TOF
253 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
257 Int_t bufsize = 32000;
258 clustersTree->Branch("TOF", &fRecPoints, bufsize);
263 Int_t detectorIndex[5];
265 Double_t cylindricalPosition[5];
269 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
272 AliTOFRawStream tofInput(rawReader);
273 rawReader->Select(5, indexDDL, indexDDL);
275 while(tofInput.Next()) {
277 detectorIndex[0] = tofInput.GetSector();
278 detectorIndex[1] = tofInput.GetPlate();
279 detectorIndex[2] = tofInput.GetStrip();
280 detectorIndex[3] = tofInput.GetPadZ();
281 detectorIndex[4] = tofInput.GetPadX();
283 for (ii=0; ii<3; ii++) position[ii] = 0.;
285 fTOFGeometry->GetPos(detectorIndex, position);
287 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
288 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
289 cylindricalPosition[2] = position[2];
290 cylindricalPosition[3] = tofInput.GetTofBin();
291 cylindricalPosition[4] = tofInput.GetADCbin();
292 tToT = tofInput.GetADCbin();
294 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
295 tofCluster->SetToT(tToT);
296 tofCluster->SetTDCND(tTdcND);
297 InsertCluster(tofCluster);
301 } // loop on DDL files
303 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
308 clustersTree->Fill();
313 //______________________________________________________________________________
315 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
318 // Converts RAW data to recpoints for TOF
321 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
323 fRunLoader->GetEvent(iEvent);
325 AliDebug(2,Form(" Event number %2i ", iEvent));
327 fTreeR = fTOFLoader->TreeR();
330 fTOFLoader->MakeTree("R");
331 fTreeR = fTOFLoader->TreeR();
334 Int_t bufsize = 32000;
335 fTreeR->Branch("TOF", &fRecPoints, bufsize);
340 Int_t detectorIndex[5];
342 Double_t cylindricalPosition[5];
346 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
349 AliTOFRawStream tofInput(rawReader);
350 rawReader->Select(5, indexDDL, indexDDL);
352 while(tofInput.Next()) {
354 detectorIndex[0] = (Int_t)tofInput.GetSector();
355 detectorIndex[1] = (Int_t)tofInput.GetPlate();
356 detectorIndex[2] = (Int_t)tofInput.GetStrip();
357 detectorIndex[3] = (Int_t)tofInput.GetPadZ();
358 detectorIndex[4] = (Int_t)tofInput.GetPadX();
360 for (ii=0; ii<3; ii++) position[ii] = 0.;
362 fTOFGeometry->GetPos(detectorIndex, position);
364 cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
365 cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
366 cylindricalPosition[2] = (Double_t)position[2];
367 cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
368 cylindricalPosition[4] = (Double_t)tofInput.GetADCbin();
369 tToT = tofInput.GetADCbin();
372 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
373 tofCluster->SetToT(tToT);
374 tofCluster->SetTDCND(tTdcND);
375 InsertCluster(tofCluster);
381 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
389 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
390 fTOFLoader->WriteRecPoints("OVERWRITE");
393 //______________________________________________________________________________
395 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
398 // Converts RAW data to MC digits for TOF
400 // (temporary solution)
403 fRunLoader->GetEvent(iEvent);
405 fTreeD = fTOFLoader->TreeD();
408 AliInfo("AliTOFClusterFinder: TreeD re-creation");
410 fTOFLoader->MakeTree("D");
411 fTreeD = fTOFLoader->TreeD();
415 fTreeR = fTOFLoader->TreeD();
418 fTOFLoader->MakeTree("D");
419 fTreeD = fTOFLoader->TreeD();
422 TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
423 Int_t bufsize = 32000;
424 fTreeD->Branch("TOF", &tofDigits, bufsize);
427 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
429 fRunLoader->GetEvent(iEvent);
431 AliDebug(2,Form(" Event number %2i ", iEvent));
435 Int_t detectorIndex[5];
438 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
441 AliTOFRawStream tofInput(rawReader);
442 rawReader->Select(5, indexDDL, indexDDL);
444 while(tofInput.Next()) {
446 detectorIndex[0] = tofInput.GetSector();
447 detectorIndex[1] = tofInput.GetPlate();
448 detectorIndex[2] = tofInput.GetStrip();
449 detectorIndex[3] = tofInput.GetPadX();
450 detectorIndex[4] = tofInput.GetPadZ();
452 digit[0] = (Float_t)tofInput.GetTofBin();
453 digit[1] = (Float_t)tofInput.GetADCbin();
455 Int_t tracknum[3]={-1,-1,-1};
457 TClonesArray &aDigits = *tofDigits;
458 Int_t last=tofDigits->GetEntriesFast();
459 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
467 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
468 fTOFLoader->WriteDigits("OVERWRITE");
471 //______________________________________________________________________________
473 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
474 //---------------------------------------------------------------------------//
475 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
476 //---------------------------------------------------------------------------//
477 if (fNumberOfTofClusters==kTofMaxCluster) {
478 AliError("Too many clusters !");
482 if (fNumberOfTofClusters==0) {
483 fTofClusters[fNumberOfTofClusters++] = tofCluster;
487 Int_t ii = FindClusterIndex(tofCluster->GetZ());
488 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
489 fTofClusters[ii] = tofCluster;
490 fNumberOfTofClusters++;
495 //_________________________________________________________________________
497 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
498 //--------------------------------------------------------------------
499 // This function returns the index of the nearest cluster
500 //--------------------------------------------------------------------
501 if (fNumberOfTofClusters==0) return 0;
502 if (z <= fTofClusters[0]->GetZ()) return 0;
503 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
504 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
505 for (; b<e; m=(b+e)/2) {
506 if (z > fTofClusters[m]->GetZ()) b=m+1;
513 //_________________________________________________________________________
515 void AliTOFClusterFinder::FillRecPoint()
518 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
519 // in Z) in the global TClonesArray of AliTOFcluster,
525 Int_t detectorIndex[5];
526 Double_t cylindricalPosition[5];
527 Int_t trackLabels[3];
528 Int_t digitIndex = -1;
532 TClonesArray &lRecPoints = *fRecPoints;
534 for (ii=0; ii<fNumberOfTofClusters; ii++) {
536 digitIndex = fTofClusters[ii]->GetIndex();
537 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
538 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
539 cylindricalPosition[0] = fTofClusters[ii]->GetR();
540 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
541 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
542 cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
543 cylindricalPosition[4] = fTofClusters[ii]->GetADC();
544 tToT = fTofClusters[ii]->GetToT();
545 tTdcND = fTofClusters[ii]->GetTDCND();
547 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND);
549 //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]));
551 } // loop on clusters
555 //_________________________________________________________________________
556 void AliTOFClusterFinder::CalibrateRecPoint()
559 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
560 // in Z) in the global TClonesArray of AliTOFcluster,
566 Int_t detectorIndex[5];
567 Int_t digitIndex = -1;
570 AliInfo(" Calibrating TOF Clusters: ")
571 AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
572 // calib->ReadParFromCDB("TOF/Calib",0); // original
573 calib->ReadParFromCDB("TOF/Calib",-1); // Use AliCDBManager's run number
574 AliTOFCal *calTOFArray = calib->GetTOFCalArray();
576 for (ii=0; ii<fNumberOfTofClusters; ii++) {
577 digitIndex = fTofClusters[ii]->GetIndex();
578 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
580 Int_t index = calib->GetIndex(detectorIndex);
582 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
584 for (Int_t j = 0; j<6; j++){
585 par[j]=calChannel->GetSlewPar(j);
587 tToT = fTofClusters[ii]->GetToT();
588 Float_t timeCorr=par[0]+par[1]*tToT+par[2]*tToT*tToT+par[3]*tToT*tToT*tToT+par[4]*tToT*tToT*tToT*tToT+par[5]*tToT*tToT*tToT*tToT*tToT;
589 tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
590 tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
591 fTofClusters[ii]->SetTDC(tdcCorr);
593 } // loop on clusters
597 //______________________________________________________________________________
599 void AliTOFClusterFinder::ResetRecpoint()
602 // Clear the list of reconstructed points
605 fNumberOfTofClusters = 0;
606 if (fRecPoints) fRecPoints->Clear();
609 //______________________________________________________________________________
611 void AliTOFClusterFinder::Load()
614 // Load TOF.Digits.root and TOF.RecPoints.root files
617 fTOFLoader->LoadDigits("READ");
618 fTOFLoader->LoadRecPoints("recreate");
621 //______________________________________________________________________________
623 void AliTOFClusterFinder::LoadClusters()
626 // Load TOF.RecPoints.root file
629 fTOFLoader->LoadRecPoints("recreate");
632 //______________________________________________________________________________
634 void AliTOFClusterFinder::UnLoad()
637 // Unload TOF.Digits.root and TOF.RecPoints.root files
640 fTOFLoader->UnloadDigits();
641 fTOFLoader->UnloadRecPoints();
644 //______________________________________________________________________________
646 void AliTOFClusterFinder::UnLoadClusters()
649 // Unload TOF.RecPoints.root file
652 fTOFLoader->UnloadRecPoints();
655 //______________________________________________________________________________