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 **************************************************************************/
17 $Log: AliTOFClusterFinder.cxx,v $
18 Revision 1.31 2007/11/24 14:53:19 zampolli
19 Status flag implemented as UChar_t
21 Revision 1.30 2007/10/04 13:08:52 arcelli
22 updates to comply with AliTOFGeometryV5 becoming AliTOFGeometry
24 Revision 1.29 2007/10/03 10:42:33 arcelli
25 updates to handle new AliTOFcluster, inheriting form AliCluster3D
27 Revision 1.28 2007/05/31 16:06:05 arcelli
28 move instance of AliRawStream outside loop on DDL
30 Revision 1.27 2007/05/02 16:31:49 arcelli
31 Add methods to handle single event reconstruction. retrieval of Calib
32 info moved to AliTOFReconstructor ctor, and passed via a pointer to
35 Revision 1.26 2007/04/30 19:02:24 arcelli
36 hopefully the last refinements for correct type conversion in calibration
38 Revision 1.25 2007/04/30 15:22:17 arcelli
39 Change TOF digit Time, Tot etc to int type
41 Revision 1.24 2007/04/27 11:19:31 arcelli
42 updates for the new decoder
44 Revision 1.23 2007/04/23 16:51:39 decaro
45 Digits-to-raw_data conversion: correction for a more real description
46 (A.De Caro, R.Preghenella)
48 Revision 1.22 2007/04/19 17:26:32 arcelli
49 Fix a bug (add some debug printout
51 Revision 1.21 2007/04/18 17:28:12 arcelli
52 Set the ToT bin width to the one actually used...
54 Revision 1.20 2007/03/09 09:57:23 arcelli
55 Remove a forgotten include of Riostrem
57 Revision 1.19 2007/03/08 15:41:20 arcelli
58 set uncorrected times when filling RecPoints
60 Revision 1.18 2007/03/06 16:31:20 arcelli
61 Add Uncorrected TOF Time signal
63 Revision 1.17 2007/02/28 18:09:11 arcelli
64 Add protection against failed retrieval of the CDB cal object,
65 now Reconstruction exits with AliFatal
67 Revision 1.16 2007/02/20 15:57:00 decaro
68 Raw data update: to read the TOF raw data defined in UNPACKED mode
71 Revision 0.03 2005/07/28 A. De Caro
72 Implement public method
73 Raw2Digits(Int_t, AliRawReader *)
74 to convert digits from raw data in MC digits
77 Revision 0.02 2005/07/27 A. De Caro
78 Implement public method
79 Digits2RecPoint(Int_t)
80 to convert digits in clusters
82 Revision 0.02 2005/07/26 A. De Caro
83 Implement private methods
84 InsertCluster(AliTOFcluster *)
85 FindClusterIndex(Double_t)
86 originally implemented in AliTOFtracker
87 by S. Arcelli and C. Zampolli
89 Revision 0.01 2005/07/25 A. De Caro
90 Implement public methods
91 Digits2RecPoint(AliRawReader *, TTree *)
92 Digits2RecPoint(Int_t, AliRawReader *)
93 to convert raw data in clusters
96 ////////////////////////////////////////////////////////////////
98 // Class for TOF cluster finder //
100 // Starting from Raw Data, create rec points, //
101 // fill TreeR for TOF, //
102 // write TOF.RecPoints.root file //
104 ////////////////////////////////////////////////////////////////
106 #include "Riostream.h"
108 #include "TClonesArray.h"
109 #include "TStopwatch.h"
111 //#include <TGeoManager.h>
112 #include <TGeoMatrix.h>
113 //#include <TGeoPhysicalNode.h>
116 #include "AliLoader.h"
118 #include "AliRawReader.h"
119 #include "AliRunLoader.h"
120 //#include "AliAlignObj.h"
121 #include <AliGeomManager.h>
123 #include "AliTOFcalib.h"
124 #include "AliTOFChannelOnlineArray.h"
125 #include "AliTOFChannelOnlineStatusArray.h"
126 #include "AliTOFChannelOffline.h"
127 #include "AliTOFClusterFinder.h"
128 #include "AliTOFcluster.h"
129 #include "AliTOFdigit.h"
130 #include "AliTOFGeometry.h"
131 #include "AliTOFrawData.h"
132 #include "AliTOFRawStream.h"
134 //extern TFile *gFile;
136 ClassImp(AliTOFClusterFinder)
138 AliTOFClusterFinder::AliTOFClusterFinder(AliTOFcalib *calib):
143 fDigits(new TClonesArray("AliTOFdigit", 4000)),
144 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
145 fNumberOfTofClusters(0),
155 //______________________________________________________________________________
157 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader, AliTOFcalib *calib):
158 fRunLoader(runLoader),
159 fTOFLoader(runLoader->GetLoader("TOFLoader")),
162 fDigits(new TClonesArray("AliTOFdigit", 4000)),
163 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
164 fNumberOfTofClusters(0),
175 //------------------------------------------------------------------------
176 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
182 fDigits(new TClonesArray("AliTOFdigit", 4000)),
183 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
184 fNumberOfTofClusters(0),
190 this->fDigits=source.fDigits;
191 this->fRecPoints=source.fRecPoints;
192 this->fDecoderVersion=source.fDecoderVersion;
193 this->fTOFcalib=source.fTOFcalib;
197 //------------------------------------------------------------------------
198 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
201 this->fDigits=source.fDigits;
202 this->fRecPoints=source.fRecPoints;
203 this->fVerbose=source.fVerbose;
204 this->fDecoderVersion=source.fDecoderVersion;
205 this->fTOFcalib=source.fTOFcalib;
209 //______________________________________________________________________________
211 AliTOFClusterFinder::~AliTOFClusterFinder()
226 fRecPoints->Delete();
232 //______________________________________________________________________________
234 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
237 // Converts digits to recpoints for TOF
240 TStopwatch stopwatch;
243 fRunLoader->GetEvent(iEvent);
245 fTreeD = fTOFLoader->TreeD();
248 AliFatal("AliTOFClusterFinder: Can not get TreeD");
251 TBranch *branch = fTreeD->GetBranch("TOF");
253 AliError("can't get the branch with the TOF digits !");
257 TClonesArray staticdigits("AliTOFdigit",10000);
258 staticdigits.Clear();
259 TClonesArray *digits =&staticdigits;
260 branch->SetAddress(&digits);
264 fTreeR = fTOFLoader->TreeR();
267 fTOFLoader->MakeTree("R");
268 fTreeR = fTOFLoader->TreeR();
271 Int_t bufsize = 32000;
272 fTreeR->Branch("TOF", &fRecPoints, bufsize);
275 Int_t nDigits = digits->GetEntriesFast();
276 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
279 Int_t dig[5]; //cluster detector indeces
280 Int_t parTOF[5]; //The TOF signal parameters
281 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
282 for (ii=0; ii<nDigits; ii++) {
283 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
284 dig[0]=d->GetSector();
285 dig[1]=d->GetPlate();
286 dig[2]=d->GetStrip();
290 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
292 parTOF[0] = d->GetTdc(); //the TDC signal
293 parTOF[1] = d->GetToT(); //the ToT signal
294 parTOF[2] = d->GetAdc(); // the adc charge
295 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
296 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
299 UShort_t volIdClus=GetClusterVolIndex(dig);
300 GetClusterPars(dig, posClus,covClus);
301 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);
302 InsertCluster(tofCluster);
306 AliInfo(Form("Number of found clusters: %i for event: %i", fNumberOfTofClusters, iEvent));
314 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
315 fTOFLoader->WriteRecPoints("OVERWRITE");
317 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
318 stopwatch.RealTime(),stopwatch.CpuTime()));
322 //______________________________________________________________________________
324 void AliTOFClusterFinder::Digits2RecPoints(TTree* digitsTree, TTree* clusterTree)
327 // Converts digits to recpoints for TOF
330 TStopwatch stopwatch;
333 /// fRunLoader->GetEvent(iEvent);
335 if (digitsTree == 0x0)
337 AliFatal("AliTOFClusterFinder: Can not get TreeD");
340 TBranch *branch = digitsTree->GetBranch("TOF");
342 AliError("can't get the branch with the TOF digits !");
346 TClonesArray staticdigits("AliTOFdigit",10000);
347 staticdigits.Clear();
348 TClonesArray *digits = & staticdigits;
349 branch->SetAddress(&digits);
354 Int_t bufsize = 32000;
355 fTreeR->Branch("TOF", &fRecPoints, bufsize);
357 digitsTree->GetEvent(0);
358 Int_t nDigits = digits->GetEntriesFast();
359 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
362 Int_t dig[5]; //cluster detector indeces
363 Int_t parTOF[5]; //The TOF signal parameters
364 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
365 for (ii=0; ii<nDigits; ii++) {
366 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
367 dig[0]=d->GetSector();
368 dig[1]=d->GetPlate();
369 dig[2]=d->GetStrip();
373 // AliDebug(2,Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
375 parTOF[0] = d->GetTdc(); //the TDC signal
376 parTOF[1] = d->GetToT(); //the ToT signal
377 parTOF[2] = d->GetAdc(); // the adc charge
378 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
379 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
383 UShort_t volIdClus=GetClusterVolIndex(dig);
384 GetClusterPars(dig,posClus,covClus);
385 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);
386 InsertCluster(tofCluster);
390 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
398 AliInfo(Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
399 stopwatch.RealTime(),stopwatch.CpuTime()));
402 //______________________________________________________________________________
404 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
408 // Converts RAW data to recpoints for TOF
411 TStopwatch stopwatch;
414 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
415 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
419 Int_t bufsize = 32000;
420 clustersTree->Branch("TOF", &fRecPoints, bufsize);
422 TClonesArray * clonesRawData;
426 Int_t detectorIndex[5];
430 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
432 AliTOFRawStream tofInput(rawReader);
435 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
438 if (fDecoderVersion) {
439 AliInfo("Using New Decoder \n");
440 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
442 else tofInput.LoadRawData(indexDDL);
444 clonesRawData = (TClonesArray*)tofInput.GetRawData();
446 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
448 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
450 //if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
451 if (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 parTOF[0] = tofRawDatum->GetTOF(); //TDC
482 parTOF[1] = tofRawDatum->GetTOT(); // TOT
483 parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
484 parTOF[3] = -1;//raw data: no track of undecalib sim time
485 parTOF[4] = tofRawDatum->GetTOF(); // RAW time
488 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
489 Int_t lab[3]={-1,-1,-1};
491 GetClusterPars(detectorIndex,posClus,covClus);
492 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);
493 InsertCluster(tofCluster);
496 if (parTOF[1]<10)ftxt << " " << parTOF[1];
497 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
498 else ftxt << " " << parTOF[1];
499 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
500 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
501 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
502 else ftxt << " " << parTOF[3] << endl;
505 } // closed loop on TOF raw data per current DDL file
507 clonesRawData->Clear();
509 } // closed loop on DDL index
511 if (fVerbose==2) ftxt.close();
513 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
518 clustersTree->Fill();
522 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
523 stopwatch.RealTime(),stopwatch.CpuTime()));
526 //______________________________________________________________________________
528 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
531 // Converts RAW data to recpoints for TOF
534 TStopwatch stopwatch;
537 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
538 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
540 fRunLoader->GetEvent(iEvent);
542 AliDebug(2,Form(" Event number %2i ", iEvent));
544 fTreeR = fTOFLoader->TreeR();
547 fTOFLoader->MakeTree("R");
548 fTreeR = fTOFLoader->TreeR();
551 Int_t bufsize = 32000;
552 fTreeR->Branch("TOF", &fRecPoints, bufsize);
554 TClonesArray * clonesRawData;
558 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
561 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
563 AliTOFRawStream tofInput(rawReader);
566 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
569 if (fDecoderVersion) {
570 AliInfo("Using New Decoder \n");
571 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
573 else tofInput.LoadRawData(indexDDL);
575 clonesRawData = (TClonesArray*)tofInput.GetRawData();
577 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
579 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
581 //if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
582 if (tofRawDatum->GetTOF()==-1) continue;
585 if (indexDDL<10) ftxt << " " << indexDDL;
586 else ftxt << " " << indexDDL;
587 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
588 else ftxt << " " << tofRawDatum->GetTRM();
589 ftxt << " " << tofRawDatum->GetTRMchain();
590 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
591 else ftxt << " " << tofRawDatum->GetTDC();
592 ftxt << " " << tofRawDatum->GetTDCchannel();
595 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
596 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
597 dummy = detectorIndex[3];
598 detectorIndex[3] = detectorIndex[4];
599 detectorIndex[4] = dummy;
602 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
603 else ftxt << " -> " << detectorIndex[0];
604 ftxt << " " << detectorIndex[1];
605 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
606 else ftxt << " " << detectorIndex[2];
607 ftxt << " " << detectorIndex[3];
608 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
609 else ftxt << " " << detectorIndex[4];
612 parTOF[0] = tofRawDatum->GetTOF(); // TDC
613 parTOF[1] = tofRawDatum->GetTOT(); // TOT
614 parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
615 parTOF[3] = -1; //raw data: no track of the undecalib sim time
616 parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
619 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
620 Int_t lab[3]={-1,-1,-1};
622 GetClusterPars(detectorIndex,posClus,covClus);
623 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);
624 InsertCluster(tofCluster);
627 if (parTOF[1]<10)ftxt << " " << parTOF[1];
628 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
629 else ftxt << " " << parTOF[1];
630 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
631 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
632 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
633 else ftxt << " " << parTOF[3] << endl;
636 } // closed loop on TOF raw data per current DDL file
638 clonesRawData->Clear();
640 } // closed loop on DDL index
642 if (fVerbose==2) ftxt.close();
644 AliInfo(Form("Number of found clusters: %i for event: %i", fNumberOfTofClusters, iEvent));
652 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
653 fTOFLoader->WriteRecPoints("OVERWRITE");
655 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
656 stopwatch.RealTime(),stopwatch.CpuTime()));
659 //______________________________________________________________________________
661 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
664 // Converts RAW data to MC digits for TOF
666 // (temporary solution)
669 TStopwatch stopwatch;
672 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
673 const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
675 fRunLoader->GetEvent(iEvent);
677 fTreeD = fTOFLoader->TreeD();
680 AliInfo("TreeD re-creation");
682 fTOFLoader->MakeTree("D");
683 fTreeD = fTOFLoader->TreeD();
686 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
687 Int_t bufsize = 32000;
688 fTreeD->Branch("TOF", &tofDigits, bufsize);
690 fRunLoader->GetEvent(iEvent);
692 AliDebug(2,Form(" Event number %2i ", iEvent));
694 TClonesArray * clonesRawData;
698 Int_t detectorIndex[5];
701 AliTOFRawStream tofInput(rawReader);
704 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
707 if (fDecoderVersion) {
708 AliInfo("Using New Decoder \n");
709 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
711 else tofInput.LoadRawData(indexDDL);
713 clonesRawData = (TClonesArray*)tofInput.GetRawData();
715 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
717 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
719 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
721 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
722 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
723 dummy = detectorIndex[3];
724 detectorIndex[3] = detectorIndex[4];
725 detectorIndex[4] = dummy;
727 digit[0] = tofInput.GetTofBin();
728 digit[1] = tofInput.GetToTbin();
729 digit[2] = tofInput.GetToTbin();
732 Int_t tracknum[3]={-1,-1,-1};
734 TClonesArray &aDigits = *tofDigits;
735 Int_t last=tofDigits->GetEntriesFast();
736 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
740 clonesRawData->Clear();
746 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
747 fTOFLoader->WriteDigits("OVERWRITE");
749 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
750 stopwatch.RealTime(),stopwatch.CpuTime()));
754 //______________________________________________________________________________
756 void AliTOFClusterFinder::Raw2Digits(AliRawReader *rawReader, TTree* digitsTree)
759 // Converts RAW data to MC digits for TOF for the current event
763 TStopwatch stopwatch;
766 const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
770 AliError("No input digits Tree");
774 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
775 Int_t bufsize = 32000;
776 digitsTree->Branch("TOF", &tofDigits, bufsize);
778 /// fRunLoader->GetEvent(iEvent);
780 /// AliDebug(2,Form(" Event number %2i ", iEvent));
782 TClonesArray * clonesRawData;
786 Int_t detectorIndex[5];
789 AliTOFRawStream tofInput(rawReader);
792 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
795 if (fDecoderVersion) {
796 AliInfo("Using New Decoder \n");
797 tofInput.LoadRawDataBuffers(indexDDL,fVerbose);
799 else tofInput.LoadRawData(indexDDL);
801 clonesRawData = (TClonesArray*)tofInput.GetRawData();
803 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
805 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
807 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
809 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
810 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
811 dummy = detectorIndex[3];
812 detectorIndex[3] = detectorIndex[4];
813 detectorIndex[4] = dummy;
815 digit[0] = tofInput.GetTofBin();
816 digit[1] = tofInput.GetToTbin();
817 digit[2] = tofInput.GetToTbin();
820 Int_t tracknum[3]={-1,-1,-1};
822 TClonesArray &aDigits = *tofDigits;
823 Int_t last=tofDigits->GetEntriesFast();
824 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
828 clonesRawData->Clear();
834 AliDebug(1, Form("Got %d digits: ", tofDigits->GetEntries()));
835 AliDebug(1, Form("Execution time to read TOF raw data and fill TOF digit tree : R:%.2fs C:%.2fs",
836 stopwatch.RealTime(),stopwatch.CpuTime()));
839 //______________________________________________________________________________
841 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
842 //---------------------------------------------------------------------------//
843 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
844 //---------------------------------------------------------------------------//
845 if (fNumberOfTofClusters==kTofMaxCluster) {
846 AliError("Too many clusters !");
850 if (fNumberOfTofClusters==0) {
851 fTofClusters[fNumberOfTofClusters++] = tofCluster;
855 Int_t ii = FindClusterIndex(tofCluster->GetZ());
856 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
857 fTofClusters[ii] = tofCluster;
858 fNumberOfTofClusters++;
863 //_________________________________________________________________________
865 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
866 //--------------------------------------------------------------------
867 // This function returns the index of the nearest cluster in z
868 //--------------------------------------------------------------------
869 if (fNumberOfTofClusters==0) return 0;
870 if (z <= fTofClusters[0]->GetZ()) return 0;
871 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
872 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
873 for (; b<e; m=(b+e)/2) {
874 if (z > fTofClusters[m]->GetZ()) b=m+1;
881 //_________________________________________________________________________
883 void AliTOFClusterFinder::FillRecPoint()
886 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
887 // in Z) in the global TClonesArray of AliTOFcluster,
893 Int_t detectorIndex[5];
895 Int_t trackLabels[3];
896 Int_t digitIndex = -1;
899 TClonesArray &lRecPoints = *fRecPoints;
901 for (ii=0; ii<fNumberOfTofClusters; ii++) {
903 digitIndex = fTofClusters[ii]->GetIndex();
904 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
905 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
906 parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
907 parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
908 parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
909 parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
910 parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
911 status=fTofClusters[ii]->GetStatus();
914 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
915 GetClusterPars(detectorIndex,posClus,covClus);
916 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);
918 } // loop on clusters
922 //_________________________________________________________________________
923 void AliTOFClusterFinder::CalibrateRecPoint()
926 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
927 // in Z) in the global TClonesArray of AliTOFcluster,
933 Int_t detectorIndex[5];
934 Int_t digitIndex = -1;
938 AliInfo(" Calibrating TOF Clusters: ")
940 AliTOFChannelOnlineArray *calDelay = fTOFcalib->GetTOFOnlineDelay();
941 AliTOFChannelOnlineStatusArray *calStatus = fTOFcalib->GetTOFOnlineStatus();
942 TObjArray *calTOFArrayOffline = fTOFcalib->GetTOFCalArrayOffline();
943 TString validity = (TString)fTOFcalib->GetOfflineValidity();
944 AliInfo(Form(" validity = %s",validity.Data()));
945 Int_t calibration = -1;
946 if (validity.CompareTo("valid")==0) {
947 AliInfo(" Using offline calibration parameters");
951 AliInfo(" Using online calibration parameters");
954 for (ii=0; ii<fNumberOfTofClusters; ii++) {
955 digitIndex = fTofClusters[ii]->GetIndex();
956 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
958 Int_t index = AliTOFGeometry::GetIndex(detectorIndex);
960 UChar_t statusPulser=calStatus->GetPulserStatus(index);
961 UChar_t statusNoise=calStatus->GetNoiseStatus(index);
962 UChar_t statusHW=calStatus->GetHWStatus(index);
963 UChar_t status=calStatus->GetStatus(index);
965 //check the status, also unknown is fine!!!!!!!
967 AliDebug(2, Form(" Status for channel %i = %i",index, (Int_t)status));
968 if((statusPulser & AliTOFChannelOnlineStatusArray::kTOFPulserBad)==(AliTOFChannelOnlineStatusArray::kTOFPulserBad)||(statusNoise & AliTOFChannelOnlineStatusArray::kTOFNoiseBad)==(AliTOFChannelOnlineStatusArray::kTOFNoiseBad)||(statusHW & AliTOFChannelOnlineStatusArray::kTOFHWBad)==(AliTOFChannelOnlineStatusArray::kTOFHWBad)){
969 AliDebug(2, Form(" Bad Status for channel %i",index));
970 fTofClusters[ii]->SetStatus(kFALSE); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
973 AliDebug(2, Form(" Good Status for channel %i",index));
975 // Get Rough channel online equalization
976 Double_t roughDelay=(Double_t)calDelay->GetDelay(index); // in ns
977 AliDebug(2,Form(" channel delay (ns) = %f", roughDelay));
978 // Get Refined channel offline calibration parameters
979 if (calibration ==1){
980 AliTOFChannelOffline * calChannelOffline = (AliTOFChannelOffline*)calTOFArrayOffline->At(index);
982 for (Int_t j = 0; j<6; j++){
983 par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
985 AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
986 AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %i , %i", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
987 tToT = (Double_t)(fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth());
988 tToT*=1.E-3; //ToT in ns
989 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
990 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)
993 timeCorr = roughDelay; // correction in ns
995 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()));
996 AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
997 timeCorr=(Double_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()*1.E-3-timeCorr;//redefine the time
999 AliDebug(2,Form(" The channel time, corr (ps)= %e",timeCorr ));
1000 tdcCorr=(Int_t)(timeCorr/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
1001 fTofClusters[ii]->SetTDC(tdcCorr);
1002 } // loop on clusters
1005 //______________________________________________________________________________
1007 void AliTOFClusterFinder::ResetRecpoint()
1010 // Clear the list of reconstructed points
1013 fNumberOfTofClusters = 0;
1014 if (fRecPoints) fRecPoints->Clear();
1017 //______________________________________________________________________________
1019 void AliTOFClusterFinder::Load()
1022 // Load TOF.Digits.root and TOF.RecPoints.root files
1025 fTOFLoader->LoadDigits("READ");
1026 fTOFLoader->LoadRecPoints("recreate");
1029 //______________________________________________________________________________
1031 void AliTOFClusterFinder::LoadClusters()
1034 // Load TOF.RecPoints.root file
1037 fTOFLoader->LoadRecPoints("recreate");
1040 //______________________________________________________________________________
1042 void AliTOFClusterFinder::UnLoad()
1045 // Unload TOF.Digits.root and TOF.RecPoints.root files
1048 fTOFLoader->UnloadDigits();
1049 fTOFLoader->UnloadRecPoints();
1052 //______________________________________________________________________________
1054 void AliTOFClusterFinder::UnLoadClusters()
1057 // Unload TOF.RecPoints.root file
1060 fTOFLoader->UnloadRecPoints();
1063 //-------------------------------------------------------------------------
1064 UShort_t AliTOFClusterFinder::GetClusterVolIndex(Int_t *ind) const {
1066 //First of all get the volume ID to retrieve the l2t transformation...
1068 // Detector numbering scheme
1075 Int_t isector =ind[0];
1076 if (isector >= nSector)
1077 AliError(Form("Wrong sector number in TOF (%d) !",isector));
1078 Int_t iplate = ind[1];
1079 if (iplate >= nPlate)
1080 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1081 Int_t istrip = ind[2];
1083 Int_t stripOffset = 0;
1089 stripOffset = nStripC;
1092 stripOffset = nStripC+nStripB;
1095 stripOffset = nStripC+nStripB+nStripA;
1098 stripOffset = nStripC+nStripB+nStripA+nStripB;
1101 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1105 Int_t index= (2*(nStripC+nStripB)+nStripA)*isector +
1109 UShort_t volIndex = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,index);
1113 //-------------------------------------------------------------------------
1114 void AliTOFClusterFinder::GetClusterPars(Int_t *ind, Double_t* pos,Double_t* cov) const {
1116 //First of all get the volume ID to retrieve the l2t transformation...
1118 UShort_t volIndex = GetClusterVolIndex(ind);
1121 //we now go in the system of the strip: determine the local coordinates
1124 // 47---------------------------------------------------0 ^ z
1125 // | | | | | | | | | | | | | | | | | | | | | | | | | | | 1 |
1126 // ----------------------------------------------------- | y going outwards
1127 // | | | | | | | | | | | | | | | | | | | | | | | | | | | 0 | par[0]=0;
1129 // ----------------------------------------------------- |
1130 // x <-----------------------------------------------------
1132 Float_t localX=(ind[4]-23.5)*2.5;
1134 Float_t localZ=(ind[3]-0.5)*3.5;
1136 //move to the tracking ref system
1143 const TGeoHMatrix *l2t= AliGeomManager::GetTracking2LocalMatrix(volIndex);
1144 // Get The position in the track ref system
1146 l2t->MasterToLocal(lpos,tpos);
1151 //Get The cluster covariance in the track ref system
1154 //cluster covariance in the local system:
1159 lcov[0]=2.5*2.5/12.;
1167 lcov[8]=3.5*3.5/12.;
1169 //cluster covariance in the tracking system:
1171 m.SetRotation(lcov);
1173 m.MultiplyLeft(&l2t->Inverse());
1174 Double_t *tcov = m.GetRotationMatrix();
1175 cov[0] = tcov[0]; cov[1] = tcov[1]; cov[2] = tcov[2];
1176 cov[3] = tcov[4]; cov[4] = tcov[5];