]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFClusterFinder.cxx
ESDtrack TOF cluster index now pointing to recPoints (S.Arcelli, C.Zampolli)
[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/*
17
18Revision 0.03 2005/07/28 A. De Caro
19 Implement public method
20 Raw2Digits(Int_t, AliRawReader *)
21 to convert digits from raw data in MC digits
22 (temporary solution)
23
24Revision 0.02 2005/07/27 A. De Caro
25 Implement public method
26 Digits2RecPoint(Int_t)
27 to convert digits in clusters
28
29Revision 0.02 2005/07/26 A. De Caro
30 Implement private methods
31 InsertCluster(AliTOFcluster *)
32 FindClusterIndex(Double_t)
33 originally implemented in AliTOFtracker
34 by S. Arcelli and C. Zampolli
35
36Revision 0.01 2005/07/25 A. De Caro
37 Implement public methods
38 Digits2RecPoint(AliRawReader *, TTree *)
39 Digits2RecPoint(Int_t, AliRawReader *)
40 to convert raw data in clusters
41 */
42
43////////////////////////////////////////////////////////////////
44// //
45// Class for TOF cluster finder //
46// //
47// Starting from Raw Data, create rec points, //
48// fill TreeR for TOF, //
49// write TOF.RecPoints.root file //
50// //
51////////////////////////////////////////////////////////////////
52
53#include <Riostream.h>
54#include <TTree.h>
55#include <TObjArray.h>
56#include <TClonesArray.h>
57
58#include "AliLog.h"
59#include "AliRunLoader.h"
60#include "AliLoader.h"
61
62#include "AliTOFdigit.h"
63#include "AliTOFcluster.h"
64#include "AliTOFGeometry.h"
65#include "AliTOFRawStream.h"
66
67#include "AliTOFClusterFinder.h"
68
69ClassImp(AliTOFClusterFinder)
70
71AliTOFClusterFinder::AliTOFClusterFinder():
72 fRunLoader(0),
73 fTOFLoader(0),
74 fTreeD(0),
75 fTreeR(0),
76 fDigits(new TClonesArray("AliTOFdigit", 4000)),
77 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
78 fNumberOfTofClusters(0)
79{
80//
81// Constructor
82//
83
84 fTOFGeometry = new AliTOFGeometry();
85
86}
87//______________________________________________________________________________
88
89AliTOFClusterFinder::AliTOFClusterFinder(AliRunLoader* runLoader):
90 fRunLoader(runLoader),
91 fTOFLoader(runLoader->GetLoader("TOFLoader")),
92 fTreeD(0),
93 fTreeR(0),
94 fDigits(new TClonesArray("AliTOFdigit", 4000)),
95 fRecPoints(new TClonesArray("AliTOFcluster", 4000)),
96 fNumberOfTofClusters(0)
97{
98//
99// Constructor
100//
101
102 fTOFGeometry = new AliTOFGeometry();
103
104}
105//______________________________________________________________________________
106
107AliTOFClusterFinder::~AliTOFClusterFinder()
108{
109
110 //
111 // Destructor
112 //
113
114 if (fDigits)
115 {
116 fDigits->Delete();
117 delete fDigits;
118 fDigits=0;
119 }
120 if (fRecPoints)
121 {
122 fRecPoints->Delete();
123 delete fRecPoints;
124 fRecPoints=0;
125 }
126
127 delete fTOFGeometry;
128
129}
130//______________________________________________________________________________
131
132void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent)
133{
134 //
135 // Converts digits to recpoints for TOF
136 //
137
138 fRunLoader->GetEvent(iEvent);
139
140 fTreeD = fTOFLoader->TreeD();
141 if (fTreeD == 0x0)
142 {
143 AliFatal("AliTOFClusterFinder: Can not get TreeD");
144 }
145
146 TBranch *branch = fTreeD->GetBranch("TOF");
147 if (!branch) {
148 AliError("can't get the branch with the TOF digits !");
149 return;
150 }
151
152 TClonesArray dummy("AliTOFdigit",10000), *digits=&dummy;
153 branch->SetAddress(&digits);
154
155 ResetRecpoint();
156
157 fTreeR = fTOFLoader->TreeR();
158 if (fTreeR == 0x0)
159 {
160 fTOFLoader->MakeTree("R");
161 fTreeR = fTOFLoader->TreeR();
162 }
163
164 Int_t bufsize = 32000;
165 fTreeR->Branch("TOF", &fRecPoints, bufsize);
166
167 fTreeD->GetEvent(0);
168 Int_t nDigits = digits->GetEntriesFast();
169 AliDebug(2,Form("Number of TOF digits: %d",nDigits));
170
171 Int_t ii, jj;
172 Int_t dig[5];
173 Float_t g[3];
174 Double_t h[5];
175 for (ii=0; ii<nDigits; ii++) {
176 AliTOFdigit *d = (AliTOFdigit*)digits->UncheckedAt(ii);
177 dig[0]=d->GetSector();
178 dig[1]=d->GetPlate();
179 dig[2]=d->GetStrip();
180 dig[3]=d->GetPadz();
181 dig[4]=d->GetPadx();
182
183 for (jj=0; jj<3; jj++) g[jj] = 0.;
184 fTOFGeometry->GetPos(dig,g);
185
186 h[0] = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
187 h[1] = TMath::ATan2(g[1],g[0]);
188 h[2] = g[2];
189 h[3] = d->GetTdc();
190 h[4] = d->GetAdc();
191
192 AliTOFcluster *tofCluster = new AliTOFcluster(h,d->GetTracks(),dig,ii);
193 InsertCluster(tofCluster);
194
195 }
196
197 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
198
199 FillRecPoint();
200
201 fTreeR->Fill();
202 ResetRecpoint();
203
204 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
205 fTOFLoader->WriteRecPoints("OVERWRITE");
206
207}
208//______________________________________________________________________________
209
210void AliTOFClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
211 TTree *clustersTree)
212{
213 //
214 // Converts RAW data to recpoints for TOF
215 //
216
217 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
218
219 ResetRecpoint();
220
221 Int_t bufsize = 32000;
222 clustersTree->Branch("TOF", &fRecPoints, bufsize);
223
224 Int_t ii = 0;
225 Int_t indexDDL = 0;
226
227 Int_t detectorIndex[5];
228 Float_t position[3];
229 Double_t cylindricalPosition[5];
230
231 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
232
233 rawReader->Reset();
234 AliTOFRawStream tofInput(rawReader);
235 rawReader->Select(5, indexDDL, indexDDL);
236
237 while(tofInput.Next()) {
238
239 detectorIndex[0] = tofInput.GetSector();
240 detectorIndex[1] = tofInput.GetPlate();
241 detectorIndex[2] = tofInput.GetStrip();
242 detectorIndex[3] = tofInput.GetPadZ();
243 detectorIndex[4] = tofInput.GetPadX();
244
245 for (ii=0; ii<3; ii++) position[ii] = 0.;
246
247 fTOFGeometry->GetPos(detectorIndex, position);
248
249 cylindricalPosition[0] = TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
250 cylindricalPosition[1] = TMath::ATan2(position[1], position[0]);
251 cylindricalPosition[2] = position[2];
252 cylindricalPosition[3] = tofInput.GetTofBin();
253 cylindricalPosition[4] = tofInput.GetADCbin();
254
255 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
256 InsertCluster(tofCluster);
257
258 } // while loop
259
260 } // loop on DDL files
261
262 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
263
264 FillRecPoint();
265
266 clustersTree->Fill();
267 ResetRecpoint();
268
269}
270//______________________________________________________________________________
271
272void AliTOFClusterFinder::Digits2RecPoints(Int_t iEvent, AliRawReader *rawReader)
273{
274 //
275 // Converts RAW data to recpoints for TOF
276 //
277
278 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
279
280 fRunLoader->GetEvent(iEvent);
281
282 AliDebug(2,Form(" Event number %2i ", iEvent));
283
284 fTreeR = fTOFLoader->TreeR();
285
286 if (fTreeR == 0x0){
287 fTOFLoader->MakeTree("R");
288 fTreeR = fTOFLoader->TreeR();
289 }
290
291 Int_t bufsize = 32000;
292 fTreeR->Branch("TOF", &fRecPoints, bufsize);
293
294 Int_t ii = 0;
295 Int_t indexDDL = 0;
296
297 Int_t detectorIndex[5];
298 Float_t position[3];
299 Double_t cylindricalPosition[5];
300
301 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
302
303 rawReader->Reset();
304 AliTOFRawStream tofInput(rawReader);
305 rawReader->Select(5, indexDDL, indexDDL);
306
307 while(tofInput.Next()) {
308
309 detectorIndex[0] = (Int_t)tofInput.GetSector();
310 detectorIndex[1] = (Int_t)tofInput.GetPlate();
311 detectorIndex[2] = (Int_t)tofInput.GetStrip();
312 detectorIndex[3] = (Int_t)tofInput.GetPadZ();
313 detectorIndex[4] = (Int_t)tofInput.GetPadX();
314
315 for (ii=0; ii<3; ii++) position[ii] = 0.;
316
317 fTOFGeometry->GetPos(detectorIndex, position);
318
319 cylindricalPosition[0] = (Double_t)TMath::Sqrt(position[0]*position[0] + position[1]*position[1]);
320 cylindricalPosition[1] = (Double_t)TMath::ATan2(position[1], position[0]);
321 cylindricalPosition[2] = (Double_t)position[2];
322 cylindricalPosition[3] = (Double_t)tofInput.GetTofBin();
323 cylindricalPosition[4] = (Double_t)tofInput.GetADCbin();
324
325 AliTOFcluster *tofCluster = new AliTOFcluster(cylindricalPosition, detectorIndex);
326 InsertCluster(tofCluster);
327
328 } // while loop
329
330 } // DDL Loop
331
332 AliInfo(Form("Number of found clusters: %i", fNumberOfTofClusters));
333
334 FillRecPoint();
335
336 fTreeR->Fill();
337 ResetRecpoint();
338
339 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
340 fTOFLoader->WriteRecPoints("OVERWRITE");
341
342}
343//______________________________________________________________________________
344
345void AliTOFClusterFinder::Raw2Digits(Int_t iEvent, AliRawReader *rawReader)
346{
347 //
348 // Converts RAW data to MC digits for TOF
349 //
350 // (temporary solution)
351 //
352
353 fRunLoader->GetEvent(iEvent);
354
355 fTreeD = fTOFLoader->TreeD();
356 if (fTreeD)
357 {
358 AliInfo("AliTOFClusterFinder: TreeD re-creation");
359 fTreeD = 0x0;
360 fTOFLoader->MakeTree("D");
361 fTreeD = fTOFLoader->TreeD();
362 }
363
364
365 fTreeR = fTOFLoader->TreeD();
366 if (fTreeD == 0x0)
367 {
368 fTOFLoader->MakeTree("D");
369 fTreeD = fTOFLoader->TreeD();
370 }
371
372 TClonesArray dummy("AliTOFdigit",10000), *tofDigits=&dummy;
373 Int_t bufsize = 32000;
374 fTreeD->Branch("TOF", &tofDigits, bufsize);
375
376
377 const Int_t kDDL = fTOFGeometry->NDDL()*fTOFGeometry->NSectors();
378
379 fRunLoader->GetEvent(iEvent);
380
381 AliDebug(2,Form(" Event number %2i ", iEvent));
382
383 Int_t indexDDL = 0;
384
385 Int_t detectorIndex[5];
386 Float_t digit[2];
387
388 for (indexDDL = 0; indexDDL < kDDL; indexDDL++) {
389
390 rawReader->Reset();
391 AliTOFRawStream tofInput(rawReader);
392 rawReader->Select(5, indexDDL, indexDDL);
393
394 while(tofInput.Next()) {
395
396 detectorIndex[0] = tofInput.GetSector();
397 detectorIndex[1] = tofInput.GetPlate();
398 detectorIndex[2] = tofInput.GetStrip();
399 detectorIndex[3] = tofInput.GetPadX();
400 detectorIndex[4] = tofInput.GetPadZ();
401
402 digit[0] = (Float_t)tofInput.GetTofBin();
403 digit[1] = (Float_t)tofInput.GetADCbin();
404
405 Int_t tracknum[3]={-1,-1,-1};
406
407 TClonesArray &aDigits = *tofDigits;
408 Int_t last=tofDigits->GetEntriesFast();
409 new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
410
411 } // while loop
412
413 } // DDL Loop
414
415 fTreeD->Fill();
416
417 fTOFLoader = fRunLoader->GetLoader("TOFLoader");
418 fTOFLoader->WriteDigits("OVERWRITE");
419
420}
421//______________________________________________________________________________
422
423Int_t AliTOFClusterFinder::InsertCluster(AliTOFcluster *tofCluster) {
424 //---------------------------------------------------------------------------//
425 // This function adds a TOF cluster to the array of TOF clusters sorted in Z //
426 //---------------------------------------------------------------------------//
427 if (fNumberOfTofClusters==kTofMaxCluster) {
428 AliError("Too many clusters !");
429 return 1;
430 }
431
432 if (fNumberOfTofClusters==0) {
433 fTofClusters[fNumberOfTofClusters++] = tofCluster;
434 return 0;
435 }
436
437 Int_t ii = FindClusterIndex(tofCluster->GetZ());
438 memmove(fTofClusters+ii+1 ,fTofClusters+ii,(fNumberOfTofClusters-ii)*sizeof(AliTOFcluster*));
439 fTofClusters[ii] = tofCluster;
440 fNumberOfTofClusters++;
441
442 return 0;
443
444}
445//_________________________________________________________________________
446
447Int_t AliTOFClusterFinder::FindClusterIndex(Double_t z) const {
448 //--------------------------------------------------------------------
449 // This function returns the index of the nearest cluster
450 //--------------------------------------------------------------------
451 if (fNumberOfTofClusters==0) return 0;
452 if (z <= fTofClusters[0]->GetZ()) return 0;
453 if (z > fTofClusters[fNumberOfTofClusters-1]->GetZ()) return fNumberOfTofClusters;
454 Int_t b = 0, e = fNumberOfTofClusters-1, m = (b+e)/2;
455 for (; b<e; m=(b+e)/2) {
456 if (z > fTofClusters[m]->GetZ()) b=m+1;
457 else e=m;
458 }
459
460 return m;
461
462}
463//_________________________________________________________________________
464
465void AliTOFClusterFinder::FillRecPoint()
466{
467 //
468 // Copy the global array of AliTOFcluster, i.e. fTofClusters (sorted
469 // in Z) in the global TClonesArray of AliTOFcluster,
470 // i.e. fRecPoints.
471 //
472
473 Int_t ii, jj;
474
475 Int_t detectorIndex[5];
476 Double_t cylindricalPosition[5];
477 Int_t trackLabels[3];
478 Int_t digitIndex = -1;
479
480 TClonesArray &lRecPoints = *fRecPoints;
481
482 for (ii=0; ii<fNumberOfTofClusters; ii++) {
483
484 digitIndex = fTofClusters[ii]->GetIndex();
485 for(jj=0; jj<5; jj++) detectorIndex[jj] = fTofClusters[ii]->GetDetInd(jj);
486 for(jj=0; jj<3; jj++) trackLabels[jj] = fTofClusters[ii]->GetLabel(jj);
487 cylindricalPosition[0] = fTofClusters[ii]->GetR();
488 cylindricalPosition[1] = fTofClusters[ii]->GetPhi();
489 cylindricalPosition[2] = fTofClusters[ii]->GetZ();
490 cylindricalPosition[3] = fTofClusters[ii]->GetTDC();
491 cylindricalPosition[4] = fTofClusters[ii]->GetADC();
492
493 new(lRecPoints[ii]) AliTOFcluster(cylindricalPosition, trackLabels, detectorIndex, digitIndex);
494
495 //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]));
496
497 } // loop on clusters
498
499}
500//______________________________________________________________________________
501
502void AliTOFClusterFinder::ResetRecpoint()
503{
504 //
505 // Clear the list of reconstructed points
506 //
507
508 fNumberOfTofClusters = 0;
509 if (fRecPoints) fRecPoints->Clear();
510
511}
512//______________________________________________________________________________
513
514void AliTOFClusterFinder::Load()
515{
516 //
517 // Load TOF.Digits.root and TOF.RecPoints.root files
518 //
519
520 fTOFLoader->LoadDigits("READ");
521 fTOFLoader->LoadRecPoints("recreate");
522
523}
524//______________________________________________________________________________
525
526void AliTOFClusterFinder::LoadClusters()
527{
528 //
529 // Load TOF.RecPoints.root file
530 //
531
532 fTOFLoader->LoadRecPoints("recreate");
533
534}
535//______________________________________________________________________________
536
537void AliTOFClusterFinder::UnLoad()
538{
539 //
540 // Unload TOF.Digits.root and TOF.RecPoints.root files
541 //
542
543 fTOFLoader->UnloadDigits();
544 fTOFLoader->UnloadRecPoints();
545
546}
547//______________________________________________________________________________
548
549void AliTOFClusterFinder::UnLoadClusters()
550{
551 //
552 // Unload TOF.RecPoints.root file
553 //
554
555 fTOFLoader->UnloadRecPoints();
556
557}
558//______________________________________________________________________________