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 1.27 2007/05/02 16:31:49 arcelli
19 Add methods to handle single event reconstruction. retrieval of Calib info moved to AliTOFReconstructor ctor, and passed via a pointer to AliTOFcalib
21 Revision 1.26 2007/04/30 19:02:24 arcelli
22 hopefully the last refinements for correct type conversion in calibration
24 Revision 1.25 2007/04/30 15:22:17 arcelli
25 Change TOF digit Time, Tot etc to int type
27 Revision 1.24 2007/04/27 11:19:31 arcelli
28 updates for the new decoder
30 Revision 1.23 2007/04/23 16:51:39 decaro
31 Digits-to-raw_data conversion: correction for a more real description (A.De Caro, R.Preghenella)
33 Revision 1.22 2007/04/19 17:26:32 arcelli
34 Fix a bug (add some debug printout
36 Revision 1.21 2007/04/18 17:28:12 arcelli
37 Set the ToT bin width to the one actually used...
39 Revision 1.20 2007/03/09 09:57:23 arcelli
40 Remove a forgotten include of Riostrem
42 Revision 1.19 2007/03/08 15:41:20 arcelli
43 set uncorrected times when filling RecPoints
45 Revision 1.18 2007/03/06 16:31:20 arcelli
46 Add Uncorrected TOF Time signal
48 Revision 1.17 2007/02/28 18:09:11 arcelli
49 Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
51 Revision 1.16 2007/02/20 15:57:00 decaro
52 Raw data update: to read the TOF raw data defined in UNPACKED mode
55 Revision 0.03 2005/07/28 A. De Caro
56 Implement public method
57 Raw2Digits(Int_t, AliRawReader *)
58 to convert digits from raw data in MC digits
61 Revision 0.02 2005/07/27 A. De Caro
62 Implement public method
63 Digits2RecPoint(Int_t)
64 to convert digits in clusters
66 Revision 0.02 2005/07/26 A. De Caro
67 Implement private methods
68 InsertCluster(AliTOFcluster *)
69 FindClusterIndex(Double_t)
70 originally implemented in AliTOFtracker
71 by S. Arcelli and C. Zampolli
73 Revision 0.01 2005/07/25 A. De Caro
74 Implement public methods
75 Digits2RecPoint(AliRawReader *, TTree *)
76 Digits2RecPoint(Int_t, AliRawReader *)
77 to convert raw data in clusters
80 ////////////////////////////////////////////////////////////////
82 // Class for TOF cluster finder //
84 // Starting from Raw Data, create rec points, //
85 // fill TreeR for TOF, //
86 // write TOF.RecPoints.root file //
88 ////////////////////////////////////////////////////////////////
91 #include "TClonesArray.h"
93 #include "TStopwatch.h"
97 #include "AliLoader.h"
99 #include "AliRawReader.h"
100 #include "AliRunLoader.h"
102 #include "AliTOFCal.h"
103 #include "AliTOFcalib.h"
104 #include "AliTOFChannel.h"
105 #include "AliTOFClusterFinder.h"
106 #include "AliTOFcluster.h"
107 #include "AliTOFdigit.h"
108 #include "AliTOFGeometry.h"
109 #include "AliTOFGeometryV5.h"
110 #include "AliTOFrawData.h"
111 #include "AliTOFRawStream.h"
112 #include "Riostream.h"
114 //extern TFile *gFile;
116 ClassImp(AliTOFClusterFinder)
118 AliTOFClusterFinder::AliTOFClusterFinder(AliTOFcalib *calib):
123 fTOFGeometry(new AliTOFGeometryV5()),
124 fDigits(new TClonesArray("AliTOFdigit", 4000)),
125 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
126 fNumberOfTofClusters(0),
135 AliInfo("V5 TOF Geometry is taken as the default");
138 //______________________________________________________________________________
140 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader, AliTOFcalib *calib):
141 fRunLoader(runLoader),
142 fTOFLoader(runLoader->GetLoader("TOFLoader")),
145 fTOFGeometry(new AliTOFGeometryV5()),
146 fDigits(new TClonesArray("AliTOFdigit", 4000)),
147 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
148 fNumberOfTofClusters(0),
159 //------------------------------------------------------------------------
160 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
166 fTOFGeometry(new AliTOFGeometryV5()),
167 fDigits(new TClonesArray("AliTOFdigit", 4000)),
168 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
169 fNumberOfTofClusters(0),
175 this->fDigits=source.fDigits;
176 this->fRecPoints=source.fRecPoints;
177 this->fTOFGeometry=source.fTOFGeometry;
178 this->fDecoderVersion=source.fDecoderVersion;
179 this->fTOFcalib=source.fTOFcalib;
183 //------------------------------------------------------------------------
184 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
187 this->fDigits=source.fDigits;
188 this->fRecPoints=source.fRecPoints;
189 this->fTOFGeometry=source.fTOFGeometry;
190 this->fVerbose=source.fVerbose;
191 this->fDecoderVersion=source.fDecoderVersion;
192 this->fTOFcalib=source.fTOFcalib;
196 //______________________________________________________________________________
198 AliTOFClusterFinder::~AliTOFClusterFinder()
213 fRecPoints->Delete();
221 //______________________________________________________________________________
223 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
226 // Converts digits to recpoints for TOF
229 TStopwatch stopwatch;
232 fRunLoader->GetEvent(iEvent);
234 fTreeD = fTOFLoader->TreeD();
237 AliFatal("AliTOFClusterFinder: Can not get TreeD");
240 TBranch *branch = fTreeD->GetBranch("TOF");
242 AliError("can't get the branch with the TOF digits !");
246 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
247 branch->SetAddress(&digits);
251 fTreeR = fTOFLoader->TreeR();
254 fTOFLoader->MakeTree("R");
255 fTreeR = fTOFLoader->TreeR();
258 Int_t bufsize = 32000;
259 fTreeR->Branch("TOF", &fRecPoints, bufsize);
262 Int_t nDigits = digits->GetEntriesFast();
263 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
266 Int_t dig[5]; //cluster detector indeces
267 Float_t g[3]; //cluster cartesian coord
268 Double_t h[3]; // the cluster spatial cyl. coordinates
269 Int_t parTOF[5]; //The TOF signal parameters
270 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
271 for (ii=0; ii<nDigits; ii++) {
272 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
273 dig[0]=d->GetSector();
274 dig[1]=d->GetPlate();
275 dig[2]=d->GetStrip();
279 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
281 for (jj=0; jj<3; jj++) g[jj] = 0.;
282 fTOFGeometry->GetPos(dig,g);
284 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
285 h[1] = TMath::ATan2(g[1],g[0]);
288 parTOF[0] = d->GetTdc(); //the TDC signal
289 parTOF[1] = d->GetToT(); //the ToT signal
290 parTOF[2] = d->GetAdc(); // the adc charge
291 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
292 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
293 AliTOFcluster *tofCluster = new AliTOFcluster(h,dig,parTOF,status,d->GetTracks(),ii);
294 InsertCluster(tofCluster);
298 AliInfo(Form("Number of found clusters: %i for event: %i", fNumberOfTofClusters, iEvent));
306 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
307 fTOFLoader->WriteRecPoints("OVERWRITE");
309 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
310 stopwatch.RealTime(),stopwatch.CpuTime()));
316 //______________________________________________________________________________
318 void AliTOFClusterFinder::Digits2RecPoints(TTree* digitsTree, TTree* clusterTree)
321 // Converts digits to recpoints for TOF
324 TStopwatch stopwatch;
327 /// fRunLoader->GetEvent(iEvent);
329 if (digitsTree == 0x0)
331 AliFatal("AliTOFClusterFinder: Can not get TreeD");
334 TBranch *branch = digitsTree->GetBranch("TOF");
336 AliError("can't get the branch with the TOF digits !");
340 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
341 branch->SetAddress(&digits);
346 Int_t bufsize = 32000;
347 fTreeR->Branch("TOF", &fRecPoints, bufsize);
349 digitsTree->GetEvent(0);
350 Int_t nDigits = digits->GetEntriesFast();
351 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
354 Int_t dig[5]; //cluster detector indeces
355 Float_t g[3]; //cluster cartesian coord
356 Double_t h[3]; // the cluster spatial cyl. coordinates
357 Int_t parTOF[5]; //The TOF signal parameters
358 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
359 for (ii=0; ii<nDigits; ii++) {
360 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
361 dig[0]=d->GetSector();
362 dig[1]=d->GetPlate();
363 dig[2]=d->GetStrip();
367 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
369 for (jj=0; jj<3; jj++) g[jj] = 0.;
370 fTOFGeometry->GetPos(dig,g);
372 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
373 h[1] = TMath::ATan2(g[1],g[0]);
376 parTOF[0] = d->GetTdc(); //the TDC signal
377 parTOF[1] = d->GetToT(); //the ToT signal
378 parTOF[2] = d->GetAdc(); // the adc charge
379 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
380 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
381 AliTOFcluster *tofCluster = new AliTOFcluster(h,dig,parTOF,status,d->GetTracks(),ii);
382 InsertCluster(tofCluster);
386 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
394 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
395 stopwatch.RealTime(),stopwatch.CpuTime()));
400 //______________________________________________________________________________
402 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
406 // Converts RAW data to recpoints for TOF
409 TStopwatch stopwatch;
412 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
413 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
417 Int_t bufsize = 32000;
418 clustersTree->Branch("TOF", &fRecPoints, bufsize);
420 TClonesArray * clonesRawData;
425 Int_t detectorIndex[5];
427 Double_t cylindricalPosition[3];
431 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
433 AliTOFRawStream tofInput(rawReader);
436 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
439 if (fDecoderVersion) {
440 AliInfo("Using New Decoder \n");
441 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
443 else tofInput.LoadRawData(indexDDL);
445 clonesRawData = (TClonesArray*)tofInput.GetRawData();
447 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
449 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
451 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
454 if (indexDDL<10) ftxt << " " << indexDDL;
455 else ftxt << " " << indexDDL;
456 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
457 else ftxt << " " << tofRawDatum->GetTRM();
458 ftxt << " " << tofRawDatum->GetTRMchain();
459 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
460 else ftxt << " " << tofRawDatum->GetTDC();
461 ftxt << " " << tofRawDatum->GetTDCchannel();
464 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
465 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
466 dummy = detectorIndex[3];
467 detectorIndex[3] = detectorIndex[4];
468 detectorIndex[4] = dummy;
471 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
472 else ftxt << " -> " << detectorIndex[0];
473 ftxt << " " << detectorIndex[1];
474 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
475 else ftxt << " " << detectorIndex[2];
476 ftxt << " " << detectorIndex[3];
477 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
478 else ftxt << " " << detectorIndex[4];
481 for (ii=0; ii<3; ii++) position[ii] = 0.;
482 fTOFGeometry->GetPos(detectorIndex, position);
484 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
485 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
486 cylindricalPosition[2] = position[2];
488 parTOF[0] = tofRawDatum->GetTOF(); //TDC
489 parTOF[1] = tofRawDatum->GetTOT(); // TOT
490 parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
491 parTOF[3] = -1;//raw data: no track of undecalib sim time
492 parTOF[4] = tofRawDatum->GetTOF(); // RAW time
493 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex, parTOF);
494 InsertCluster(tofCluster);
497 if (parTOF[1]<10)ftxt << " " << parTOF[1];
498 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
499 else ftxt << " " << parTOF[1];
500 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
501 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
502 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
503 else ftxt << " " << parTOF[3] << endl;
506 } // closed loop on TOF raw data per current DDL file
508 clonesRawData->Clear();
510 } // closed loop on DDL index
512 if (fVerbose==2) ftxt.close();
514 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
519 clustersTree->Fill();
523 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
524 stopwatch.RealTime(),stopwatch.CpuTime()));
527 //______________________________________________________________________________
529 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
532 // Converts RAW data to recpoints for TOF
535 TStopwatch stopwatch;
538 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
539 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
541 fRunLoader->GetEvent(iEvent);
543 AliDebug(2,Form(" Event number %2i ", iEvent));
545 fTreeR = fTOFLoader->TreeR();
548 fTOFLoader->MakeTree("R");
549 fTreeR = fTOFLoader->TreeR();
552 Int_t bufsize = 32000;
553 fTreeR->Branch("TOF", &fRecPoints, bufsize);
555 TClonesArray * clonesRawData;
560 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
562 Double_t cylindricalPosition[5];
565 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
567 AliTOFRawStream tofInput(rawReader);
570 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
573 if (fDecoderVersion) {
574 AliInfo("Using New Decoder \n");
575 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
577 else tofInput.LoadRawData(indexDDL);
579 clonesRawData = (TClonesArray*)tofInput.GetRawData();
581 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
583 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
585 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
588 if (indexDDL<10) ftxt << " " << indexDDL;
589 else ftxt << " " << indexDDL;
590 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
591 else ftxt << " " << tofRawDatum->GetTRM();
592 ftxt << " " << tofRawDatum->GetTRMchain();
593 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
594 else ftxt << " " << tofRawDatum->GetTDC();
595 ftxt << " " << tofRawDatum->GetTDCchannel();
598 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
599 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
600 dummy = detectorIndex[3];
601 detectorIndex[3] = detectorIndex[4];
602 detectorIndex[4] = dummy;
605 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
606 else ftxt << " -> " << detectorIndex[0];
607 ftxt << " " << detectorIndex[1];
608 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
609 else ftxt << " " << detectorIndex[2];
610 ftxt << " " << detectorIndex[3];
611 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
612 else ftxt << " " << detectorIndex[4];
615 for (ii=0; ii<3; ii++) position[ii] = 0.;
616 fTOFGeometry->GetPos(detectorIndex, position);
618 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
619 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
620 cylindricalPosition[2] = position[2];
621 parTOF[0] = tofRawDatum->GetTOF(); // TDC
622 parTOF[1] = tofRawDatum->GetTOT(); // TOT
623 parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
624 parTOF[3] = -1; //raw data: no track of the undecalib sim time
625 parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
626 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex, parTOF);
627 InsertCluster(tofCluster);
630 if (parTOF[1]<10)ftxt << " " << parTOF[1];
631 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
632 else ftxt << " " << parTOF[1];
633 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
634 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
635 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
636 else ftxt << " " << parTOF[3] << endl;
639 } // closed loop on TOF raw data per current DDL file
641 clonesRawData->Clear();
643 } // closed loop on DDL index
645 if (fVerbose==2) ftxt.close();
647 AliInfo(Form("Number of found clusters: %i for event: %i", fNumberOfTofClusters, iEvent));
655 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
656 fTOFLoader->WriteRecPoints("OVERWRITE");
658 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
659 stopwatch.RealTime(),stopwatch.CpuTime()));
662 //______________________________________________________________________________
664 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
667 // Converts RAW data to MC digits for TOF
669 // (temporary solution)
672 TStopwatch stopwatch;
675 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
676 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
678 fRunLoader->GetEvent(iEvent);
680 fTreeD = fTOFLoader->TreeD();
683 AliInfo("TreeD re-creation");
685 fTOFLoader->MakeTree("D");
686 fTreeD = fTOFLoader->TreeD();
689 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
690 Int_t bufsize = 32000;
691 fTreeD->Branch("TOF", &tofDigits, bufsize);
693 fRunLoader->GetEvent(iEvent);
695 AliDebug(2,Form(" Event number %2i ", iEvent));
697 TClonesArray * clonesRawData;
701 Int_t detectorIndex[5];
704 AliTOFRawStream tofInput(rawReader);
707 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
710 if (fDecoderVersion) {
711 AliInfo("Using New Decoder \n");
712 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
714 else tofInput.LoadRawData(indexDDL);
716 clonesRawData = (TClonesArray*)tofInput.GetRawData();
718 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
720 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
722 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
724 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
725 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
726 dummy = detectorIndex[3];
727 detectorIndex[3] = detectorIndex[4];
728 detectorIndex[4] = dummy;
730 digit[0] = tofInput.GetTofBin();
731 digit[1] = tofInput.GetToTbin();
732 digit[2] = tofInput.GetToTbin();
735 Int_t tracknum[3]={-1,-1,-1};
737 TClonesArray &aDigits = *tofDigits;
738 Int_t last=tofDigits->GetEntriesFast();
739 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
743 clonesRawData->Clear();
749 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
750 fTOFLoader->WriteDigits("OVERWRITE");
752 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
753 stopwatch.RealTime(),stopwatch.CpuTime()));
757 //______________________________________________________________________________
759 void AliTOFClusterFinder::Raw2Digits(AliRawReader *rawReader, TTree* digitsTree)
762 // Converts RAW data to MC digits for TOF for the current event
766 TStopwatch stopwatch;
769 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
773 AliError("No input digits Tree");
777 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
778 Int_t bufsize = 32000;
779 digitsTree->Branch("TOF", &tofDigits, bufsize);
781 /// fRunLoader->GetEvent(iEvent);
783 /// AliDebug(2,Form(" Event number %2i ", iEvent));
785 TClonesArray * clonesRawData;
789 Int_t detectorIndex[5];
792 AliTOFRawStream tofInput(rawReader);
795 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
798 if (fDecoderVersion) {
799 AliInfo("Using New Decoder \n");
800 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
802 else tofInput.LoadRawData(indexDDL);
804 clonesRawData = (TClonesArray*)tofInput.GetRawData();
806 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
808 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
810 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
812 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
813 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
814 dummy = detectorIndex[3];
815 detectorIndex[3] = detectorIndex[4];
816 detectorIndex[4] = dummy;
818 digit[0] = tofInput.GetTofBin();
819 digit[1] = tofInput.GetToTbin();
820 digit[2] = tofInput.GetToTbin();
823 Int_t tracknum[3]={-1,-1,-1};
825 TClonesArray &aDigits = *tofDigits;
826 Int_t last=tofDigits->GetEntriesFast();
827 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
831 clonesRawData->Clear();
837 AliDebug(1, Form("Got %d digits: ", tofDigits->GetEntries()));
838 AliDebug(1, Form("Execution time to read TOF raw data and fill TOF digit tree : R:%.2fs C:%.2fs",
839 stopwatch.RealTime(),stopwatch.CpuTime()));
842 //______________________________________________________________________________
844 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
845 //---------------------------------------------------------------------------//
846 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
847 //---------------------------------------------------------------------------//
848 if (fNumberOfTofClusters==kTofMaxCluster) {
849 AliError("Too many clusters !");
853 if (fNumberOfTofClusters==0) {
854 fTofClusters[fNumberOfTofClusters++] = tofCluster;
858 Int_t ii = FindClusterIndex(tofCluster->GetZ());
859 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
860 fTofClusters[ii] = tofCluster;
861 fNumberOfTofClusters++;
866 //_________________________________________________________________________
868 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
869 //--------------------------------------------------------------------
870 // This function returns the index of the nearest cluster
871 //--------------------------------------------------------------------
872 if (fNumberOfTofClusters==0) return 0;
873 if (z <= fTofClusters[0]->GetZ()) return 0;
874 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
875 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
876 for (; b<e; m=(b+e)/2) {
877 if (z > fTofClusters[m]->GetZ()) b=m+1;
884 //_________________________________________________________________________
886 void AliTOFClusterFinder::FillRecPoint()
889 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
890 // in Z) in the global TClonesArray of AliTOFcluster,
896 Int_t detectorIndex[5];
897 Double_t cylindricalPosition[3];
899 Int_t trackLabels[3];
900 Int_t digitIndex = -1;
903 TClonesArray &lRecPoints = *fRecPoints;
905 for (ii=0; ii<fNumberOfTofClusters; ii++) {
907 digitIndex = fTofClusters[ii]->GetIndex();
908 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
909 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
910 cylindricalPosition[0] = fTofClusters[ii]->GetR();
911 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
912 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
913 parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
914 parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
915 parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
916 parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
917 parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
918 status=fTofClusters[ii]->GetStatus();
919 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, detectorIndex, parTOF,status,trackLabels,digitIndex);
921 } // loop on clusters
925 //_________________________________________________________________________
926 void AliTOFClusterFinder::CalibrateRecPoint()
929 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
930 // in Z) in the global TClonesArray of AliTOFcluster,
936 Int_t detectorIndex[5];
937 Int_t digitIndex = -1;
941 AliInfo(" Calibrating TOF Clusters: ")
943 AliTOFCal *calTOFArray = fTOFcalib->GetTOFCalArray();
945 for (ii=0; ii<fNumberOfTofClusters; ii++) {
946 digitIndex = fTofClusters[ii]->GetIndex();
947 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
949 Int_t index = fTOFcalib->GetIndex(detectorIndex);
951 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
953 // Get channel status
954 Bool_t status=calChannel->GetStatus();
955 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
956 // Get Rough channel online equalization
957 Float_t roughDelay=calChannel->GetDelay();
958 AliDebug(2,Form(" channel delay = %f", roughDelay));
959 // Get Refined channel offline calibration parameters
961 for (Int_t j = 0; j<6; j++){
962 par[j]=calChannel->GetSlewPar(j);
963 AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
965 AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %i , %i", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
966 tToT = (Double_t)(fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()); tToT*=1.E-3; //ToT in ns
967 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
968 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+roughDelay; // the time correction
969 AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
970 timeCorr=(Double_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()*1.E-3-timeCorr;//redefine the time
972 AliDebug(2,Form(" The channel time, corr (ps)= %e",timeCorr ));
973 tdcCorr=(Int_t)(timeCorr/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
974 fTofClusters[ii]->SetTDC(tdcCorr);
975 AliDebug(2,Form(" The channel time, corr (counts) = %i",fTofClusters[ii]->GetTDC()));
977 } // loop on clusters
980 //______________________________________________________________________________
982 void AliTOFClusterFinder::ResetRecpoint()
985 // Clear the list of reconstructed points
988 fNumberOfTofClusters = 0;
989 if (fRecPoints) fRecPoints->Clear();
992 //______________________________________________________________________________
994 void AliTOFClusterFinder::Load()
997 // Load TOF.Digits.root and TOF.RecPoints.root files
1000 fTOFLoader->LoadDigits("READ");
1001 fTOFLoader->LoadRecPoints("recreate");
1004 //______________________________________________________________________________
1006 void AliTOFClusterFinder::LoadClusters()
1009 // Load TOF.RecPoints.root file
1012 fTOFLoader->LoadRecPoints("recreate");
1015 //______________________________________________________________________________
1017 void AliTOFClusterFinder::UnLoad()
1020 // Unload TOF.Digits.root and TOF.RecPoints.root files
1023 fTOFLoader->UnloadDigits();
1024 fTOFLoader->UnloadRecPoints();
1027 //______________________________________________________________________________
1029 void AliTOFClusterFinder::UnLoadClusters()
1032 // Unload TOF.RecPoints.root file
1035 fTOFLoader->UnloadRecPoints();
1038 //______________________________________________________________________________