]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFClusterFinder.cxx
Include Methods to derive TOF AlignObjs from Survey Data
[u/mrichter/AliRoot.git] / TOF / AliTOFClusterFinder.cxx
CommitLineData
d08a92dd 1/***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16/*
15ec34b9 17$Log$
c77771e2 18Revision 1.19 2007/03/08 15:41:20 arcelli
19set uncorrected times when filling RecPoints
20
3432ebfa 21Revision 1.18 2007/03/06 16:31:20 arcelli
22Add Uncorrected TOF Time signal
23
aa5476d8 24Revision 1.17 2007/02/28 18:09:11 arcelli
25Add protection against failed retrieval of the CDB cal object, now Reconstruction exits with AliFatal
26
a1cab04c 27Revision 1.16 2007/02/20 15:57:00 decaro
28Raw data update: to read the TOF raw data defined in UNPACKED mode
29
d08a92dd 30
31Revision 0.03 2005/07/28 A. De Caro
32 Implement public method
33 Raw2Digits(Int_t, AliRawReader *)
34 to convert digits from raw data in MC digits
35 (temporary solution)
36
37Revision 0.02 2005/07/27 A. De Caro
38 Implement public method
39 Digits2RecPoint(Int_t)
40 to convert digits in clusters
41
42Revision 0.02 2005/07/26 A. De Caro
43 Implement private methods
44 InsertCluster(AliTOFcluster *)
45 FindClusterIndex(Double_t)
46 originally implemented in AliTOFtracker
47 by S. Arcelli and C. Zampolli
48
49Revision 0.01 2005/07/25 A. De Caro
50 Implement public methods
51 Digits2RecPoint(AliRawReader *, TTree *)
52 Digits2RecPoint(Int_t, AliRawReader *)
53 to convert raw data in clusters
54 */
55
56////////////////////////////////////////////////////////////////
57// //
58// Class for TOF cluster finder //
59// //
60// Starting from Raw Data, create rec points, //
61// fill TreeR for TOF, //
62// write TOF.RecPoints.root file //
63// //
64////////////////////////////////////////////////////////////////
65
15ec34b9 66
0e46b9ae 67#include "TClonesArray.h"
15ec34b9 68//#include "TFile.h"
0e46b9ae 69#include "TTree.h"
d08a92dd 70
d0eb8f39 71#include "AliDAQ.h"
0e46b9ae 72#include "AliLoader.h"
d08a92dd 73#include "AliLog.h"
0e46b9ae 74#include "AliRawReader.h"
d08a92dd 75#include "AliRunLoader.h"
d08a92dd 76
37879eed 77#include "AliTOFCal.h"
15ec34b9 78#include "AliTOFcalib.h"
37879eed 79#include "AliTOFChannel.h"
d08a92dd 80#include "AliTOFClusterFinder.h"
0e46b9ae 81#include "AliTOFcluster.h"
82#include "AliTOFdigit.h"
0e46b9ae 83#include "AliTOFGeometry.h"
15ec34b9 84#include "AliTOFGeometryV5.h"
85#include "AliTOFrawData.h"
0e46b9ae 86#include "AliTOFRawStream.h"
3432ebfa 87#include "Riostream.h"
0e46b9ae 88
15ec34b9 89//extern TFile *gFile;
d08a92dd 90
91ClassImp(AliTOFClusterFinder)
92
93AliTOFClusterFinder::AliTOFClusterFinder():
94 fRunLoader(0),
95 fTOFLoader(0),
96 fTreeD(0),
97 fTreeR(0),
58d8d9a3 98 fTOFGeometry(new AliTOFGeometryV5()),
d08a92dd 99 fDigits(new TClonesArray("AliTOFdigit", 4000)),
100 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
15ec34b9 101 fNumberOfTofClusters(0),
102 fVerbose(0)
d08a92dd 103{
104//
105// Constructor
106//
107
3bbdf45d 108 AliInfo("V5 TOF Geometry is taken as the default");
37879eed 109
d08a92dd 110}
111//______________________________________________________________________________
112
113AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
114 fRunLoader(runLoader),
115 fTOFLoader(runLoader->GetLoader("TOFLoader")),
116 fTreeD(0),
117 fTreeR(0),
58d8d9a3 118 fTOFGeometry(new AliTOFGeometryV5()),
d08a92dd 119 fDigits(new TClonesArray("AliTOFdigit", 4000)),
120 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
15ec34b9 121 fNumberOfTofClusters(0),
122 fVerbose(0)
d08a92dd 123{
124//
125// Constructor
126//
127
0f49d1bc 128// runLoader->CdGAFile();
129// TFile *in=(TFile*)gFile;
130// in->cd();
131// fTOFGeometry = (AliTOFGeometry*)in->Get("TOFgeometry");
d08a92dd 132
7aeeaf38 133}
134
135//------------------------------------------------------------------------
136AliTOFClusterFinder::AliTOFClusterFinder(const AliTOFClusterFinder &source)
58d8d9a3 137 :TObject(),
138 fRunLoader(0),
139 fTOFLoader(0),
140 fTreeD(0),
141 fTreeR(0),
142 fTOFGeometry(new AliTOFGeometryV5()),
143 fDigits(new TClonesArray("AliTOFdigit", 4000)),
144 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
15ec34b9 145 fNumberOfTofClusters(0),
146 fVerbose(0)
7aeeaf38 147{
148 // copy constructor
149 this->fDigits=source.fDigits;
150 this->fRecPoints=source.fRecPoints;
151 this->fTOFGeometry=source.fTOFGeometry;
152
153}
154
155//------------------------------------------------------------------------
d0eb8f39 156AliTOFClusterFinder& AliTOFClusterFinder::operator=(const AliTOFClusterFinder &source)
7aeeaf38 157{
158 // ass. op.
159 this->fDigits=source.fDigits;
160 this->fRecPoints=source.fRecPoints;
161 this->fTOFGeometry=source.fTOFGeometry;
15ec34b9 162 this->fVerbose=source.fVerbose;
7aeeaf38 163 return *this;
164
d08a92dd 165}
166//______________________________________________________________________________
167
168AliTOFClusterFinder::~AliTOFClusterFinder()
169{
170
171 //
172 // Destructor
173 //
174
175 if (fDigits)
176 {
177 fDigits->Delete();
178 delete fDigits;
179 fDigits=0;
180 }
181 if (fRecPoints)
182 {
183 fRecPoints->Delete();
184 delete fRecPoints;
185 fRecPoints=0;
186 }
187
188 delete fTOFGeometry;
189
190}
191//______________________________________________________________________________
192
193void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
194{
195 //
196 // Converts digits to recpoints for TOF
197 //
198
199 fRunLoader->GetEvent(iEvent);
200
201 fTreeD = fTOFLoader->TreeD();
202 if (fTreeD == 0x0)
203 {
204 AliFatal("AliTOFClusterFinder: Can not get TreeD");
205 }
206
207 TBranch *branch = fTreeD->GetBranch("TOF");
208 if (!branch) {
209 AliError("can't get the branch with the TOF digits !");
210 return;
211 }
212
15ec34b9 213 TClonesArray *digits = new TClonesArray("AliTOFdigit",10000);
d08a92dd 214 branch->SetAddress(&digits);
215
216 ResetRecpoint();
217
218 fTreeR = fTOFLoader->TreeR();
219 if (fTreeR == 0x0)
220 {
221 fTOFLoader->MakeTree("R");
222 fTreeR = fTOFLoader->TreeR();
223 }
224
225 Int_t bufsize = 32000;
226 fTreeR->Branch("TOF", &fRecPoints, bufsize);
227
228 fTreeD->GetEvent(0);
229 Int_t nDigits = digits->GetEntriesFast();
230 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
231
232 Int_t ii, jj;
233 Int_t dig[5];
234 Float_t g[3];
235 Double_t h[5];
340693af 236 Float_t tToT;
237 Double_t tTdcND;
d08a92dd 238 for (ii=0; ii<nDigits; ii++) {
239 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
240 dig[0]=d->GetSector();
241 dig[1]=d->GetPlate();
242 dig[2]=d->GetStrip();
243 dig[3]=d->GetPadz();
244 dig[4]=d->GetPadx();
245
d0eb8f39 246 //AliInfo(Form(" %2i %1i %2i %1i %2i ",dig[0],dig[1],dig[2],dig[3],dig[4]));
247
d08a92dd 248 for (jj=0; jj<3; jj++) g[jj] = 0.;
249 fTOFGeometry->GetPos(dig,g);
250
251 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
252 h[1] = TMath::ATan2(g[1],g[0]);
253 h[2] = g[2];
254 h[3] = d->GetTdc();
255 h[4] = d->GetAdc();
340693af 256 tToT = d->GetToT();
257 tTdcND = d->GetTdcND();
d08a92dd 258
340693af 259 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii,tToT, tTdcND);
aa5476d8 260 tofCluster->SetTDCRAW(d->GetTdc());
d08a92dd 261 InsertCluster(tofCluster);
262
263 }
264
265 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
266
37879eed 267 CalibrateRecPoint();
d08a92dd 268 FillRecPoint();
269
270 fTreeR->Fill();
271 ResetRecpoint();
272
273 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
274 fTOFLoader->WriteRecPoints("OVERWRITE");
275
276}
277//______________________________________________________________________________
278
279void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
280 TTree *clustersTree)
281{
282 //
283 // Converts RAW data to recpoints for TOF
284 //
285
d0eb8f39 286 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
287 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
d08a92dd 288
289 ResetRecpoint();
290
291 Int_t bufsize = 32000;
292 clustersTree->Branch("TOF", &fRecPoints, bufsize);
293
15ec34b9 294 TClonesArray * clonesRawData;
295
d08a92dd 296 Int_t ii = 0;
15ec34b9 297 Int_t dummy = -1;
d08a92dd 298
299 Int_t detectorIndex[5];
300 Float_t position[3];
301 Double_t cylindricalPosition[5];
340693af 302 Float_t tToT;
303 Double_t tTdcND;
d08a92dd 304
15ec34b9 305 ofstream ftxt;
306 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
307
308 Int_t indexDDL = 0;
309 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
310
311 rawReader->Reset();
312 AliTOFRawStream tofInput(rawReader);
313 tofInput.LoadRawData(indexDDL);
314
315 clonesRawData = (TClonesArray*)tofInput.GetRawData();
316
317 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
318
319 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
320
321 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
322
323 if (fVerbose==2) {
324 if (indexDDL<10) ftxt << " " << indexDDL;
325 else ftxt << " " << indexDDL;
326 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
327 else ftxt << " " << tofRawDatum->GetTRM();
328 ftxt << " " << tofRawDatum->GetTRMchain();
329 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
330 else ftxt << " " << tofRawDatum->GetTDC();
331 ftxt << " " << tofRawDatum->GetTDCchannel();
332 }
333
334 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
335 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
336 dummy = detectorIndex[3];
337 detectorIndex[3] = detectorIndex[4];
338 detectorIndex[4] = dummy;
339
340 if (fVerbose==2) {
341 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
342 else ftxt << " -> " << detectorIndex[0];
343 ftxt << " " << detectorIndex[1];
344 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
345 else ftxt << " " << detectorIndex[2];
346 ftxt << " " << detectorIndex[3];
347 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
348 else ftxt << " " << detectorIndex[4];
349 }
350
351 for (ii=0; ii<3; ii++) position[ii] = 0.;
352 fTOFGeometry->GetPos(detectorIndex, position);
353
354 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
355 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
356 cylindricalPosition[2] = position[2];
357 cylindricalPosition[3] = tofRawDatum->GetTOF();
358 cylindricalPosition[4] = tofRawDatum->GetTOT();
359 tToT = tofRawDatum->GetTOT();
360 tTdcND = -1.;
361 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
362 tofCluster->SetToT(tToT);
363 tofCluster->SetTDCND(tTdcND);
aa5476d8 364 tofCluster->SetTDCRAW(tofRawDatum->GetTOF());
15ec34b9 365 InsertCluster(tofCluster);
366
367 if (fVerbose==2) {
368 if (cylindricalPosition[4]<10) ftxt << " " << cylindricalPosition[4];
369 else if (cylindricalPosition[4]>=10 && cylindricalPosition[4]<100) ftxt << " " << cylindricalPosition[4];
370 else ftxt << " " << cylindricalPosition[4];
371 if (cylindricalPosition[3]<10) ftxt << " " << cylindricalPosition[3] << endl;
372 else if (cylindricalPosition[3]>=10 && cylindricalPosition[3]<100) ftxt << " " << cylindricalPosition[3] << endl;
373 else if (cylindricalPosition[3]>=100 && cylindricalPosition[3]<1000) ftxt << " " << cylindricalPosition[3] << endl;
374 else ftxt << " " << cylindricalPosition[3] << endl;
375 }
376
377 } // closed loop on TOF raw data per current DDL file
378
379 clonesRawData->Clear();
380
381 } // closed loop on DDL index
382
383 /*
384 Int_t indexDDL = 0;
d08a92dd 385 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
386
387 rawReader->Reset();
388 AliTOFRawStream tofInput(rawReader);
362c9d61 389 rawReader->Select("TOF", indexDDL, indexDDL);
d08a92dd 390
391 while(tofInput.Next()) {
392
d0eb8f39 393 for (ii=0; ii<5; ii++) detectorIndex[ii] = -1;
394
d08a92dd 395 detectorIndex[0] = tofInput.GetSector();
396 detectorIndex[1] = tofInput.GetPlate();
397 detectorIndex[2] = tofInput.GetStrip();
398 detectorIndex[3] = tofInput.GetPadZ();
399 detectorIndex[4] = tofInput.GetPadX();
400
d0eb8f39 401 //AliInfo(Form(" %2i %1i %2i %1i %2i ",detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
402
403 if (detectorIndex[0]==-1 ||
404 detectorIndex[1]==-1 ||
405 detectorIndex[2]==-1 ||
406 detectorIndex[3]==-1 ||
407 detectorIndex[4]==-1) continue;
408
d08a92dd 409 for (ii=0; ii<3; ii++) position[ii] = 0.;
410
411 fTOFGeometry->GetPos(detectorIndex, position);
412
413 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
414 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
415 cylindricalPosition[2] = position[2];
416 cylindricalPosition[3] = tofInput.GetTofBin();
d0eb8f39 417 cylindricalPosition[4] = tofInput.GetToTbin();
418 tToT = tofInput.GetToTbin();
340693af 419 tTdcND = -1.;
d08a92dd 420 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
340693af 421 tofCluster->SetToT(tToT);
422 tofCluster->SetTDCND(tTdcND);
d08a92dd 423 InsertCluster(tofCluster);
424
425 } // while loop
426
427 } // loop on DDL files
15ec34b9 428 */
429
430 if (fVerbose==2) ftxt.close();
d08a92dd 431
432 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
433
37879eed 434 CalibrateRecPoint();
d08a92dd 435 FillRecPoint();
436
437 clustersTree->Fill();
37879eed 438
d08a92dd 439 ResetRecpoint();
440
441}
442//______________________________________________________________________________
443
444void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
445{
446 //
447 // Converts RAW data to recpoints for TOF
448 //
449
d0eb8f39 450 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
451 const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
d08a92dd 452
453 fRunLoader->GetEvent(iEvent);
454
455 AliDebug(2,Form(" Event number %2i ", iEvent));
456
457 fTreeR = fTOFLoader->TreeR();
458
459 if (fTreeR == 0x0){
460 fTOFLoader->MakeTree("R");
461 fTreeR = fTOFLoader->TreeR();
462 }
463
464 Int_t bufsize = 32000;
465 fTreeR->Branch("TOF", &fRecPoints, bufsize);
466
15ec34b9 467 TClonesArray * clonesRawData;
468
d08a92dd 469 Int_t ii = 0;
15ec34b9 470 Int_t dummy = -1;
d08a92dd 471
d0eb8f39 472 Int_t detectorIndex[5] = {-1, -1, -1, -1, -1};
d08a92dd 473 Float_t position[3];
474 Double_t cylindricalPosition[5];
340693af 475 Float_t tToT;
476 Double_t tTdcND;
d08a92dd 477
15ec34b9 478 ofstream ftxt;
479 if (fVerbose==2) ftxt.open("TOFdigitsRead.txt",ios::app);
480
481 Int_t indexDDL = 0;
d08a92dd 482 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
483
484 rawReader->Reset();
485 AliTOFRawStream tofInput(rawReader);
15ec34b9 486 tofInput.LoadRawData(indexDDL);
487
488 clonesRawData = (TClonesArray*)tofInput.GetRawData();
489
490 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
491
492 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
493
494 if (tofRawDatum->GetTOT()==-1 || tofRawDatum->GetTOF()==-1) continue;
495
496 if (fVerbose==2) {
497 if (indexDDL<10) ftxt << " " << indexDDL;
498 else ftxt << " " << indexDDL;
499 if (tofRawDatum->GetTRM()<10) ftxt << " " << tofRawDatum->GetTRM();
500 else ftxt << " " << tofRawDatum->GetTRM();
501 ftxt << " " << tofRawDatum->GetTRMchain();
502 if (tofRawDatum->GetTDC()<10) ftxt << " " << tofRawDatum->GetTDC();
503 else ftxt << " " << tofRawDatum->GetTDC();
504 ftxt << " " << tofRawDatum->GetTDCchannel();
505 }
506
507 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
508 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
509 dummy = detectorIndex[3];
510 detectorIndex[3] = detectorIndex[4];
511 detectorIndex[4] = dummy;
512
513 if (fVerbose==2) {
514 if (detectorIndex[0]<10) ftxt << " -> " << detectorIndex[0];
515 else ftxt << " -> " << detectorIndex[0];
516 ftxt << " " << detectorIndex[1];
517 if (detectorIndex[2]<10) ftxt << " " << detectorIndex[2];
518 else ftxt << " " << detectorIndex[2];
519 ftxt << " " << detectorIndex[3];
520 if (detectorIndex[4]<10) ftxt << " " << detectorIndex[4];
521 else ftxt << " " << detectorIndex[4];
522 }
d0eb8f39 523
d08a92dd 524 for (ii=0; ii<3; ii++) position[ii] = 0.;
d08a92dd 525 fTOFGeometry->GetPos(detectorIndex, position);
15ec34b9 526
527 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
528 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
529 cylindricalPosition[2] = position[2];
530 cylindricalPosition[3] = tofRawDatum->GetTOF();
531 cylindricalPosition[4] = tofRawDatum->GetTOT();
532 tToT = tofRawDatum->GetTOT();
340693af 533 tTdcND = -1.;
d08a92dd 534 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
340693af 535 tofCluster->SetToT(tToT);
536 tofCluster->SetTDCND(tTdcND);
aa5476d8 537 tofCluster->SetTDCRAW(tofRawDatum->GetTOF());
d08a92dd 538 InsertCluster(tofCluster);
539
15ec34b9 540 if (fVerbose==2) {
541 if (cylindricalPosition[4]<10) ftxt << " " << cylindricalPosition[4];
542 else if (cylindricalPosition[4]>=10 && cylindricalPosition[4]<100) ftxt << " " << cylindricalPosition[4];
543 else ftxt << " " << cylindricalPosition[4];
544 if (cylindricalPosition[3]<10) ftxt << " " << cylindricalPosition[3] << endl;
545 else if (cylindricalPosition[3]>=10 && cylindricalPosition[3]<100) ftxt << " " << cylindricalPosition[3] << endl;
546 else if (cylindricalPosition[3]>=100 && cylindricalPosition[3]<1000) ftxt << " " << cylindricalPosition[3] << endl;
547 else ftxt << " " << cylindricalPosition[3] << endl;
548 }
d08a92dd 549
15ec34b9 550 } // closed loop on TOF raw data per current DDL file
551
552 clonesRawData->Clear();
553
554 } // closed loop on DDL index
555
556 if (fVerbose==2) ftxt.close();
d08a92dd 557
558 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
559
37879eed 560 CalibrateRecPoint();
d08a92dd 561 FillRecPoint();
562
563 fTreeR->Fill();
564 ResetRecpoint();
565
566 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
567 fTOFLoader->WriteRecPoints("OVERWRITE");
568
569}
570//______________________________________________________________________________
571
572void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
573{
574 //
575 // Converts RAW data to MC digits for TOF
576 //
577 // (temporary solution)
578 //
579
15ec34b9 580 //const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
581 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
582
d08a92dd 583 fRunLoader->GetEvent(iEvent);
584
585 fTreeD = fTOFLoader->TreeD();
586 if (fTreeD)
587 {
d0eb8f39 588 AliInfo("TreeD re-creation");
d08a92dd 589 fTreeD = 0x0;
590 fTOFLoader->MakeTree("D");
591 fTreeD = fTOFLoader->TreeD();
592 }
593
15ec34b9 594 TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
d08a92dd 595 Int_t bufsize = 32000;
596 fTreeD->Branch("TOF", &tofDigits, bufsize);
597
d08a92dd 598 fRunLoader->GetEvent(iEvent);
599
600 AliDebug(2,Form(" Event number %2i ", iEvent));
601
15ec34b9 602 TClonesArray * clonesRawData;
603
604 Int_t dummy = -1;
d08a92dd 605
606 Int_t detectorIndex[5];
d0eb8f39 607 Float_t digit[4];
d08a92dd 608
15ec34b9 609 Int_t indexDDL = 0;
d08a92dd 610 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
611
612 rawReader->Reset();
613 AliTOFRawStream tofInput(rawReader);
15ec34b9 614 tofInput.LoadRawData(indexDDL);
d08a92dd 615
15ec34b9 616 clonesRawData = (TClonesArray*)tofInput.GetRawData();
d08a92dd 617
15ec34b9 618 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
d0eb8f39 619
15ec34b9 620 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
d0eb8f39 621
15ec34b9 622 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
623
624 tofInput.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
625 tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
626 dummy = detectorIndex[3];
627 detectorIndex[3] = detectorIndex[4];
628 detectorIndex[4] = dummy;
d0eb8f39 629
d08a92dd 630 digit[0] = (Float_t)tofInput.GetTofBin();
d0eb8f39 631 digit[1] = (Float_t)tofInput.GetToTbin();
632 digit[2] = (Float_t)tofInput.GetToTbin();
633 digit[3] = -1.;
d08a92dd 634
635 Int_t tracknum[3]={-1,-1,-1};
636
637 TClonesArray &aDigits = *tofDigits;
638 Int_t last=tofDigits->GetEntriesFast();
639 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
640
641 } // while loop
642
15ec34b9 643 clonesRawData->Clear();
644
d08a92dd 645 } // DDL Loop
646
647 fTreeD->Fill();
648
649 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
650 fTOFLoader->WriteDigits("OVERWRITE");
651
652}
653//______________________________________________________________________________
654
655Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
656 //---------------------------------------------------------------------------//
657 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
658 //---------------------------------------------------------------------------//
659 if (fNumberOfTofClusters==kTofMaxCluster) {
660 AliError("Too many clusters !");
661 return 1;
662 }
663
664 if (fNumberOfTofClusters==0) {
665 fTofClusters[fNumberOfTofClusters++] = tofCluster;
666 return 0;
667 }
668
669 Int_t ii = FindClusterIndex(tofCluster->GetZ());
670 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
671 fTofClusters[ii] = tofCluster;
672 fNumberOfTofClusters++;
673
674 return 0;
675
676}
677//_________________________________________________________________________
678
679Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
680 //--------------------------------------------------------------------
681 // This function returns the index of the nearest cluster
682 //--------------------------------------------------------------------
683 if (fNumberOfTofClusters==0) return 0;
684 if (z <= fTofClusters[0]->GetZ()) return 0;
685 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
686 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
687 for (; b<e; m=(b+e)/2) {
688 if (z > fTofClusters[m]->GetZ()) b=m+1;
689 else e=m;
690 }
691
692 return m;
693
694}
695//_________________________________________________________________________
696
697void AliTOFClusterFinder::FillRecPoint()
698{
699 //
700 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
701 // in Z) in the global TClonesArray of AliTOFcluster,
702 // i.e. fRecPoints.
703 //
704
705 Int_t ii, jj;
706
707 Int_t detectorIndex[5];
708 Double_t cylindricalPosition[5];
709 Int_t trackLabels[3];
710 Int_t digitIndex = -1;
58d8d9a3 711 Float_t tToT=0.;
712 Double_t tTdcND=0.;
3432ebfa 713 Double_t tTdcRAW=0.;
58d8d9a3 714 Bool_t cStatus = kTRUE;
d08a92dd 715
716 TClonesArray &lRecPoints = *fRecPoints;
717
718 for (ii=0; ii<fNumberOfTofClusters; ii++) {
719
720 digitIndex = fTofClusters[ii]->GetIndex();
721 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
722 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
723 cylindricalPosition[0] = fTofClusters[ii]->GetR();
724 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
725 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
726 cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
727 cylindricalPosition[4] = fTofClusters[ii]->GetADC();
340693af 728 tToT = fTofClusters[ii]->GetToT();
729 tTdcND = fTofClusters[ii]->GetTDCND();
58d8d9a3 730 cStatus=fTofClusters[ii]->GetStatus();
3432ebfa 731 tTdcRAW=fTofClusters[ii]->GetTDCRAW();
732 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex, tToT, tTdcND, tTdcRAW,cStatus);
d08a92dd 733
734 //AliInfo(Form("%3i %3i %f %f %f %f %f %2i %2i %2i %1i %2i",ii,digitIndex, cylindricalPosition[2],cylindricalPosition[0],cylindricalPosition[1],cylindricalPosition[3],cylindricalPosition[4],detectorIndex[0],detectorIndex[1],detectorIndex[2],detectorIndex[3],detectorIndex[4]));
735
736 } // loop on clusters
737
738}
37879eed 739
740//_________________________________________________________________________
741void AliTOFClusterFinder::CalibrateRecPoint()
742{
743 //
744 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
745 // in Z) in the global TClonesArray of AliTOFcluster,
746 // i.e. fRecPoints.
747 //
748
749 Int_t ii, jj;
750
751 Int_t detectorIndex[5];
752 Int_t digitIndex = -1;
340693af 753 Float_t tToT;
754 Float_t tdcCorr;
37879eed 755 AliInfo(" Calibrating TOF Clusters: ")
756 AliTOFcalib *calib = new AliTOFcalib(fTOFGeometry);
3f2bddfa 757 // calib->ReadParFromCDB("TOF/Calib",0); // original
a1cab04c 758 // Use AliCDBManager's run number
759 if(!calib->ReadParFromCDB("TOF/Calib",-1)) {AliFatal("Exiting, no CDB object found!!!");exit(0);}
760
340693af 761 AliTOFCal *calTOFArray = calib->GetTOFCalArray();
37879eed 762
763 for (ii=0; ii<fNumberOfTofClusters; ii++) {
764 digitIndex = fTofClusters[ii]->GetIndex();
765 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
766
767 Int_t index = calib->GetIndex(detectorIndex);
768
340693af 769 AliTOFChannel * calChannel = calTOFArray->GetChannel(index);
58d8d9a3 770
771 // Get channel status
772 Bool_t status=calChannel->GetStatus();
773 if(status)fTofClusters[ii]->SetStatus(!status); //odd convention, to avoid conflict with calibration objects currently in the db (temporary solution).
774
775 // Get Rough channel online equalization
776 Float_t roughDelay=calChannel->GetDelay();
777
778 // Get Refined channel offline calibration parameters
37879eed 779 Float_t par[6];
780 for (Int_t j = 0; j<6; j++){
340693af 781 par[j]=calChannel->GetSlewPar(j);
37879eed 782 }
340693af 783 tToT = fTofClusters[ii]->GetToT();
58d8d9a3 784 Float_t timeCorr=par[0]+par[1]*tToT+par[2]*tToT*tToT+par[3]*tToT*tToT*tToT+par[4]*tToT*tToT*tToT*tToT+par[5]*tToT*tToT*tToT*tToT*tToT+roughDelay;
340693af 785 tdcCorr=(fTofClusters[ii]->GetTDC()*AliTOFGeometry::TdcBinWidth()+32)*1.E-3-timeCorr;
786 tdcCorr=(tdcCorr*1E3-32)/AliTOFGeometry::TdcBinWidth();
787 fTofClusters[ii]->SetTDC(tdcCorr);
37879eed 788
789 } // loop on clusters
790
791 delete calib;
792}
d08a92dd 793//______________________________________________________________________________
794
795void AliTOFClusterFinder::ResetRecpoint()
796{
797 //
798 // Clear the list of reconstructed points
799 //
800
801 fNumberOfTofClusters = 0;
802 if (fRecPoints) fRecPoints->Clear();
803
804}
805//______________________________________________________________________________
806
807void AliTOFClusterFinder::Load()
808{
809 //
810 // Load TOF.Digits.root and TOF.RecPoints.root files
811 //
812
813 fTOFLoader->LoadDigits("READ");
814 fTOFLoader->LoadRecPoints("recreate");
815
816}
817//______________________________________________________________________________
818
819void AliTOFClusterFinder::LoadClusters()
820{
821 //
822 // Load TOF.RecPoints.root file
823 //
824
825 fTOFLoader->LoadRecPoints("recreate");
826
827}
828//______________________________________________________________________________
829
830void AliTOFClusterFinder::UnLoad()
831{
832 //
833 // Unload TOF.Digits.root and TOF.RecPoints.root files
834 //
835
836 fTOFLoader->UnloadDigits();
837 fTOFLoader->UnloadRecPoints();
838
839}
840//______________________________________________________________________________
841
842void AliTOFClusterFinder::UnLoadClusters()
843{
844 //
845 // Unload TOF.RecPoints.root file
846 //
847
848 fTOFLoader->UnloadRecPoints();
849
850}
851//______________________________________________________________________________