]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TOF/AliTOFReconstructionerV2.cxx
Path to ROOT and AliRoot libraries added (by F. Carminati)
[u/mrichter/AliRoot.git] / TOF / AliTOFReconstructionerV2.cxx
CommitLineData
bf6bf84c 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// This is a TTask for reconstruction V2 in TOF
18// Description of the algorithm
19//-- Author: F. Pierella | pierella@bo.infn.it
20//////////////////////////////////////////////////////////////////////////////
21
22#include "TBenchmark.h"
23#include "TFile.h"
24#include "TFolder.h"
25#include "TParticle.h"
26#include "TTask.h"
27#include "TTree.h"
28#include "TClonesArray.h"
29#include "TROOT.h"
30#include "TSystem.h"
31#include <stdlib.h>
32#include <iostream.h>
33#include <fstream.h>
34#include "../TGeant3/TGeant3.h"
35
36#include "AliDetector.h"
37#include "AliMC.h"
38#include "AliRun.h"
39#include "AliTOF.h"
40#include "AliTOFDigitMap.h"
41#include "AliTOFHitMap.h"
42#include "AliTOFhitT0.h"
43#include "AliTOFdigit.h"
44#include "AliTOFConstants.h"
45#include "AliTOFReconstructionerV2.h"
46#include "AliTOFTrackV2.h"
47
48#include "../STEER/AliKalmanTrack.h"
49#include "../TPC/AliTPCtrack.h"
50#include "../TRD/AliTRDtrack.h"
51
52ClassImp(AliTOFReconstructionerV2)
53
54//____________________________________________________________________________
55 AliTOFReconstructionerV2::AliTOFReconstructionerV2():TTask("AliTOFReconstructionerV2","")
56{
57 //
58 // std ctor
59 // set all member vars to zero
60 //
61 fdbg =0;
62 fDigitsMap=0x0;
63 fField =0;
64 fG3Geom =0x0;
65 fNDummyTracks=0;
66 fScaleSigmaFactor=0.;
67 fStep =0;
68 fTOFDigits=0x0;
69 fTOFTracks=0x0;
70 fTOFDigitsFile ="digits.root";
71 fTPCBackTracksFile="AliTPCBackTracks.root";
72 fKalmanTree =0x0;
73 fBranchWithTracks =0x0;
74}
75
76//____________________________________________________________________________
77AliTOFReconstructionerV2::AliTOFReconstructionerV2(char* tpcBackTracks, char* tofDigits):TTask("AliTOFReconstructionerV2","")
78{
79 //
80 // par ctor
81 // default arguments are specified only in
82 // the header file
83 // Parameters:
84 // tpcBackTracks -> file name with backpropagated tracks in TPC
85 // tofDigits -> file name with TOF digits
86 //
87
88 fdbg =0;
89 fDigitsMap=0x0;
90 fField =0.2; // default value 0.2 [Tesla]
91 fNDummyTracks=20; // by default 20 test tracks used
92 fScaleSigmaFactor=1.;
93 fStep =0.01; // [cm]
94 fTOFDigits=0x0;
95 fTOFTracks=0x0;
96 fTOFDigitsFile =tofDigits;
97 fTPCBackTracksFile=tpcBackTracks;
98 fKalmanTree =0x0;
99 fBranchWithTracks =0x0;
100
101 // initialize the G3 geometry
102 gAlice->Init();
103 gAlice->Print();
104 // set the fg3 pointer to geometry used by IsInsideThePad method
105 fG3Geom = (TGeant3*) gMC;
106
107 // add Task to //root/Tasks folder
108 TTask * roottasks = (TTask*)gROOT->GetRootFolder()->FindObject("Tasks") ;
109 roottasks->Add(this) ;
110}
111
112//____________________________________________________________________________
113AliTOFReconstructionerV2::AliTOFReconstructionerV2(const AliTOFReconstructionerV2 & rec)
114{
115 //
116 // Dummy copy constructor
117 // required by coding conventions
118 ;
119}
120
121
122
123//____________________________________________________________________________
124 AliTOFReconstructionerV2::~AliTOFReconstructionerV2()
125{
126 //
127 // dtor
128 // some delete has to be moved
129 //
130
131 if (fDigitsMap)
132 {
133 delete fDigitsMap;
134 fDigitsMap = 0;
135 }
136
137 if (fG3Geom)
138 {
139 delete fG3Geom;
140 fG3Geom = 0;
141 }
142
143 if (fTOFDigits)
144 {
145 delete fTOFDigits;
146 fTOFDigits = 0;
147 }
148
149 if (fTOFTracks)
150 {
151 delete fTOFTracks;
152 fTOFTracks = 0;
153 }
154
155 if (fTOFDigitsFile)
156 {
157 delete fTOFDigitsFile;
158 }
159
160 if (fTPCBackTracksFile)
161 {
162 delete fTPCBackTracksFile;
163 }
164
165 if (fKalmanTree)
166 {
167 delete fKalmanTree;
168 fKalmanTree = 0;
169 }
170
171 if (fBranchWithTracks)
172 {
173 delete fBranchWithTracks;
174 fBranchWithTracks = 0;
175 }
176}
177
178
179//____________________________________________________________________________
180void AliTOFReconstructionerV2::Exec(Option_t* option)
181{
182 //
183 // Description of the algorithm:
184 //
185 //
186 //
187 //
188 //
189 //
190
191 // load TOF digits and fill the digit map
192 Int_t tofDigitsLoading=LoadTOFDigits();
193
194 // load back-propagated tracks in TPC
195 Int_t tpcTracksLoading=LoadTPCTracks();
196
197 if(tofDigitsLoading || tpcTracksLoading) {
198 cout<<" Couldn't start reconstruction V2. Exit."<<endl;
199 exit(1);
200 }
201 // or load TRD reconstructed tracks
202 // Int_t trdTracksLoading=LoadTRDTracks();
203
204 // create a TObjArray to store reconstructed tracks
205 // and reject fake tracks
206 const Int_t maxRecTracks = 10000; // max number of reconstructed tracks
207 // per event
208
209 const Float_t stripRegionHeight = 15.3; // [cm] height in radial direction
210 // of the volume where strips are placed
211 TObjArray trackArray(maxRecTracks);
212
213 // buffer for output
214 fTOFTracks= new TClonesArray("AliTOFTrackV2");
215 // create a reference to fill the TClonesArray
216 TClonesArray &aTOFTracks = *fTOFTracks;
217
218 const Int_t maxIndex = 100000; // max number of primary tracks to be analysed
219 // the index of the rtIndex array is the track label
220 // the content of the array is -1 for fake tracks
221 // and the track index for good tracks
222 Float_t dEdXarray[maxRecTracks];
223 Int_t rtIndex[maxIndex];
224 for(Int_t i = 0; i < maxIndex; i++) rtIndex[i] = -1;
225
226 AliKalmanTrack::SetConvConst(100/0.299792458/fField);
227
228 Int_t nRecTracks = (Int_t) fKalmanTree->GetEntries();
229 cout<<"Found "<<nRecTracks<<" entries in the track tree "<<endl;
230
231 // load the tracks into the array
232 for (Int_t i=0; i<nRecTracks; i++) {
233 AliTPCtrack *iotrack=new AliTPCtrack();
234 fBranchWithTracks->SetAddress(&iotrack);
235 fKalmanTree->GetEvent(i);
236 trackArray.AddLast(iotrack);
237 Int_t trackLabel = iotrack->GetLabel();
238 dEdXarray[i]=iotrack->GetdEdx(); // usefull for PID
239 //
240 // start filling the TClonesArray of AliTOFTrackV2
241 Float_t trdXYZ[3]={0.,0.,0.};
242 Float_t trdPxPyPz[3]={0.,0.,0.};
243
244 // tpc outer wall positions
245 Double_t xk=iotrack->GetX();
246 // get the running coordinates in the lrf
247 // and alpha
248 Double_t x=xk;
249 Double_t y=iotrack->GetY();
250 Double_t z=iotrack->GetZ();
251 Double_t alpha=iotrack->GetAlpha();
252 GetGlobalXYZ(alpha, x, y, z);
253 Float_t tpcXYZ[3]={x,y,z};
254
255 // momentum at the end of TPC
256 Float_t lambda=TMath::ATan(iotrack->GetTgl());
257 Float_t invpt=TMath::Abs(iotrack->Get1Pt());
258 Float_t pt=-99.;
259 if (invpt) pt= 1./invpt; // pt
260 Float_t tpcmom=1./(invpt*TMath::Cos(lambda));
261 Float_t pz=tpcmom*TMath::Sin(lambda);
262 Float_t tpcPtPz[2]={pt,pz};
263
264 Int_t matchingStatus=-1;
265 if(trackLabel < 0) matchingStatus=0;
266 new(aTOFTracks[i]) AliTOFTrackV2(trackLabel,matchingStatus,tpcmom,dEdXarray[i],tpcXYZ,tpcPtPz,trdXYZ,trdPxPyPz);
267 // printf("rt with %d clusters and label %d \n",
268 // iotrack->GetNumberOfClusters(), trackLabel);
269
270 if(trackLabel < 0) continue;
271 if(trackLabel >= maxIndex) continue;
272 rtIndex[trackLabel] = i;
273 delete iotrack;
274 }
275
276 if(strstr(option,"MC")) Comparison(rtIndex);
277
278 // start loop on tracks
279 // and backpropagate them from TPC to TOF
280 // backpropagation is performed only
281 // for good tracks (fake tracks rejected)
282 AliTPCtrack *rt;
283
284 for (Int_t i=0; i<nRecTracks; i++) {
285
286 //******* tracking: extract track coordinates, momentum, etc.
287 rt = (AliTPCtrack*) trackArray.UncheckedAt(i);
288 // track length to be implemented
289 // Double_t tr_length = rt->GetLength();
290
291 Int_t tpcTrackLabel=rt->GetLabel();
292 // reject fake tracks
293 //if(tpcTrackLabel< 0) continue;
294
295 // starting backpropagation to TOF
296 // here we assume to have backpropagated tracks in TPC
297 // starting point xk=246.055
298 // starting back propagation
299 // outer wall of the TPC
300 Int_t iOuterTPCWall=rt->PropagateTo(261.53,40.,0.06124);
301 // frame with air just to the beginning of the TOF
302 Int_t iFrameWithAir=rt->PropagateTo(370.,36.66,1.2931e-3);
303 // trough the wall of the TOF plate
304 // thickness has changed according to the
305 // last value
306 Int_t iTOFWall=rt->PropagateTo(370.11,24.01,2.7);
307
308 /*
309 // outer wall of the TPC
310 Int_t iOuterTPCWall=rt->PropagateTo(261.53,40.,0.06124);
311 // air in between TPC and TRD
312 Int_t iFrameWithAir=rt->PropagateTo(294.5,36.66,1.2931e-3);
313 // TRD itself
314 // mean density for the TRD calculated from
315 // TRD Technical Design Report
316 // page 11 -> global thickness
317 // page 23 -> different main layers thickness (Radiator Air/ Drift Chamber Gas /G10)
318 // page 139 -> material budget and radiation lengths
319 Int_t iTRD=rt->PropagateTo(369.1,171.7,0.33);
320 // air in between TRD and TOF
321 Int_t iFrameWithAirbis=rt->PropagateTo(370.,36.66,1.2931e-3);
322 // trough the wall of the TOF plate
323 Int_t iTOFWall=rt->PropagateTo(370.11,24.01,2.7);
324 */
325
326 // select only cases when
327 // backpropagation succeded
328 // and particle is in the geometrical TOF acceptance along Z
329 // |z|<380 [cm]
330 AliTOFTrackV2* oTOFtracks=(AliTOFTrackV2*)fTOFTracks->UncheckedAt(i);
331 Bool_t outOfZacceptance=(rt->GetZ()<=380.);
332 if(outOfZacceptance) oTOFtracks->SetMatchingStatus(-2);
333
334 if(iOuterTPCWall==1 && iFrameWithAir==1 && iTOFWall==1 && (!outOfZacceptance)){
335 Double_t cc[15];
336 // get sigmaY and sigmaZ
337 rt->GetExternalCovariance(cc);
338 Double_t sigmaY =TMath::Sqrt(cc[0]); // [cm]
339 Double_t sigmaZ =TMath::Sqrt(cc[2]); // [cm]
340
341 // arrays used by the DigitFinder
342 Int_t nSlot=1;
343 TArrayI *secArray= new TArrayI(nSlot);
344 TArrayI *plaArray= new TArrayI(nSlot);
345 TArrayI *strArray= new TArrayI(nSlot);
346 TArrayI *pdzArray= new TArrayI(nSlot);
347 TArrayI *pdxArray= new TArrayI(nSlot);
348
349 // make fNDummyTracks clones of the current backpropagated track
350 // make a copy of the current track
351 // smear according to the backpropagated area
352 for (Int_t j=0; j<fNDummyTracks; i++) {
353 AliTPCtrack *dummyrt=new AliTPCtrack(*rt);
354 // get ylrf and zlrf
355 Double_t ylrf= dummyrt->GetY(); // P0
356 Double_t zlrf= dummyrt->GetZ(); // P1
357
358 // smear according to sigmaY and sigmaZ
359 Double_t ylrfNew=gRandom->Gaus(ylrf,fScaleSigmaFactor*sigmaY);
360 Double_t zlrfNew=gRandom->Gaus(zlrf,fScaleSigmaFactor*sigmaZ);
361
362 // set Y and Z accordingly
363 // setter to be added in the class AliTPCtrack
364 // here I have to modify the AliTPCtrack class
365 // adding the setters for Y and Z
366 //dummyrt->SetY(ylrfNew);
367 //dummyrt->SetZ(zlrfNew);
368
369 // start fine-backpropagation inside the TOF
370 Bool_t padNotFound =kTRUE;
371 Bool_t isInStripsRegion=kTRUE;
372 Double_t xk=dummyrt->GetX();
373
374 while (padNotFound && isInStripsRegion){
375 xk+=fStep;
376 // here we assume a frame with air
377 dummyrt->PropagateTo(xk,36.66,1.2931e-3);
378 // get the running coordinates in the lrf
379 // and alpha
380 Double_t x=xk;
381 Double_t y=dummyrt->GetY();
382 Double_t z=dummyrt->GetZ();
383 Double_t alpha=dummyrt->GetAlpha();
384 GetGlobalXYZ(alpha, x, y, z);
385
386 // check if the point falls into a pad
387 // using the G3 geometry
388 Int_t* volumeID = new Int_t[AliTOFConstants::fgkmaxtoftree];
389 // volumeID[0] -> TOF Sector range [1-18]
390 // volumeID[1] -> TOF Plate range [1- 5]
391 // volumeID[2] -> TOF Strip max range [1-20]
392 // volumeID[3] -> TOF Pad along Z range [1- 2]
393 // volumeID[4] -> TOF Pad along X range [1-48]
394
395 Float_t zInPadFrame=0.;
396 Float_t xInPadFrame=0.;
397 IsInsideThePad((Float_t)x,(Float_t)y,(Float_t)z, volumeID, zInPadFrame, xInPadFrame);
398 // adding protection versus wrong volume numbering
399 // to be released in the next release after debugging
400 if(volumeID[4]){
401 padNotFound=kFALSE;
402 nSlot++;
403 secArray->Set(nSlot-1);
404 plaArray->Set(nSlot-1);
405 strArray->Set(nSlot-1);
406 pdzArray->Set(nSlot-1);
407 pdxArray->Set(nSlot-1);
408
409 (*secArray)[nSlot-1]=volumeID[0];
410 (*plaArray)[nSlot-1]=volumeID[1];
411 (*strArray)[nSlot-1]=volumeID[2];
412 (*pdzArray)[nSlot-1]=volumeID[3];
413 (*pdxArray)[nSlot-1]=volumeID[4];
414
415 } // track falls into a pad volume
416
417 delete [] volumeID;
418
419 // check on xk to stop the fine-propagation
420 if(xk>=(370.+stripRegionHeight)) isInStripsRegion=kFALSE;
421
422 } // close the while for fine-propagation
423
424 delete dummyrt;
425 } // end loop on test tracks
426
427 // start TOF digit finder
428 Int_t assignedVol[5]={0,0,0,0,0};
429 Float_t tdc=-1.;
430 Int_t* digitTrackArray=0x0;
431 Bool_t assignedDigit=DigitFinder(secArray, plaArray, strArray, pdzArray, pdxArray, assignedVol, digitTrackArray, tdc);
432
433
434 if(assignedDigit){
435 // fill the tree for tracks with time of flight
436 // tof is given in tdc bin
437 // conversion to [ns]
438 Float_t binWidth=50.; // [ps]
439 Float_t timeOfFlight=tdc*binWidth/1000.;
440
441 // only the first track number contributing
442 // to the assigned digit
443 Int_t tofDigitTrackLabel=digitTrackArray[0];
444
445 // matching status for the current track
446 Int_t matching=4;
447 if(tpcTrackLabel==digitTrackArray[0] || tpcTrackLabel==digitTrackArray[1] || tpcTrackLabel==digitTrackArray[0]) matching=3;
448 oTOFtracks->UpdateTrack(tofDigitTrackLabel, matching, timeOfFlight);
449 } else {
450 // fill the TClonesArray for tracks with no time of flight
451 Int_t matching=1;
452 oTOFtracks->SetMatchingStatus(matching);
453 }
454
455 // delete used memory for tmp arrays used by DigitFinder
456 delete secArray;
457 delete plaArray;
458 delete strArray;
459 delete pdzArray;
460 delete pdxArray;
461
462 } // close the if for succeded backpropagation in TOF acceptance along z
463
464 } // end loop on reconstructed tracks
465
466 // free used memory for digitmap
467 delete fDigitsMap;
468
469 // save array with TOF tracks
470 Int_t output=SaveTracks();
471 if(output) cout << "Error writing TOF tracks " << endl;
472}
473
474
475//__________________________________________________________________
476void AliTOFReconstructionerV2::Init(Option_t* opt)
477{
478 //
479 // Initialize the AliTOFReconstructionerV2
480 //
481 //
482}
483
484
485//__________________________________________________________________
486Int_t AliTOFReconstructionerV2::LoadTPCTracks()
487{
488 //
489 // Connect the tree and the branch
490 // with reconstructed tracks
491 // backpropagated in the TPC
492
493 gBenchmark->Start("LoadTPCTracks");
494
495 TFile *kalFile = TFile::Open(fTPCBackTracksFile.Data());
496 if (!kalFile->IsOpen()) {cerr<<"Can't open AliTPCBackTracks.root !\n"; return 3;}
497
498 // tracks from Kalman
499 Int_t event=0;
500 char treename[100]; sprintf(treename,"TreeT_TPCb_%d",event);
501 fKalmanTree=(TTree*)kalFile->Get(treename);
502 if (!fKalmanTree) {cerr<<"Can't get a tree with TPC back tracks !\n"; return 4;}
503
504 // otherwise you get always 0 for 1/pt
505 AliKalmanTrack::SetConvConst(100/0.299792458/fField);
506
507 fBranchWithTracks=fKalmanTree->GetBranch("tracks");
508 Int_t kalEntries =(Int_t)fKalmanTree->GetEntries();
509 cout<<"Number of loaded Tracks :"<< kalEntries <<endl;
510
511 gBenchmark->Stop("LoadTPCTracks");
512 gBenchmark->Show("LoadTPCTracks");
513 return 0;
514}
515
516//__________________________________________________________________
517Int_t AliTOFReconstructionerV2::LoadTRDTracks()
518{
519 //
520 // Connect the tree and the branch
521 // with reconstructed tracks in TRD
522
523 gBenchmark->Start("LoadTRDTracks");
524
525 Int_t nEvent = 0;
526 const Int_t nPrimaries = 84210/16;
527 const Int_t maxIndex = nPrimaries;
528 Int_t rtIndex[maxIndex];
529
530 TFile *tf=TFile::Open("AliTRDtracks.root");
531
532 if (!tf->IsOpen()) {cerr<<"Can't open AliTRDtracks.root !\n"; return 3;}
533 TObjArray tarray(2000);
534 char tname[100];
535 sprintf(tname,"TRDb_%d",nEvent);
536 TTree *tracktree=(TTree*)tf->Get(tname);
537
538 TBranch *tbranch=tracktree->GetBranch("tracks");
539
540 Int_t nRecTracks = (Int_t) tracktree->GetEntries();
541 cerr<<"Found "<<nRecTracks<<" entries in the track tree"<<endl;
542
543 for (Int_t i=0; i<nRecTracks; i++) {
544 AliTRDtrack *iotrack=new AliTRDtrack();
545 tbranch->SetAddress(&iotrack);
546 tracktree->GetEvent(i);
547 tarray.AddLast(iotrack);
548 Int_t trackLabel = iotrack->GetLabel();
549
550 // printf("rt with %d clusters and label %d \n",
551 // iotrack->GetNumberOfClusters(), trackLabel);
552
553 if(trackLabel < 0) continue;
554 if(trackLabel >= maxIndex) continue;
555 rtIndex[trackLabel] = i;
556 }
557 tf->Close();
558 gBenchmark->Stop("LoadTRDTracks");
559 gBenchmark->Show("LoadTRDTracks");
560 return 0;
561
562}
563
564//__________________________________________________________________
565Int_t AliTOFReconstructionerV2::LoadTOFDigits()
566{
567 //
568 // Connect the TClonesArray with TOF
569 // digits and fill the digit map
570 // used by the DigitFinder
571
572 Int_t rc=0;
573
574 gBenchmark->Start("LoadTOFDigits");
575
576 TFile *file = (TFile*)gROOT->GetListOfFiles()->FindObject(fTOFDigitsFile.Data());
577 if(file){
578 cout<<"headerFile already open \n";
579 }
580 else {
581 if(!file)file=TFile::Open(fTOFDigitsFile.Data());
582 }
583
584 // Get AliRun object from file
585 if (!gAlice) {
586 gAlice = (AliRun*)file->Get("gAlice");
587 if (gAlice) printf("AliRun object found on file\n");
588 }
589
590 Int_t iEvNum = 0;
591 if (iEvNum == 0) iEvNum = (Int_t) gAlice->TreeE()->GetEntries();
592
593
594 AliTOFdigit *tofdigit;
595
596 AliTOF * tof = (AliTOF *) gAlice->GetDetector("TOF") ;
597
598 if (!tof) {
599 cout << "<LoadTOFDigits> No TOF detector found" << endl;
600 rc = 2;
601 return rc;
602 }
603
604 for (Int_t ievent = 0; ievent < iEvNum; ievent++) {
605
606 gAlice->GetEvent(ievent) ;
607 if(gAlice->TreeD()==0) {
608 cout << "<LoadTOFDigits> No TreeD found" << endl;
609 rc = 4;
610 return rc;
611 }
612
613
614 Int_t ndig;
615 gAlice->ResetDigits();
616 gAlice->TreeD()->GetEvent(ievent);
617 fTOFDigits = tof->Digits();
618
619 ndig=fTOFDigits->GetEntries();
620
621 // create the digit map
622 fDigitsMap = new AliTOFDigitMap(fTOFDigits);
623
624
625 cout << "<LoadTOFDigits> found " << ndig
626 << " TOF digits for event " << ievent << endl;
627
628 for (Int_t k=0; k<ndig; k++) {
629 tofdigit= (AliTOFdigit*) fTOFDigits->UncheckedAt(k);
630 Float_t tdc=tofdigit->GetTdc();
631 // adc value can be used for weighting
632 //Float_t adc=tofdigit->GetAdc();
633
634 // TOF digit volumes
635 Int_t vol[5]; // location for a digit
636 Int_t sector = tofdigit->GetSector(); // range [1-18]
637 Int_t plate = tofdigit->GetPlate(); // range [1- 5]
638 Int_t strip = tofdigit->GetStrip(); // range [1-20]
639 Int_t padx = tofdigit->GetPadx(); // range [1-48]
640 Int_t padz = tofdigit->GetPadz(); // range [1- 2]
641
642 vol[0] = sector;
643 vol[1] = plate;
644 vol[2] = strip;
645 vol[3] = padx;
646 vol[4] = padz;
647
648 // QA
649 Bool_t isDigitBad = (sector<1 || sector>18 || plate<1 || plate >5 || padz<1 || padz>2 || padx<1 || padx>48);
650
651 if (isDigitBad) {
652 cout << "<LoadTOFDigits> strange digit found" << endl;
653 rc = 3;
654 return rc;
655 } // if (isDigitBad)
656
657 // Fill the digit map checking if the location is already used
658 // in this case we take the earliest signal
659 if (fDigitsMap->TestHit(vol) != kEmpty) {
660 // start comparison in between the 2 digit
661 AliTOFdigit *dig = static_cast<AliTOFdigit*>(fDigitsMap->GetHit(vol));
662 if(tdc < (dig->GetTdc())) fDigitsMap->SetHit(vol,k);
663 // we can add also the check on adc value
664 // by selecting the largest adc value
665 } else {
666 fDigitsMap->SetHit(vol,k);
667 }
668 // to be added protection versus 2-digit on the same pad
669 // we have to have memory also of the second digit
670
671 } // for (k=0; k<ndig; k++)
672
673 } // end loop on events
674
675 gBenchmark->Stop("LoadTOFDigits");
676 gBenchmark->Show("LoadTOFDigits");
677 return rc;
678}
679
680//__________________________________________________________________
681void AliTOFReconstructionerV2::IsInsideThePad(Float_t x, Float_t y, Float_t z, Int_t *nGeom, Float_t& zPad, Float_t& xPad)
682{
683 // input: x,y,z - coordinates of a point in the mrf [cm]
684 // output: array nGeom[]
685 // nGeom[0] - the TOF sector number, 1,2,...,18 along azimuthal direction starting from -90 deg.
686 // nGeom[1] - the TOF module number, 1,2,3,4,5=C,B,A,B,C along z-direction
687 // nGeom[2] - the TOF strip number, 1,2,... along z-direction
688 // nGeom[3] - the TOF padz number, 1,2=NPZ across a strip
689 // nGeom[4] - the TOF padx number, 1,2,...,48=NPX along a strip
690 // zPad, xPad - coordinates of the hit in the pad frame
691 // numbering is adopted for the version 3.08 of AliRoot
692 // example:
693 // from Hits: sec,pla,str,padz,padx=4,2,14,2,35
694 // Vol. n.0: ALIC, copy number 1
695 // Vol. n.1: B077, copy number 1
696 // Vol. n.2: B074, copy number 5
697 // Vol. n.3: BTO2, copy number 1
698 // Vol. n.4: FTOB, copy number 2
699 // Vol. n.5: FLTB, copy number 0
700 // Vol. n.6: FSTR, copy number 14
701 // Vol. n.7: FSEN, copy number 0
702 // Vol. n.8: FSEZ, copy number 2
703 // Vol. n.9: FSEX, copy number 35
704 // Vol. n.10: FPAD, copy number 0
705
706
707 Float_t xTOF[3];
708 Int_t sector=0,module=0,strip=0,padz=0,padx=0;
709 Int_t i,numed,nLevel,copyNumber;
710 Gcvolu_t* gcvolu;
711 char name[5];
712 name[4]=0;
713
714 for (i=0; i<AliTOFConstants::fgkmaxtoftree; i++) nGeom[i]=0;
715 zPad=100.;
716 xPad=100.;
717
718 xTOF[0]=x;
719 xTOF[1]=y;
720 xTOF[2]=z;
721
722 fG3Geom->Gmedia(xTOF, numed);
723 gcvolu=fG3Geom->Gcvolu();
724 nLevel=gcvolu->nlevel;
725 if(fdbg) {
726 for (Int_t i=0; i<nLevel; i++) {
727 strncpy(name,(char*) (&gcvolu->names[i]),4);
728 cout<<"Vol. n."<<i<<": "<<name<<", copy number "<<gcvolu->number[i]<<endl;
729 }
730 }
731 if(nLevel>=2) {
732 // sector type name: B071(1,2,...,10),B074(1,2,3,4,5-PHOS),B075(1,2,3-RICH)
733 strncpy(name,(char*) (&gcvolu->names[2]),4);
734 // volume copy: 1,2,...,10 for B071, 1,2,3,4,5 for B074, 1,2,3 for B075
735 copyNumber=gcvolu->number[2];
736 if(!strcmp(name,"B071")) {
737 if (copyNumber>=6 && copyNumber<=8) {
738 sector=copyNumber+10;
739 } else if (copyNumber>=1 && copyNumber<=5){
740 sector=copyNumber+7;
741 } else {
742 sector=copyNumber-8;
743 }
744 } else if(!strcmp(name,"B075")) {
745 sector=copyNumber+12;
746 } else if(!strcmp(name,"B074")) {
747 if (copyNumber>=1 && copyNumber<=3){
748 sector=copyNumber+4;
749 } else {
750 sector=copyNumber-1;
751 }
752 }
753 }
754 if(sector) {
755 nGeom[0]=sector;
756 if(nLevel>=4) {
757 // we'll use the module value in z-direction:
758 // 1 2 3 4 5
759 // the module order in z-direction: FTOC,FTOB,FTOA,FTOB,FTOC
760 // the module copy: 2 2 0 1 1
761 // module type name: FTOA, FTOB, FTOC
762 strncpy(name,(char*) (&gcvolu->names[4]),4);
763 // module copy:
764 copyNumber=gcvolu->number[4];
765 if(!strcmp(name,"FTOC")) {
766 if (copyNumber==2) {
767 module=1;
768 } else {
769 module=5;
770 }
771 } else if(!strcmp(name,"FTOB")) {
772 if (copyNumber==2) {
773 module=2;
774 } else {
775 module=4;
776 }
777 } else if(!strcmp(name,"FTOA")) {
778 module=3;
779 }
780 }
781 }
782
783 if(module) {
784 nGeom[1]=module;
785 if(nLevel>=6) {
786 // strip type name: FSTR
787 strncpy(name,(char*) (&gcvolu->names[6]),4);
788 // strip copy:
789 copyNumber=gcvolu->number[6];
790 if(!strcmp(name,"FSTR")) strip=copyNumber;
791 }
792 }
793
794 if(strip) {
795 nGeom[2]=strip;
796 if(nLevel>=8) {
797 // padz type name: FSEZ
798 strncpy(name,(char*) (&gcvolu->names[8]),4);
799 // padz copy:
800 copyNumber=gcvolu->number[8];
801 if(!strcmp(name,"FSEZ")) padz=copyNumber;
802 }
803 }
804 if(padz) {
805 nGeom[3]=padz;
806 if(nLevel>=9) {
807 // padx type name: FSEX
808 strncpy(name,(char*) (&gcvolu->names[9]),4);
809 // padx copy:
810 copyNumber=gcvolu->number[9];
811 if(!strcmp(name,"FSEX")) padx=copyNumber;
812 }
813 }
814
815 if(padx) {
816 nGeom[4]=padx;
817 zPad=gcvolu->glx[2];
818 xPad=gcvolu->glx[0];
819 }
820
821}
822
823//__________________________________________________________________
824void AliTOFReconstructionerV2::GetGlobalXYZ(Double_t alpha, Double_t& x, Double_t& y, Double_t& z)
825{
826 //
827 // return the current running coordinates of
828 // the track in the global reference frame
829 // x, y and z have to initialized to the
830 // local frame coordinates by the caller
831 // alpha is the alpha coordinate in the TPC Kalman
832 // reference frame
833
834 // it take into account differences in between
835 // TPC and TRD local coordinates frames
836 if (alpha < -TMath::Pi()) alpha += 2*TMath::Pi();
837 else if (alpha >= TMath::Pi()) alpha -= 2*TMath::Pi();
838
839 // conversion
840 Double_t tmp=x*TMath::Cos(alpha) - y*TMath::Sin(alpha);
841 y=x*TMath::Sin(alpha) + y*TMath::Cos(alpha);
842 x=tmp;
843}
844
845//__________________________________________________________________
846Bool_t AliTOFReconstructionerV2::DigitFinder(TArrayI *secArray, TArrayI *plaArray, TArrayI *strArray, TArrayI *pdzArray, TArrayI *pdxArray, Int_t* assignedVol, Int_t* digitTrackArray, Float_t& tdc)
847{
848 //
849 // Description
850 // input: arrays with sectors, plates, strips, padz, padx
851 // found during fine-propagation of probe tracks
852 //
853 // output kFALSE if signal is not found
854 // kTRUE if signal is found
855 // in this case the assignedVol array contains the digit volume numbers
856 // and digitTrackArray the track numbers (max 3) contributing to the
857 // digit
858
859 Bool_t dummy=kFALSE;
860 Int_t nFilledSlot=secArray->GetSize();
861
862 // start loop
863 Float_t maxWeight=-1.;
864 Int_t indexOfMaxWeight=-1;
865 for (Int_t i = 0; i < nFilledSlot; i++) {
866 Int_t vol[5]; // location for a digit
867 vol[0] = (*secArray)[i];
868 vol[1] = (*plaArray)[i];
869 vol[2] = (*strArray)[i];
870 vol[3] = (*pdxArray)[i];
871 vol[4] = (*pdzArray)[i];
872
873 // check for digit in the current location
874 if (fDigitsMap->TestHit(vol) != kEmpty) {
875
876 AliTOFdigit *dig = static_cast<AliTOFdigit*>(fDigitsMap->GetHit(vol));
877 Float_t adcWeight=dig->GetAdc();
878 if(adcWeight > maxWeight){
879 maxWeight=adcWeight;
880 indexOfMaxWeight=i;
881 tdc=dig->GetTdc();
882 digitTrackArray=dig->GetTracks();
883
884 } // if(adcWeight > maxWeight)
885 } // close if (fDigitsMap->TestHit(vol) != kEmpty)
886
887 } // end loop
888
889 if(indexOfMaxWeight!=-1){
890 assignedVol[0]=(*secArray)[indexOfMaxWeight];
891 assignedVol[1]=(*plaArray)[indexOfMaxWeight];
892 assignedVol[2]=(*strArray)[indexOfMaxWeight];
893 assignedVol[3]=(*pdxArray)[indexOfMaxWeight];
894 assignedVol[4]=(*pdzArray)[indexOfMaxWeight];
895 dummy=kTRUE;
896 }
897
898 return dummy;
899}
900
901//__________________________________________________________________
902Int_t AliTOFReconstructionerV2::SaveTracks(const Char_t *outname, const Int_t split)
903{
904 //
905 // save reconstructed tracks into
906 // outname file
907 //
908 TDirectory *savedir=gDirectory;
909 const Char_t *name="Writing Output";
910 cerr<<'\n'<<name<<"...\n";
911 gBenchmark->Start(name);
912
913 TFile *out=TFile::Open(outname,"RECREATE");
914 if (!out->IsOpen()) {
915 cerr<<"AliTOFReconstructionerV2::SaveTracks(): ";
916 cerr<<"file for TOF tracks is not open !\n";
917 return 2;
918 }
919
920 out->cd();
921 TTree T("T","tree with TOF tracks");
922 T.Branch("tracks",&fTOFTracks,256000,split);
923
924
925 T.Fill();
926 T.Write();
927 savedir->cd();
928 out->Close();
929 return 0;
930 gBenchmark->Stop(name);
931 gBenchmark->Show(name);
932}
933
934//__________________________________________________________________
935void AliTOFReconstructionerV2::Comparison(Int_t* rtIndex)
936{
937 //
938 // perform MC comparison
939 // used also for track length
940 // for the time being
941 // Connect the AliRoot file containing Geometry, Kine, Hits, and Digits
942
943
944 Int_t nEvent = 0;
945 Char_t *alifile = "galice.root";
946
947 TFile *gafl = (TFile*) gROOT->GetListOfFiles()->FindObject(alifile);
948 if (!gafl) {
949 cout << "Open the ALIROOT-file " << alifile << endl;
950 gafl = new TFile(alifile);
951 }
952 else {
953 cout << alifile << " is already open" << endl;
954 }
955
956 // Get AliRun object from file or create it if not on file
957 gAlice = (AliRun*) gafl->Get("gAlice");
958 if (gAlice)
959 cout << "AliRun object found on file" << endl;
960 else
961 gAlice = new AliRun("gAlice","Alice test program");
962
963
964 // Import the Trees for the event nEvent in the file
965 const Int_t nparticles = gAlice->GetEvent(nEvent);
966 if (nparticles <= 0) return;
967
968 // Get pointers to Alice detectors and Hits containers
969 AliDetector *TOF = gAlice->GetDetector("TOF");
970 Int_t ntracks = (Int_t) gAlice->TreeH()->GetEntries();
971
972 // Start loop on tracks in the hits containers
973 for (Int_t track=0; track < ntracks; track++) {
974
975 if(TOF) {
976 for(AliTOFhitT0* tofHit = (AliTOFhitT0*)TOF->FirstHit(track);
977 tofHit;
978 tofHit=(AliTOFhitT0*)TOF->NextHit()) {
979
980 Int_t ipart = tofHit->GetTrack();
981 if(ipart >= 80000) continue;
982 if(rtIndex[ipart] < 0) continue;
983
984 TParticle *part = gAlice->Particle(ipart);
985
986 // get first the pdg code
987 Int_t pdgCode=part->GetPdgCode();
988
989 // then track length
990 Float_t trackLength=tofHit->GetLen(); // [cm]
991
992 // update the tof TClonesArray with TOF tracks
993 AliTOFTrackV2* oTOFtracks=(AliTOFTrackV2*)fTOFTracks->UncheckedAt(rtIndex[ipart]);
994 oTOFtracks->UpdateTrack(pdgCode,trackLength);
995
996 } // loop on hits connected to the current track
997 } // if(TOF)
998 } // end loop on primary tracks
999}