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 #include "AliTOFDeltaBCOffset.h"
135 #include "AliTOFCTPLatency.h"
136 #include "AliTOFRunParams.h"
138 //extern TFile *gFile;
140 ClassImp(AliTOFClusterFinder)
142 AliTOFClusterFinder::AliTOFClusterFinder(AliTOFcalib *calib):
147 fDigits(new TClonesArray("AliTOFdigit", 4000)),
148 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
149 fNumberOfTofClusters(0),
153 fTOFRawStream(AliTOFRawStream())
159 TString validity = (TString)fTOFcalib->GetOfflineValidity();
160 if (validity.CompareTo("valid")==0) {
161 AliInfo(Form(" validity = %s - Using offline calibration parameters",validity.Data()));
163 AliInfo(Form(" validity = %s - Using online calibration parameters",validity.Data()));
168 //______________________________________________________________________________
170 AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader, AliTOFcalib *calib):
171 fRunLoader(runLoader),
172 fTOFLoader(runLoader->GetLoader("TOFLoader")),
175 fDigits(new TClonesArray("AliTOFdigit", 4000)),
176 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
177 fNumberOfTofClusters(0),
181 fTOFRawStream(AliTOFRawStream())
187 TString validity = (TString)fTOFcalib->GetOfflineValidity();
188 if (validity.CompareTo("valid")==0) {
189 AliInfo(Form(" validity = %s - Using offline calibration parameters",validity.Data()));
191 AliInfo(Form(" validity = %s - Using online calibration parameters",validity.Data()));
196 //------------------------------------------------------------------------
197 AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
203 fDigits(source.fDigits),
204 fRecPoints(source.fRecPoints),
205 fNumberOfTofClusters(0),
207 fDecoderVersion(source.fDecoderVersion),
208 fTOFcalib(source.fTOFcalib),
209 fTOFRawStream(source.fTOFRawStream)
214 //------------------------------------------------------------------------
215 AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
222 TObject::operator=(source);
223 fDigits=source.fDigits;
224 fRecPoints=source.fRecPoints;
225 fVerbose=source.fVerbose;
226 fDecoderVersion=source.fDecoderVersion;
227 fTOFcalib=source.fTOFcalib;
228 fTOFRawStream=source.fTOFRawStream;
232 //______________________________________________________________________________
234 AliTOFClusterFinder::~AliTOFClusterFinder()
249 fRecPoints->Delete();
254 if (fTofClusters || fNumberOfTofClusters) {
255 for (Int_t ii=0; ii<fNumberOfTofClusters; ii++)
256 if (fTofClusters[ii]) fTofClusters[ii]->Delete();
257 fNumberOfTofClusters=0;
261 //______________________________________________________________________________
263 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
266 // Converts digits to recpoints for TOF
269 TStopwatch stopwatch;
274 fRunLoader->GetEvent(iEvent);
276 fTreeD = fTOFLoader->TreeD();
279 AliFatal("AliTOFClusterFinder: Can not get TreeD");
282 TBranch *branch = fTreeD->GetBranch("TOF");
284 AliError("can't get the branch with the TOF digits !");
288 TClonesArray staticdigits("AliTOFdigit",10000);
289 staticdigits.Clear();
290 TClonesArray *digits =&staticdigits;
291 branch->SetAddress(&digits);
295 fTreeR = fTOFLoader->TreeR();
298 fTOFLoader->MakeTree("R");
299 fTreeR = fTOFLoader->TreeR();
302 Int_t bufsize = 32000;
303 fTreeR->Branch("TOF", &fRecPoints, bufsize);
306 Int_t nDigits = digits->GetEntriesFast();
307 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
310 Int_t dig[5]; //cluster detector indeces
311 Int_t parTOF[7]; //The TOF signal parameters
312 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
313 for (ii=0; ii<nDigits; ii++) {
314 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
315 dig[0]=d->GetSector();
316 dig[1]=d->GetPlate();
317 dig[2]=d->GetStrip();
321 /* check valid index */
322 if (dig[0]==-1||dig[1]==-1||dig[2]==-1||dig[3]==-1||dig[4]==-1) continue;
324 // Do not reconstruct anything in the holes
325 if (dig[0]==13 || dig[0]==14 || dig[0]==15 ) { // sectors with holes
326 if (dig[1]==2) { // plate with holes
332 AliDebug(2,Form(" %2d %1d %2d %1d %2d ",dig[0],dig[1],dig[2],dig[3],dig[4]));
334 parTOF[0] = d->GetTdc(); //the TDC signal
335 parTOF[1] = d->GetToT(); //the ToT signal
336 parTOF[2] = d->GetAdc(); // the adc charge
337 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
338 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
339 parTOF[5] = 0; // deltaBC
340 parTOF[6] = 0; // L0-L1 latency
343 UShort_t volIdClus=GetClusterVolIndex(dig);
344 GetClusterPars(dig, posClus,covClus);
345 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);
346 InsertCluster(tofCluster);
350 AliDebug(1,Form("Number of found clusters: %d for event: %d", fNumberOfTofClusters, iEvent));
358 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
359 fTOFLoader->WriteRecPoints("OVERWRITE");
361 AliDebug(1,Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
362 stopwatch.RealTime(),stopwatch.CpuTime()));
363 if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
367 //______________________________________________________________________________
369 void AliTOFClusterFinder::Digits2RecPoints(TTree* digitsTree, TTree* clusterTree)
372 // Converts digits to recpoints for TOF
375 TStopwatch stopwatch;
380 /// fRunLoader->GetEvent(iEvent);
382 if (digitsTree == 0x0)
384 AliFatal("AliTOFClusterFinder: Can not get TreeD");
387 TBranch *branch = digitsTree->GetBranch("TOF");
389 AliError("can't get the branch with the TOF digits !");
393 TClonesArray staticdigits("AliTOFdigit",10000);
394 staticdigits.Clear();
395 TClonesArray *digits = & staticdigits;
396 branch->SetAddress(&digits);
401 Int_t bufsize = 32000;
402 fTreeR->Branch("TOF", &fRecPoints, bufsize);
404 digitsTree->GetEvent(0);
405 Int_t nDigits = digits->GetEntriesFast();
406 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
409 Int_t dig[5]; //cluster detector indeces
410 Int_t parTOF[7]; //The TOF signal parameters
411 Bool_t status=kTRUE; // assume all sim channels ok in the beginning...
412 for (ii=0; ii<nDigits; ii++) {
413 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
414 dig[0]=d->GetSector();
415 dig[1]=d->GetPlate();
416 dig[2]=d->GetStrip();
420 /* check valid index */
421 if (dig[0]==-1||dig[1]==-1||dig[2]==-1||dig[3]==-1||dig[4]==-1) continue;
423 // Do not reconstruct anything in the holes
424 if (dig[0]==13 || dig[0]==14 || dig[0]==15 ) { // sectors with holes
425 if (dig[1]==2) { // plate with holes
431 // AliDebug(2,Form(" %2d %1d %2d %1d %2d ",dig[0],dig[1],dig[2],dig[3],dig[4]));
433 parTOF[0] = d->GetTdc(); //the TDC signal
434 parTOF[1] = d->GetToT(); //the ToT signal
435 parTOF[2] = d->GetAdc(); // the adc charge
436 parTOF[3] = d->GetTdcND(); // non decalibrated sim time
437 parTOF[4] = d->GetTdc(); // raw time, == Tdc time for the moment
438 parTOF[5] = 0; // deltaBC
439 parTOF[6] = 0; // L0-L1 latency
443 UShort_t volIdClus=GetClusterVolIndex(dig);
444 GetClusterPars(dig,posClus,covClus);
445 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);
446 InsertCluster(tofCluster);
450 AliDebug(1,Form("Number of found clusters: %d", fNumberOfTofClusters));
458 AliDebug(1,Form("Execution time to read TOF digits and to write TOF clusters : R:%.4fs C:%.4fs",
459 stopwatch.RealTime(),stopwatch.CpuTime()));
460 if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
463 //______________________________________________________________________________
465 void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
469 // Converts RAW data to recpoints for TOF
472 TStopwatch stopwatch;
477 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
481 Int_t bufsize = 32000;
482 clustersTree->Branch("TOF", &fRecPoints, bufsize);
484 TClonesArray * clonesRawData;
488 Int_t detectorIndex[5];
492 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
494 //AliTOFRawStream tofInput(rawReader);
495 fTOFRawStream.Clear();
496 fTOFRawStream.SetRawReader(rawReader);
499 AliInfo("Using New Decoder");
502 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
506 fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
507 else fTOFRawStream.LoadRawData(indexDDL);
509 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
511 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
513 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
515 //if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
516 if (tofRawDatum->GetTOF()==-1) continue;
519 if (indexDDL<10) ftxt << " " << indexDDL;
520 else ftxt << " " << indexDDL;
521 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
522 else ftxt << " " << tofRawDatum->GetTRM();
523 ftxt << " " << tofRawDatum->GetTRMchain();
524 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
525 else ftxt << " " << tofRawDatum->GetTDC();
526 ftxt << " " << tofRawDatum->GetTDCchannel();
529 fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
530 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
531 dummy = detectorIndex[3];
532 detectorIndex[3] = detectorIndex[4];
533 detectorIndex[4] = dummy;
536 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
537 else ftxt << " -> " << detectorIndex[0];
538 ftxt << " " << detectorIndex[1];
539 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
540 else ftxt << " " << detectorIndex[2];
541 ftxt << " " << detectorIndex[3];
542 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
543 else ftxt << " " << detectorIndex[4];
546 /* check valid index */
547 if (detectorIndex[0]==-1||detectorIndex[1]==-1||detectorIndex[2]==-1||detectorIndex[3]==-1||detectorIndex[4]==-1) continue;
549 // Do not reconstruct anything in the holes
550 if (detectorIndex[0]==13 || detectorIndex[0]==14 || detectorIndex[0]==15 ) { // sectors with holes
551 if (detectorIndex[1]==2) { // plate with holes
557 parTOF[0] = tofRawDatum->GetTOF(); //TDC
558 parTOF[1] = tofRawDatum->GetTOT(); // TOT
559 parTOF[2] = tofRawDatum->GetTOT(); //ADC==TOF
560 parTOF[3] = -1;//raw data: no track of undecalib sim time
561 parTOF[4] = tofRawDatum->GetTOF(); // RAW time
562 parTOF[5] = tofRawDatum->GetDeltaBC(); // deltaBC
563 parTOF[6] = tofRawDatum->GetL0L1Latency(); // L0-L1 latency
566 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
567 Int_t lab[3]={-1,-1,-1};
569 GetClusterPars(detectorIndex,posClus,covClus);
570 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);
571 InsertCluster(tofCluster);
574 if (parTOF[1]<10)ftxt << " " << parTOF[1];
575 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
576 else ftxt << " " << parTOF[1];
577 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
578 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
579 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
580 else ftxt << " " << parTOF[3] << endl;
583 } // closed loop on TOF raw data per current DDL file
585 clonesRawData->Clear();
587 } // closed loop on DDL index
589 if (fVerbose==2) ftxt.close();
591 AliDebug(1,Form("Number of found clusters: %d", fNumberOfTofClusters));
593 CalibrateRecPoint(rawReader->GetTimestamp());
596 clustersTree->Fill();
600 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
601 stopwatch.RealTime(),stopwatch.CpuTime()));
602 if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
605 //______________________________________________________________________________
607 void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
610 // Converts RAW data to recpoints for TOF
613 TStopwatch stopwatch;
618 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
620 fRunLoader->GetEvent(iEvent);
622 AliDebug(2,Form(" Event number %2d ", iEvent));
624 fTreeR = fTOFLoader->TreeR();
627 fTOFLoader->MakeTree("R");
628 fTreeR = fTOFLoader->TreeR();
631 Int_t bufsize = 32000;
632 fTreeR->Branch("TOF", &fRecPoints, bufsize);
634 TClonesArray * clonesRawData;
638 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
641 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
643 //AliTOFRawStream tofInput(rawReader);
644 fTOFRawStream.Clear();
645 fTOFRawStream.SetRawReader(rawReader);
648 AliInfo("Using New Decoder");
651 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
655 fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
656 else fTOFRawStream.LoadRawData(indexDDL);
658 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
660 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
662 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
664 //if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
665 if (tofRawDatum->GetTOF()==-1) continue;
668 if (indexDDL<10) ftxt << " " << indexDDL;
669 else ftxt << " " << indexDDL;
670 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
671 else ftxt << " " << tofRawDatum->GetTRM();
672 ftxt << " " << tofRawDatum->GetTRMchain();
673 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
674 else ftxt << " " << tofRawDatum->GetTDC();
675 ftxt << " " << tofRawDatum->GetTDCchannel();
678 fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
679 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
680 dummy = detectorIndex[3];
681 detectorIndex[3] = detectorIndex[4];
682 detectorIndex[4] = dummy;
685 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
686 else ftxt << " -> " << detectorIndex[0];
687 ftxt << " " << detectorIndex[1];
688 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
689 else ftxt << " " << detectorIndex[2];
690 ftxt << " " << detectorIndex[3];
691 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
692 else ftxt << " " << detectorIndex[4];
695 /* check valid index */
696 if (detectorIndex[0]==-1||detectorIndex[1]==-1||detectorIndex[2]==-1||detectorIndex[3]==-1||detectorIndex[4]==-1) continue;
698 // Do not reconstruct anything in the holes
699 if (detectorIndex[0]==13 || detectorIndex[0]==14 || detectorIndex[0]==15 ) { // sectors with holes
700 if (detectorIndex[1]==2) { // plate with holes
706 parTOF[0] = tofRawDatum->GetTOF(); // TDC
707 parTOF[1] = tofRawDatum->GetTOT(); // TOT
708 parTOF[2] = tofRawDatum->GetTOT(); // raw data have ADC=TOT
709 parTOF[3] = -1; //raw data: no track of the undecalib sim time
710 parTOF[4] = tofRawDatum->GetTOF(); // Raw time == TDC
711 parTOF[5] = tofRawDatum->GetDeltaBC(); // deltaBC
712 parTOF[6] = tofRawDatum->GetL0L1Latency(); // L0-L1 latency
715 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
716 Int_t lab[3]={-1,-1,-1};
718 GetClusterPars(detectorIndex,posClus,covClus);
719 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);
720 InsertCluster(tofCluster);
723 if (parTOF[1]<10)ftxt << " " << parTOF[1];
724 else if (parTOF[1]>=10 && parTOF[1]<100) ftxt << " " << parTOF[1];
725 else ftxt << " " << parTOF[1];
726 if (parTOF[0]<10) ftxt << " " << parTOF[0] << endl;
727 else if (parTOF[0]>=10 && parTOF[0]<100) ftxt << " " << parTOF[0] << endl;
728 else if (parTOF[0]>=100 && parTOF[0]<1000) ftxt << " " << parTOF[0] << endl;
729 else ftxt << " " << parTOF[3] << endl;
732 } // closed loop on TOF raw data per current DDL file
734 clonesRawData->Clear();
736 } // closed loop on DDL index
738 if (fVerbose==2) ftxt.close();
740 AliDebug(1,Form("Number of found clusters: %d for event: %d", fNumberOfTofClusters, iEvent));
742 CalibrateRecPoint(rawReader->GetTimestamp());
748 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
749 fTOFLoader->WriteRecPoints("OVERWRITE");
751 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.4fs C:%.4fs",
752 stopwatch.RealTime(),stopwatch.CpuTime()));
753 if (inholes) AliWarning(Form("Clusters in the TOF holes: %d",inholes));
756 //______________________________________________________________________________
758 void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
761 // Converts RAW data to MC digits for TOF
763 // (temporary solution)
766 TStopwatch stopwatch;
769 const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
771 fRunLoader->GetEvent(iEvent);
773 fTreeD = fTOFLoader->TreeD();
776 AliInfo("TreeD re-creation");
778 fTOFLoader->MakeTree("D");
779 fTreeD = fTOFLoader->TreeD();
782 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
783 Int_t bufsize = 32000;
784 fTreeD->Branch("TOF", &tofDigits, bufsize);
786 fRunLoader->GetEvent(iEvent);
788 AliDebug(2,Form(" Event number %2d ", iEvent));
790 TClonesArray * clonesRawData;
794 Int_t detectorIndex[5];
797 //AliTOFRawStream tofInput(rawReader);
798 fTOFRawStream.Clear();
799 fTOFRawStream.SetRawReader(rawReader);
802 AliInfo("Using New Decoder");
805 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
809 fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
810 else fTOFRawStream.LoadRawData(indexDDL);
812 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
814 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
816 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
818 //if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
819 if (tofRawDatum->GetTOF()==-1) continue;
821 fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
822 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
823 dummy = detectorIndex[3];
824 detectorIndex[3] = detectorIndex[4];
825 detectorIndex[4] = dummy;
827 digit[0] = fTOFRawStream.GetTofBin();
828 digit[1] = fTOFRawStream.GetToTbin();
829 digit[2] = fTOFRawStream.GetToTbin();
832 Int_t tracknum[3]={-1,-1,-1};
834 TClonesArray &aDigits = *tofDigits;
835 Int_t last=tofDigits->GetEntriesFast();
836 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
840 clonesRawData->Clear();
846 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
847 fTOFLoader->WriteDigits("OVERWRITE");
851 AliDebug(1, Form("Execution time to read TOF raw data and to write TOF clusters : R:%.2fs C:%.2fs",
852 stopwatch.RealTime(),stopwatch.CpuTime()));
856 //______________________________________________________________________________
858 void AliTOFClusterFinder::Raw2Digits(AliRawReader *rawReader, TTree* digitsTree)
861 // Converts RAW data to MC digits for TOF for the current event
865 TStopwatch stopwatch;
868 const Int_t kDDL = AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors();
872 AliError("No input digits Tree");
876 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
877 Int_t bufsize = 32000;
878 digitsTree->Branch("TOF", &tofDigits, bufsize);
880 TClonesArray * clonesRawData;
884 Int_t detectorIndex[5];
887 //AliTOFRawStream tofInput(rawReader);
888 fTOFRawStream.Clear();
889 fTOFRawStream.SetRawReader(rawReader);
892 AliInfo("Using New Decoder");
895 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
899 fTOFRawStream.LoadRawDataBuffers(indexDDL,fVerbose);
900 else fTOFRawStream.LoadRawData(indexDDL);
902 clonesRawData = (TClonesArray*)fTOFRawStream.GetRawData();
904 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
906 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
908 //if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
909 if (tofRawDatum->GetTOF()==-1) continue;
911 fTOFRawStream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
912 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
913 dummy = detectorIndex[3];
914 detectorIndex[3] = detectorIndex[4];
915 detectorIndex[4] = dummy;
917 digit[0] = fTOFRawStream.GetTofBin();
918 digit[1] = fTOFRawStream.GetToTbin();
919 digit[2] = fTOFRawStream.GetToTbin();
922 Int_t tracknum[3]={-1,-1,-1};
924 TClonesArray &aDigits = *tofDigits;
925 Int_t last=tofDigits->GetEntriesFast();
926 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
930 clonesRawData->Clear();
938 AliDebug(1, Form("Got %d digits: ", fDigits->GetEntries()));
939 AliDebug(1, Form("Execution time to read TOF raw data and fill TOF digit tree : R:%.2fs C:%.2fs",
940 stopwatch.RealTime(),stopwatch.CpuTime()));
943 //______________________________________________________________________________
945 Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
946 //---------------------------------------------------------------------------//
947 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
948 //---------------------------------------------------------------------------//
949 if (fNumberOfTofClusters==kTofMaxCluster) {
950 AliError("Too many clusters !");
954 if (fNumberOfTofClusters==0) {
955 fTofClusters[fNumberOfTofClusters++] = tofCluster;
959 Int_t ii = FindClusterIndex(tofCluster->GetZ());
960 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
961 fTofClusters[ii] = tofCluster;
962 fNumberOfTofClusters++;
967 //_________________________________________________________________________
969 Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
970 //--------------------------------------------------------------------
971 // This function returns the index of the nearest cluster in z
972 //--------------------------------------------------------------------
973 if (fNumberOfTofClusters==0) return 0;
974 if (z <= fTofClusters[0]->GetZ()) return 0;
975 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
976 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
977 for (; b<e; m=(b+e)/2) {
978 if (z > fTofClusters[m]->GetZ()) b=m+1;
985 //_________________________________________________________________________
987 void AliTOFClusterFinder::FillRecPoint()
990 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
991 // in Z) in the global TClonesArray of AliTOFcluster,
997 Int_t detectorIndex[5];
999 Int_t trackLabels[3];
1000 Int_t digitIndex = -1;
1001 Bool_t status=kTRUE;
1003 TClonesArray &lRecPoints = *fRecPoints;
1005 for (ii=0; ii<fNumberOfTofClusters; ii++) {
1007 digitIndex = fTofClusters[ii]->GetIndex();
1008 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
1009 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
1010 parTOF[0] = fTofClusters[ii]->GetTDC(); // TDC
1011 parTOF[1] = fTofClusters[ii]->GetToT(); // TOT
1012 parTOF[2] = fTofClusters[ii]->GetADC(); // ADC=TOT
1013 parTOF[3] = fTofClusters[ii]->GetTDCND(); // TDCND
1014 parTOF[4] = fTofClusters[ii]->GetTDCRAW();//RAW
1015 parTOF[5] = fTofClusters[ii]->GetDeltaBC();//deltaBC
1016 parTOF[6] = fTofClusters[ii]->GetL0L1Latency();//L0-L1 latency
1017 status=fTofClusters[ii]->GetStatus();
1018 Double_t posClus[3];
1019 Double_t covClus[6];
1020 UShort_t volIdClus=GetClusterVolIndex(detectorIndex);
1021 GetClusterPars(detectorIndex,posClus,covClus);
1022 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);
1024 AliDebug(2, Form(" %4d %4d %f %f %f %f %f %f %f %f %f %3d %3d %3d %2d %1d %2d %1d %2d %4d %3d %3d %4d %4d %1d %4d",
1025 ii, volIdClus, posClus[0], posClus[1], posClus[2],
1026 fTofClusters[ii]->GetSigmaX2(),
1027 fTofClusters[ii]->GetSigmaXY(),
1028 fTofClusters[ii]->GetSigmaXZ(),
1029 fTofClusters[ii]->GetSigmaY2(),
1030 fTofClusters[ii]->GetSigmaYZ(),
1031 fTofClusters[ii]->GetSigmaZ2(),
1032 trackLabels[0], trackLabels[1], trackLabels[2],
1033 detectorIndex[0], detectorIndex[1], detectorIndex[2], detectorIndex[3], detectorIndex[4],
1034 parTOF[0], parTOF[1], parTOF[2], parTOF[3], parTOF[4],
1035 status, digitIndex));
1037 } // loop on clusters
1041 //_________________________________________________________________________
1042 void AliTOFClusterFinder::CalibrateRecPoint(UInt_t timestamp)
1045 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
1046 // in Z) in the global TClonesArray of AliTOFcluster,
1052 Int_t detectorIndex[5];
1053 Int_t digitIndex = -1;
1057 Float_t tdcLatencyWindow;
1058 AliDebug(1," Calibrating TOF Clusters");
1060 AliTOFChannelOnlineArray *calDelay = fTOFcalib->GetTOFOnlineDelay();
1061 AliTOFChannelOnlineStatusArray *calStatus = fTOFcalib->GetTOFOnlineStatus();
1062 TObjArray *calTOFArrayOffline = fTOFcalib->GetTOFCalArrayOffline();
1064 AliTOFDeltaBCOffset *deltaBCOffsetObj = fTOFcalib->GetDeltaBCOffset();
1065 Int_t deltaBCOffset = deltaBCOffsetObj->GetDeltaBCOffset();
1066 AliTOFCTPLatency *ctpLatencyObj = fTOFcalib->GetCTPLatency();
1067 Float_t ctpLatency = ctpLatencyObj->GetCTPLatency();
1068 AliTOFRunParams *runParamsObj = fTOFcalib->GetRunParams();
1069 Float_t t0 = runParamsObj->EvalT0(timestamp);
1071 TString validity = (TString)fTOFcalib->GetOfflineValidity();
1072 Int_t calibration = -1;
1073 if (validity.CompareTo("valid")==0) {
1074 //AliInfo(Form(" validity = %s - Using offline calibration parameters",validity.Data()));
1077 //AliInfo(Form(" validity = %s - Using online calibration parameters",validity.Data()));
1081 for (ii=0; ii<fNumberOfTofClusters; ii++) {
1082 digitIndex = fTofClusters[ii]->GetIndex();
1083 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
1085 Int_t index = AliTOFGeometry::GetIndex(detectorIndex);
1087 UChar_t statusPulser=calStatus->GetPulserStatus(index);
1088 UChar_t statusNoise=calStatus->GetNoiseStatus(index);
1089 UChar_t statusHW=calStatus->GetHWStatus(index);
1090 UChar_t status=calStatus->GetStatus(index);
1091 tdcLatencyWindow = calStatus->GetLatencyWindow(index) * 1.e3; /* ns -> ps */
1093 //check the status, also unknown is fine!!!!!!!
1095 AliDebug(2, Form(" Status for channel %d = %d",index, (Int_t)status));
1096 if((statusPulser & AliTOFChannelOnlineStatusArray::kTOFPulserBad)==(AliTOFChannelOnlineStatusArray::kTOFPulserBad)||(statusNoise & AliTOFChannelOnlineStatusArray::kTOFNoiseBad)==(AliTOFChannelOnlineStatusArray::kTOFNoiseBad)||(statusHW & AliTOFChannelOnlineStatusArray::kTOFHWBad)==(AliTOFChannelOnlineStatusArray::kTOFHWBad)){
1097 AliDebug(2, Form(" Bad Status for channel %d",index));
1098 fTofClusters[ii]->SetStatus(kFALSE); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
1101 AliDebug(2, Form(" Good Status for channel %d",index));
1103 // Get Rough channel online equalization
1104 Double_t roughDelay=(Double_t)calDelay->GetDelay(index); // in ns
1105 AliDebug(2,Form(" channel delay (ns) = %f", roughDelay));
1106 // Get Refined channel offline calibration parameters
1107 if (calibration ==1){
1108 AliTOFChannelOffline * calChannelOffline = (AliTOFChannelOffline*)calTOFArrayOffline->At(index);
1110 for (Int_t j = 0; j<6; j++){
1111 par[j]=(Double_t)calChannelOffline->GetSlewPar(j);
1113 AliDebug(2,Form(" Calib Pars = %f, %f, %f, %f, %f, %f ",par[0],par[1],par[2],par[3],par[4],par[5]));
1114 AliDebug(2,Form(" The ToT and Time, uncorr (counts) = %d , %d", fTofClusters[ii]->GetToT(),fTofClusters[ii]->GetTDC()));
1115 tToT = (Double_t)(fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth());
1116 tToT*=1.E-3; //ToT in ns
1118 /* check TOT limits and set new TOT in case */
1119 if (tToT < AliTOFGeometry::SlewTOTMin()) tToT = AliTOFGeometry::SlewTOTMin();
1120 if (tToT > AliTOFGeometry::SlewTOTMax()) tToT = AliTOFGeometry::SlewTOTMax();
1122 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,tToT));
1123 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)
1126 timeCorr = roughDelay; // correction in ns
1128 AliDebug(2,Form(" The ToT and Time, uncorr (ns)= %e, %e",fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()*1.E-3,fTofClusters[ii]->GetToT()*AliTOFGeometry::ToTBinWidth()));
1129 AliDebug(2,Form(" The time correction (ns) = %f", timeCorr));
1130 timeCorr=(Double_t)(fTofClusters[ii]->GetTDC())*AliTOFGeometry::TdcBinWidth()*1.E-3-timeCorr;//redefine the time
1132 AliDebug(2,Form(" The channel time, corr (ps)= %e",timeCorr ));
1134 /* here timeCorr should be already corrected for calibration.
1135 * we now go into further corrections keeping in mind that timeCorr
1138 * the following corrections are performed in this way:
1140 * time = timeRaw - deltaBC + L0L1Latency + CTPLatency - TDCLatencyWindow - T0
1144 AliDebug(2, Form("applying further corrections (DeltaBC): DeltaBC=%d (BC bins), DeltaBCoffset=%d (BC bins)", fTofClusters[ii]->GetDeltaBC(), deltaBCOffset));
1145 AliDebug(2, Form("applying further corrections (L0L1Latency): L0L1Latency=%d (BC bins)", fTofClusters[ii]->GetL0L1Latency()));
1146 AliDebug(2, Form("applying further corrections (CTPLatency): CTPLatency=%f (ps)", ctpLatency));
1147 AliDebug(2, Form("applying further corrections (TDCLatencyWindow): TDCLatencyWindow=%f (ps)", tdcLatencyWindow));
1148 AliDebug(2, Form("applying further corrections (T0): T0=%f (ps)", t0));
1150 /* deltaBC correction (inhibited for the time being) */
1151 // timeCorr -= (fTofClusters[ii]->GetDeltaBC() - deltaBCOffset) * AliTOFGeometry::BunchCrossingBinWidth();
1152 /* L0L1-latency correction */
1153 timeCorr += fTofClusters[ii]->GetL0L1Latency() * AliTOFGeometry::BunchCrossingBinWidth();
1154 /* CTP-latency correction (from OCDB) */
1155 timeCorr += ctpLatency;
1156 /* TDC latency-window correction (from OCDB) */
1157 timeCorr -= tdcLatencyWindow;
1158 /* T0 correction (from OCDB) */
1162 * end of furhter corrections
1165 tdcCorr=(Int_t)(timeCorr/AliTOFGeometry::TdcBinWidth()); //the corrected time (tdc counts)
1166 fTofClusters[ii]->SetTDC(tdcCorr);
1167 } // loop on clusters
1170 //______________________________________________________________________________
1172 void AliTOFClusterFinder::ResetRecpoint()
1175 // Clear the list of reconstructed points
1178 fNumberOfTofClusters = 0;
1179 if (fRecPoints) fRecPoints->Clear();
1182 //______________________________________________________________________________
1184 void AliTOFClusterFinder::Load()
1187 // Load TOF.Digits.root and TOF.RecPoints.root files
1190 fTOFLoader->LoadDigits("READ");
1191 fTOFLoader->LoadRecPoints("recreate");
1194 //______________________________________________________________________________
1196 void AliTOFClusterFinder::LoadClusters()
1199 // Load TOF.RecPoints.root file
1202 fTOFLoader->LoadRecPoints("recreate");
1205 //______________________________________________________________________________
1207 void AliTOFClusterFinder::UnLoad()
1210 // Unload TOF.Digits.root and TOF.RecPoints.root files
1213 fTOFLoader->UnloadDigits();
1214 fTOFLoader->UnloadRecPoints();
1217 //______________________________________________________________________________
1219 void AliTOFClusterFinder::UnLoadClusters()
1222 // Unload TOF.RecPoints.root file
1225 fTOFLoader->UnloadRecPoints();
1228 //-------------------------------------------------------------------------
1229 UShort_t AliTOFClusterFinder::GetClusterVolIndex(const Int_t * const ind) const {
1231 //First of all get the volume ID to retrieve the l2t transformation...
1233 // Detector numbering scheme
1240 Int_t isector =ind[0];
1241 if (isector >= nSector)
1242 AliError(Form("Wrong sector number in TOF (%d) !",isector));
1243 Int_t iplate = ind[1];
1244 if (iplate >= nPlate)
1245 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1246 Int_t istrip = ind[2];
1248 Int_t stripOffset = 0;
1254 stripOffset = nStripC;
1257 stripOffset = nStripC+nStripB;
1260 stripOffset = nStripC+nStripB+nStripA;
1263 stripOffset = nStripC+nStripB+nStripA+nStripB;
1266 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
1270 Int_t index= (2*(nStripC+nStripB)+nStripA)*isector +
1274 UShort_t volIndex = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,index);
1278 //-------------------------------------------------------------------------
1279 void AliTOFClusterFinder::GetClusterPars(Int_t *ind, Double_t* pos,Double_t* cov) const {
1281 //First of all get the volume ID to retrieve the l2t transformation...
1283 UShort_t volIndex = GetClusterVolIndex(ind);
1286 //we now go in the system of the strip: determine the local coordinates
1289 // 47---------------------------------------------------0 ^ z
1290 // | | | | | | | | | | | | | | | | | | | | | | | | | | | 1 |
1291 // ----------------------------------------------------- | y going outwards
1292 // | | | | | | | | | | | | | | | | | | | | | | | | | | | 0 | par[0]=0;
1294 // ----------------------------------------------------- |
1295 // x <-----------------------------------------------------
1298 Float_t localX = (ind[4]-23.5)*AliTOFGeometry::XPad();
1300 Float_t localZ = (ind[3]- 0.5)*AliTOFGeometry::ZPad();
1302 Float_t localX = (ind[4]-AliTOFGeometry::NpadX()/2)*AliTOFGeometry::XPad()+AliTOFGeometry::XPad()/2.;
1304 Float_t localZ = (ind[3]-AliTOFGeometry::NpadZ()/2)*AliTOFGeometry::ZPad()+AliTOFGeometry::ZPad()/2.;
1305 //move to the tracking ref system
1312 const TGeoHMatrix *l2t= AliGeomManager::GetTracking2LocalMatrix(volIndex);
1313 // Get The position in the track ref system
1315 l2t->MasterToLocal(lpos,tpos);
1320 //Get The cluster covariance in the track ref system
1323 //cluster covariance in the local system:
1328 lcov[0] = AliTOFGeometry::XPad()*AliTOFGeometry::XPad()/12.;
1336 lcov[8] = AliTOFGeometry::ZPad()*AliTOFGeometry::ZPad()/12.;
1338 //cluster covariance in the tracking system:
1340 m.SetRotation(lcov);
1342 m.MultiplyLeft(&l2t->Inverse());
1343 Double_t *tcov = m.GetRotationMatrix();
1344 cov[0] = tcov[0]; cov[1] = tcov[1]; cov[2] = tcov[2];
1345 cov[3] = tcov[4]; cov[4] = tcov[5];