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.29 2007/10/03 10:42:33 arcelli
19 updates to handle new AliTOFcluster, inheriting form AliCluster3D
21 Revision 1.28 2007/05/31 16:06:05 arcelli
22 move instance of AliRawStream outside loop on DDL
24 Revision 1.27 2007/05/02 16:31:49 arcelli
25 Add methods to handle single event reconstruction. retrieval of Calib info moved to AliTOFReconstructor ctor, and passed via a pointer to AliTOFcalib
27 Revision 1.26 2007/04/30 19:02:24 arcelli
28 hopefully the last refinements for correct type conversion in calibration
30 Revision 1.25 2007/04/30 15:22:17 arcelli
31 Change TOF digit Time, Tot etc to int type
33 Revision 1.24 2007/04/27 11:19:31 arcelli
34 updates for the new decoder
36 Revision 1.23 2007/04/23 16:51:39 decaro
37 Digits-to-raw_data conversion: correction for a more real description (A.De Caro, R.Preghenella)
39 Revision 1.22 2007/04/19 17:26:32 arcelli
40 Fix a bug (add some debug printout
42 Revision 1.21 2007/04/18 17:28:12 arcelli
43 Set the ToT bin width to the one actually used...
45 Revision 1.20 2007/03/09 09:57:23 arcelli
46 Remove a forgotten include of Riostrem
48 Revision 1.19 2007/03/08 15:41:20 arcelli
49 set uncorrected times when filling RecPoints
51 Revision 1.18 2007/03/06 16:31:20 arcelli
52 Add Uncorrected TOF Time signal
54 Revision 1.17 2007/02/28 18:09:11 arcelli
55 Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
57 Revision 1.16 2007/02/20 15:57:00 decaro
58 Raw data update: to read the TOF raw data defined in UNPACKED mode
61 Revision 0.03 2005/07/28 A. De Caro
62 Implement public method
63 Raw2Digits(Int_t, AliRawReader *)
64 to convert digits from raw data in MC digits
67 Revision 0.02 2005/07/27 A. De Caro
68 Implement public method
69 Digits2RecPoint(Int_t)
70 to convert digits in clusters
72 Revision 0.02 2005/07/26 A. De Caro
73 Implement private methods
74 InsertCluster(AliTOFcluster *)
75 FindClusterIndex(Double_t)
76 originally implemented in AliTOFtracker
77 by S. Arcelli and C. Zampolli
79 Revision 0.01 2005/07/25 A. De Caro
80 Implement public methods
81 Digits2RecPoint(AliRawReader *, TTree *)
82 Digits2RecPoint(Int_t, AliRawReader *)
83 to convert raw data in clusters
86 ////////////////////////////////////////////////////////////////
88 // Class for TOF cluster finder //
90 // Starting from Raw Data, create rec points, //
91 // fill TreeR for TOF, //
92 // write TOF.RecPoints.root file //
94 ////////////////////////////////////////////////////////////////
97 #include "TClonesArray.h"
98 #include "TStopwatch.h"
100 #include <TGeoManager.h>
101 #include <AliGeomManager.h>
102 #include <TGeoMatrix.h>
103 #include <TGeoPhysicalNode.h>
104 #include "AliAlignObj.h"
107 #include "AliLoader.h"
109 #include "AliRawReader.h"
110 #include "AliRunLoader.h"
112 #include "AliTOFcalib.h"
113 #include "AliTOFChannelOnline.h"
114 #include "AliTOFChannelOffline.h"
115 #include "AliTOFClusterFinder.h"
116 #include "AliTOFcluster.h"
117 #include "AliTOFdigit.h"
118 #include "AliTOFGeometry.h"
119 #include "AliTOFrawData.h"
120 #include "AliTOFRawStream.h"
121 #include "Riostream.h"
123 //extern TFile *gFile;
125 ClassImp(AliTOFClusterFinder)
127 AliTOFClusterFinder::AliTOFClusterFinder(AliTOFcalib *calib):
132 fDigits(new TClonesArray("AliTOFdigit", 4000)),
133 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
134 fNumberOfTofClusters(0),
144 //______________________________________________________________________________
146 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader, AliTOFcalib *calib):
147 fRunLoader(runLoader),
148 fTOFLoader(runLoader->GetLoader("TOFLoader")),
151 fDigits(new TClonesArray("AliTOFdigit", 4000)),
152 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
153 fNumberOfTofClusters(0),
164 //------------------------------------------------------------------------
165 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
171 fDigits(new TClonesArray("AliTOFdigit", 4000)),
172 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
173 fNumberOfTofClusters(0),
179 this->fDigits=source.fDigits;
180 this->fRecPoints=source.fRecPoints;
181 this->fDecoderVersion=source.fDecoderVersion;
182 this->fTOFcalib=source.fTOFcalib;
186 //------------------------------------------------------------------------
187 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
190 this->fDigits=source.fDigits;
191 this->fRecPoints=source.fRecPoints;
192 this->fVerbose=source.fVerbose;
193 this->fDecoderVersion=source.fDecoderVersion;
194 this->fTOFcalib=source.fTOFcalib;
198 //______________________________________________________________________________
200 AliTOFClusterFinder::~AliTOFClusterFinder()
215 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 Int_t parTOF[5]; //The TOF signal parameters
268 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
269 for (ii=0; ii<nDigits; ii++) {
270 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
271 dig[0]=d->GetSector();
272 dig[1]=d->GetPlate();
273 dig[2]=d->GetStrip();
277 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
279 parTOF[0] = d->GetTdc(); //the TDC signal
280 parTOF[1] = d->GetToT(); //the ToT signal
281 parTOF[2] = d->GetAdc(); // the adc charge
282 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
283 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
286 UShort_t volIdClus=GetClusterVolIndex(dig);
287 GetClusterPars(dig, posClus,covClus);
288 AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],d->GetTracks(),dig,parTOF,status,ii);
289 InsertCluster(tofCluster);
293 AliInfo(Form("Number of found clusters: %i for event: %i", fNumberOfTofClusters, iEvent));
301 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
302 fTOFLoader->WriteRecPoints("OVERWRITE");
304 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
305 stopwatch.RealTime(),stopwatch.CpuTime()));
311 //______________________________________________________________________________
313 void AliTOFClusterFinder::Digits2RecPoints(TTree* digitsTree, TTree* clusterTree)
316 // Converts digits to recpoints for TOF
319 TStopwatch stopwatch;
322 /// fRunLoader->GetEvent(iEvent);
324 if (digitsTree == 0x0)
326 AliFatal("AliTOFClusterFinder: Can not get TreeD");
329 TBranch *branch = digitsTree->GetBranch("TOF");
331 AliError("can't get the branch with the TOF digits !");
335 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
336 branch->SetAddress(&digits);
341 Int_t bufsize = 32000;
342 fTreeR->Branch("TOF", &fRecPoints, bufsize);
344 digitsTree->GetEvent(0);
345 Int_t nDigits = digits->GetEntriesFast();
346 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
349 Int_t dig[5]; //cluster detector indeces
350 Int_t parTOF[5]; //The TOF signal parameters
351 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
352 for (ii=0; ii<nDigits; ii++) {
353 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
354 dig[0]=d->GetSector();
355 dig[1]=d->GetPlate();
356 dig[2]=d->GetStrip();
360 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
362 parTOF[0] = d->GetTdc(); //the TDC signal
363 parTOF[1] = d->GetToT(); //the ToT signal
364 parTOF[2] = d->GetAdc(); // the adc charge
365 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
366 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
370 UShort_t volIdClus=GetClusterVolIndex(dig);
371 GetClusterPars(dig,posClus,covClus);
372 AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],d->GetTracks(),dig,parTOF,status,ii);
373 InsertCluster(tofCluster);
377 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
385 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
386 stopwatch.RealTime(),stopwatch.CpuTime()));
391 //______________________________________________________________________________
393 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
397 // Converts RAW data to recpoints for TOF
400 TStopwatch stopwatch;
403 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
404 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
408 Int_t bufsize = 32000;
409 clustersTree->Branch("TOF", &fRecPoints, bufsize);
411 TClonesArray * clonesRawData;
415 Int_t detectorIndex[5];
419 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
421 AliTOFRawStream tofInput(rawReader);
424 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
427 if (fDecoderVersion) {
428 AliInfo("Using New Decoder \n");
429 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
431 else tofInput.LoadRawData(indexDDL);
433 clonesRawData = (TClonesArray*)tofInput.GetRawData();
435 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
437 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
439 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
442 if (indexDDL<10) ftxt << " " << indexDDL;
443 else ftxt << " " << indexDDL;
444 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
445 else ftxt << " " << tofRawDatum->GetTRM();
446 ftxt << " " << tofRawDatum->GetTRMchain();
447 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
448 else ftxt << " " << tofRawDatum->GetTDC();
449 ftxt << " " << tofRawDatum->GetTDCchannel();
452 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
453 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
454 dummy = detectorIndex[3];
455 detectorIndex[3] = detectorIndex[4];
456 detectorIndex[4] = dummy;
459 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
460 else ftxt << " -> " << detectorIndex[0];
461 ftxt << " " << detectorIndex[1];
462 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
463 else ftxt << " " << detectorIndex[2];
464 ftxt << " " << detectorIndex[3];
465 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
466 else ftxt << " " << detectorIndex[4];
469 parTOF[0] = tofRawDatum->GetTOF(); //TDC
470 parTOF[1] = tofRawDatum->GetTOT(); // TOT
471 parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
472 parTOF[3] = -1;//raw data: no track of undecalib sim time
473 parTOF[4] = tofRawDatum->GetTOF(); // RAW time
476 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
477 Int_t lab[3]={-1,-1,-1};
479 GetClusterPars(detectorIndex,posClus,covClus);
480 AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],lab,detectorIndex,parTOF,status,-1);
481 InsertCluster(tofCluster);
484 if (parTOF[1]<10)ftxt << " " << parTOF[1];
485 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
486 else ftxt << " " << parTOF[1];
487 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
488 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
489 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
490 else ftxt << " " << parTOF[3] << endl;
493 } // closed loop on TOF raw data per current DDL file
495 clonesRawData->Clear();
497 } // closed loop on DDL index
499 if (fVerbose==2) ftxt.close();
501 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
506 clustersTree->Fill();
510 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
511 stopwatch.RealTime(),stopwatch.CpuTime()));
514 //______________________________________________________________________________
516 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
519 // Converts RAW data to recpoints for TOF
522 TStopwatch stopwatch;
525 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
526 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
528 fRunLoader->GetEvent(iEvent);
530 AliDebug(2,Form(" Event number %2i ", iEvent));
532 fTreeR = fTOFLoader->TreeR();
535 fTOFLoader->MakeTree("R");
536 fTreeR = fTOFLoader->TreeR();
539 Int_t bufsize = 32000;
540 fTreeR->Branch("TOF", &fRecPoints, bufsize);
542 TClonesArray * clonesRawData;
546 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
549 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
551 AliTOFRawStream tofInput(rawReader);
554 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
557 if (fDecoderVersion) {
558 AliInfo("Using New Decoder \n");
559 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
561 else tofInput.LoadRawData(indexDDL);
563 clonesRawData = (TClonesArray*)tofInput.GetRawData();
565 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
567 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
569 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
572 if (indexDDL<10) ftxt << " " << indexDDL;
573 else ftxt << " " << indexDDL;
574 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
575 else ftxt << " " << tofRawDatum->GetTRM();
576 ftxt << " " << tofRawDatum->GetTRMchain();
577 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
578 else ftxt << " " << tofRawDatum->GetTDC();
579 ftxt << " " << tofRawDatum->GetTDCchannel();
582 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
583 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
584 dummy = detectorIndex[3];
585 detectorIndex[3] = detectorIndex[4];
586 detectorIndex[4] = dummy;
589 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
590 else ftxt << " -> " << detectorIndex[0];
591 ftxt << " " << detectorIndex[1];
592 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
593 else ftxt << " " << detectorIndex[2];
594 ftxt << " " << detectorIndex[3];
595 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
596 else ftxt << " " << detectorIndex[4];
599 parTOF[0] = tofRawDatum->GetTOF(); // TDC
600 parTOF[1] = tofRawDatum->GetTOT(); // TOT
601 parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
602 parTOF[3] = -1; //raw data: no track of the undecalib sim time
603 parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
606 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
607 Int_t lab[3]={-1,-1,-1};
609 GetClusterPars(detectorIndex,posClus,covClus);
610 AliTOFcluster *tofCluster = new AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],lab,detectorIndex,parTOF,status,-1);
611 InsertCluster(tofCluster);
614 if (parTOF[1]<10)ftxt << " " << parTOF[1];
615 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
616 else ftxt << " " << parTOF[1];
617 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
618 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
619 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
620 else ftxt << " " << parTOF[3] << endl;
623 } // closed loop on TOF raw data per current DDL file
625 clonesRawData->Clear();
627 } // closed loop on DDL index
629 if (fVerbose==2) ftxt.close();
631 AliInfo(Form("Number of found clusters: %i for event: %i", fNumberOfTofClusters, iEvent));
639 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
640 fTOFLoader->WriteRecPoints("OVERWRITE");
642 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
643 stopwatch.RealTime(),stopwatch.CpuTime()));
646 //______________________________________________________________________________
648 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
651 // Converts RAW data to MC digits for TOF
653 // (temporary solution)
656 TStopwatch stopwatch;
659 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
660 const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
662 fRunLoader->GetEvent(iEvent);
664 fTreeD = fTOFLoader->TreeD();
667 AliInfo("TreeD re-creation");
669 fTOFLoader->MakeTree("D");
670 fTreeD = fTOFLoader->TreeD();
673 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
674 Int_t bufsize = 32000;
675 fTreeD->Branch("TOF", &tofDigits, bufsize);
677 fRunLoader->GetEvent(iEvent);
679 AliDebug(2,Form(" Event number %2i ", iEvent));
681 TClonesArray * clonesRawData;
685 Int_t detectorIndex[5];
688 AliTOFRawStream tofInput(rawReader);
691 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
694 if (fDecoderVersion) {
695 AliInfo("Using New Decoder \n");
696 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
698 else tofInput.LoadRawData(indexDDL);
700 clonesRawData = (TClonesArray*)tofInput.GetRawData();
702 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
704 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
706 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
708 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
709 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
710 dummy = detectorIndex[3];
711 detectorIndex[3] = detectorIndex[4];
712 detectorIndex[4] = dummy;
714 digit[0] = tofInput.GetTofBin();
715 digit[1] = tofInput.GetToTbin();
716 digit[2] = tofInput.GetToTbin();
719 Int_t tracknum[3]={-1,-1,-1};
721 TClonesArray &aDigits = *tofDigits;
722 Int_t last=tofDigits->GetEntriesFast();
723 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
727 clonesRawData->Clear();
733 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
734 fTOFLoader->WriteDigits("OVERWRITE");
736 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
737 stopwatch.RealTime(),stopwatch.CpuTime()));
741 //______________________________________________________________________________
743 void AliTOFClusterFinder::Raw2Digits(AliRawReader *rawReader, TTree* digitsTree)
746 // Converts RAW data to MC digits for TOF for the current event
750 TStopwatch stopwatch;
753 const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
757 AliError("No input digits Tree");
761 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
762 Int_t bufsize = 32000;
763 digitsTree->Branch("TOF", &tofDigits, bufsize);
765 /// fRunLoader->GetEvent(iEvent);
767 /// AliDebug(2,Form(" Event number %2i ", iEvent));
769 TClonesArray * clonesRawData;
773 Int_t detectorIndex[5];
776 AliTOFRawStream tofInput(rawReader);
779 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
782 if (fDecoderVersion) {
783 AliInfo("Using New Decoder \n");
784 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
786 else tofInput.LoadRawData(indexDDL);
788 clonesRawData = (TClonesArray*)tofInput.GetRawData();
790 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
792 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
794 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
796 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
797 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
798 dummy = detectorIndex[3];
799 detectorIndex[3] = detectorIndex[4];
800 detectorIndex[4] = dummy;
802 digit[0] = tofInput.GetTofBin();
803 digit[1] = tofInput.GetToTbin();
804 digit[2] = tofInput.GetToTbin();
807 Int_t tracknum[3]={-1,-1,-1};
809 TClonesArray &aDigits = *tofDigits;
810 Int_t last=tofDigits->GetEntriesFast();
811 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
815 clonesRawData->Clear();
821 AliDebug(1, Form("Got %d digits: ", tofDigits->GetEntries()));
822 AliDebug(1, Form("Execution time to read TOF raw data and fill TOF digit tree : R:%.2fs C:%.2fs",
823 stopwatch.RealTime(),stopwatch.CpuTime()));
826 //______________________________________________________________________________
828 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
829 //---------------------------------------------------------------------------//
830 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
831 //---------------------------------------------------------------------------//
832 if (fNumberOfTofClusters==kTofMaxCluster) {
833 AliError("Too many clusters !");
837 if (fNumberOfTofClusters==0) {
838 fTofClusters[fNumberOfTofClusters++] = tofCluster;
842 Int_t ii = FindClusterIndex(tofCluster->GetZ());
843 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
844 fTofClusters[ii] = tofCluster;
845 fNumberOfTofClusters++;
850 //_________________________________________________________________________
852 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
853 //--------------------------------------------------------------------
854 // This function returns the index of the nearest cluster in z
855 //--------------------------------------------------------------------
856 if (fNumberOfTofClusters==0) return 0;
857 if (z <= fTofClusters[0]->GetZ()) return 0;
858 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
859 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
860 for (; b<e; m=(b+e)/2) {
861 if (z > fTofClusters[m]->GetZ()) b=m+1;
868 //_________________________________________________________________________
870 void AliTOFClusterFinder::FillRecPoint()
873 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
874 // in Z) in the global TClonesArray of AliTOFcluster,
880 Int_t detectorIndex[5];
882 Int_t trackLabels[3];
883 Int_t digitIndex = -1;
886 TClonesArray &lRecPoints = *fRecPoints;
888 for (ii=0; ii<fNumberOfTofClusters; ii++) {
890 digitIndex = fTofClusters[ii]->GetIndex();
891 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
892 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
893 parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
894 parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
895 parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
896 parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
897 parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
898 status=fTofClusters[ii]->GetStatus();
901 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
902 GetClusterPars(detectorIndex,posClus,covClus);
903 new(lRecPoints[ii]) AliTOFcluster(volIdClus,posClus[0],posClus[1],posClus[2],covClus[0],covClus[1],covClus[2],covClus[3],covClus[4],covClus[5],trackLabels,detectorIndex, parTOF,status,digitIndex);
905 } // loop on clusters
909 //_________________________________________________________________________
910 void AliTOFClusterFinder::CalibrateRecPoint()
913 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
914 // in Z) in the global TClonesArray of AliTOFcluster,
920 Int_t detectorIndex[5];
921 Int_t digitIndex = -1;
925 AliInfo(" Calibrating TOF Clusters: ")
927 TObjArray *calTOFArrayOnline = fTOFcalib->GetTOFCalArrayOnline();
928 TObjArray *calTOFArrayOffline = fTOFcalib->GetTOFCalArrayOffline();
929 TString validity = (TString)fTOFcalib->GetOfflineValidity();
930 AliInfo(Form(" validity = %s",validity.Data()));
931 Int_t calibration = -1;
932 if (validity.CompareTo("valid")==0) {
933 AliInfo(" Using offline calibration parameters");
937 AliInfo(" Using online calibration parameters");
940 for (ii=0; ii<fNumberOfTofClusters; ii++) {
941 digitIndex = fTofClusters[ii]->GetIndex();
942 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
944 Int_t index = AliTOFGeometry::GetIndex(detectorIndex);
946 AliTOFChannelOnline * calChannelOnline = (AliTOFChannelOnline* )calTOFArrayOnline->At(index);
948 // Get channel status
949 Bool_t status=calChannelOnline->GetStatus();
950 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
951 // Get Rough channel online equalization
952 Double_t roughDelay=(Double_t)calChannelOnline->GetDelay(); // in ns
953 AliDebug(2,Form(" channel delay (ns) = %f", roughDelay));
954 // Get Refined channel offline calibration parameters
955 if (calibration ==1){
956 AliTOFChannelOffline * calChannelOffline = (AliTOFChannelOffline*)calTOFArrayOffline->At(index);
958 for (Int_t j = 0; j<6; j++){
959 par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
961 AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
962 AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %i , %i", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
963 tToT = (Double_t)(fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth());
964 tToT*=1.E-3; //ToT in ns
965 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
966 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; // the time correction (ns)
969 timeCorr = roughDelay; // correction in ns
971 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()));
972 AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
973 timeCorr=(Double_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()*1.E-3-timeCorr;//redefine the time
975 AliDebug(2,Form(" The channel time, corr (ps)= %e",timeCorr ));
976 tdcCorr=(Int_t)(timeCorr/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
977 fTofClusters[ii]->SetTDC(tdcCorr);
978 } // loop on clusters
981 //______________________________________________________________________________
983 void AliTOFClusterFinder::ResetRecpoint()
986 // Clear the list of reconstructed points
989 fNumberOfTofClusters = 0;
990 if (fRecPoints) fRecPoints->Clear();
993 //______________________________________________________________________________
995 void AliTOFClusterFinder::Load()
998 // Load TOF.Digits.root and TOF.RecPoints.root files
1001 fTOFLoader->LoadDigits("READ");
1002 fTOFLoader->LoadRecPoints("recreate");
1005 //______________________________________________________________________________
1007 void AliTOFClusterFinder::LoadClusters()
1010 // Load TOF.RecPoints.root file
1013 fTOFLoader->LoadRecPoints("recreate");
1016 //______________________________________________________________________________
1018 void AliTOFClusterFinder::UnLoad()
1021 // Unload TOF.Digits.root and TOF.RecPoints.root files
1024 fTOFLoader->UnloadDigits();
1025 fTOFLoader->UnloadRecPoints();
1028 //______________________________________________________________________________
1030 void AliTOFClusterFinder::UnLoadClusters()
1033 // Unload TOF.RecPoints.root file
1036 fTOFLoader->UnloadRecPoints();
1039 //-------------------------------------------------------------------------
1040 UShort_t AliTOFClusterFinder::GetClusterVolIndex(Int_t *ind) const {
1042 //First of all get the volume ID to retrieve the l2t transformation...
1044 // Detector numbering scheme
1051 Int_t isector =ind[0];
1052 if (isector >= nSector)
1053 AliError(Form("Wrong sector number in TOF (%d) !",isector));
1054 Int_t iplate = ind[1];
1055 if (iplate >= nPlate)
1056 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1057 Int_t istrip = ind[2];
1059 Int_t stripOffset = 0;
1065 stripOffset = nStripC;
1068 stripOffset = nStripC+nStripB;
1071 stripOffset = nStripC+nStripB+nStripA;
1074 stripOffset = nStripC+nStripB+nStripA+nStripB;
1077 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1081 Int_t index= (2*(nStripC+nStripB)+nStripA)*isector +
1085 UShort_t volIndex = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,index);
1089 //-------------------------------------------------------------------------
1090 void AliTOFClusterFinder::GetClusterPars(Int_t *ind, Double_t* pos,Double_t* cov) const {
1092 //First of all get the volume ID to retrieve the l2t transformation...
1094 UShort_t volIndex = GetClusterVolIndex(ind);
1097 //we now go in the system of the strip: determine the local coordinates
1100 // 47---------------------------------------------------0 ^ z
1101 // | | | | | | | | | | | | | | | | | | | | | | | | | | | 1 |
1102 // ----------------------------------------------------- | y going outwards
1103 // | | | | | | | | | | | | | | | | | | | | | | | | | | | 0 | par[0]=0;
1105 // ----------------------------------------------------- |
1106 // x <-----------------------------------------------------
1108 Float_t localX=(ind[4]-23.5)*2.5;
1110 Float_t localZ=(ind[3]-0.5)*3.5;
1112 //move to the tracking ref system
1119 const TGeoHMatrix *l2t= AliGeomManager::GetTracking2LocalMatrix(volIndex);
1120 // Get The position in the track ref system
1122 l2t->MasterToLocal(lpos,tpos);
1127 //Get The cluster covariance in the track ref system
1130 //cluster covariance in the local system:
1135 lcov[0]=2.5*2.5/12.;
1143 lcov[8]=3.5*3.5/12.;
1145 //cluster covariance in the tracking system:
1147 m.SetRotation(lcov);
1149 m.MultiplyLeft(&l2t->Inverse());
1150 Double_t *tcov = m.GetRotationMatrix();
1151 cov[0] = tcov[0]; cov[1] = tcov[1]; cov[2] = tcov[2];
1152 cov[3] = tcov[4]; cov[4] = tcov[5];