Reconstruction and PID using transition radiation photons: first implementation ...
[u/mrichter/AliRoot.git] / TRD / AliTRDtracker.cxx
CommitLineData
46d29e70 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
0fa7dfa7 16/* $Id$ */
bbf92647 17
029cd327 18///////////////////////////////////////////////////////////////////////////////
19// //
20// The standard TRD tracker //
21// //
22///////////////////////////////////////////////////////////////////////////////
23
a2cb5b3d 24#include <Riostream.h>
46d29e70 25#include <TFile.h>
46d29e70 26#include <TBranch.h>
5443e65e 27#include <TTree.h>
c630aafd 28#include <TObjArray.h>
46d29e70 29
46d29e70 30#include "AliTRDgeometry.h"
5443e65e 31#include "AliTRDparameter.h"
3c625a9b 32#include "AliTRDgeometryHole.h"
46d29e70 33#include "AliTRDcluster.h"
34#include "AliTRDtrack.h"
b7a0917f 35#include "AliBarrelTrack.h"
36#include "AliESD.h"
46d29e70 37
38#include "AliTRDtracker.h"
39
40ClassImp(AliTRDtracker)
41
029cd327 42 const Float_t AliTRDtracker::fgkSeedDepth = 0.5;
43 const Float_t AliTRDtracker::fgkSeedStep = 0.10;
44 const Float_t AliTRDtracker::fgkSeedGap = 0.25;
5443e65e 45
029cd327 46 const Float_t AliTRDtracker::fgkMaxSeedDeltaZ12 = 40.;
47 const Float_t AliTRDtracker::fgkMaxSeedDeltaZ = 25.;
48 const Float_t AliTRDtracker::fgkMaxSeedC = 0.0052;
49 const Float_t AliTRDtracker::fgkMaxSeedTan = 1.2;
50 const Float_t AliTRDtracker::fgkMaxSeedVertexZ = 150.;
a819a5f7 51
029cd327 52 const Double_t AliTRDtracker::fgkSeedErrorSY = 0.2;
53 const Double_t AliTRDtracker::fgkSeedErrorSY3 = 2.5;
54 const Double_t AliTRDtracker::fgkSeedErrorSZ = 0.1;
bbf92647 55
029cd327 56 const Float_t AliTRDtracker::fgkMinClustersInSeed = 0.7;
bbf92647 57
029cd327 58 const Float_t AliTRDtracker::fgkMinClustersInTrack = 0.5;
59 const Float_t AliTRDtracker::fgkMinFractionOfFoundClusters = 0.8;
bbf92647 60
029cd327 61 const Float_t AliTRDtracker::fgkSkipDepth = 0.3;
62 const Float_t AliTRDtracker::fgkLabelFraction = 0.8;
63 const Float_t AliTRDtracker::fgkWideRoad = 20.;
5443e65e 64
029cd327 65 const Double_t AliTRDtracker::fgkMaxChi2 = 12.;
a819a5f7 66
029cd327 67const Int_t AliTRDtracker::fgkFirstPlane = 5;
68const Int_t AliTRDtracker::fgkLastPlane = 17;
9c9d2487 69
bbf92647 70
89f05372 71//____________________________________________________________________
72AliTRDtracker::AliTRDtracker():AliTracker(),
73 fGeom(0),
74 fPar(0),
75 fNclusters(0),
76 fClusters(0),
77 fNseeds(0),
78 fSeeds(0),
79 fNtracks(0),
80 fTracks(0),
81 fSY2corr(0),
82 fSZ2corr(0),
83 fTimeBinsPerPlane(0),
84 fMaxGap(0),
85 fVocal(kFALSE),
86 fAddTRDseeds(kFALSE),
87 fNoTilt(kFALSE)
88{
b7a0917f 89 // Default constructor
90
89f05372 91 for(Int_t i=0;i<kTrackingSectors;i++) fTrSec[i]=0;
92 for(Int_t j=0;j<5;j++)
93 for(Int_t k=0;k<18;k++) fHoles[j][k]=kFALSE;
94}
46d29e70 95//____________________________________________________________________
c630aafd 96AliTRDtracker::AliTRDtracker(const TFile *geomfile):AliTracker()
46d29e70 97{
5443e65e 98 //
99 // Main constructor
100 //
46d29e70 101
c630aafd 102 //Float_t fTzero = 0;
b8dc2353 103
5443e65e 104 fAddTRDseeds = kFALSE;
5443e65e 105 fGeom = NULL;
b8dc2353 106 fNoTilt = kFALSE;
5443e65e 107
108 TDirectory *savedir=gDirectory;
109 TFile *in=(TFile*)geomfile;
110 if (!in->IsOpen()) {
111 printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
112 printf(" DETAIL TRD geometry and DEFAULT TRD parameter will be used\n");
113 }
114 else {
115 in->cd();
c630aafd 116// in->ls();
5443e65e 117 fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
118 fPar = (AliTRDparameter*) in->Get("TRDparameter");
c630aafd 119// fGeom->Dump();
5443e65e 120 }
46d29e70 121
5443e65e 122 if(fGeom) {
123 // fTzero = geo->GetT0();
b8dc2353 124 printf("Found geometry version %d on file \n", fGeom->IsVersion());
5443e65e 125 }
126 else {
c630aafd 127 printf("AliTRDtracker::AliTRDtracker(): can't find TRD geometry!\n");
3c625a9b 128 //printf("The DETAIL TRD geometry will be used\n");
129 //fGeom = new AliTRDgeometryDetail();
130 fGeom = new AliTRDgeometryHole();
131 fGeom->SetPHOShole();
132 fGeom->SetRICHhole();
c630aafd 133 }
134
135 if (!fPar) {
136 printf("AliTRDtracker::AliTRDtracker(): can't find TRD parameter!\n");
137 printf("The DEFAULT TRD parameter will be used\n");
5443e65e 138 fPar = new AliTRDparameter();
139 }
c630aafd 140 fPar->ReInit();
5443e65e 141
142 savedir->cd();
46d29e70 143
5443e65e 144
145 // fGeom->SetT0(fTzero);
0a29d0f1 146
46d29e70 147 fNclusters = 0;
148 fClusters = new TObjArray(2000);
149 fNseeds = 0;
5443e65e 150 fSeeds = new TObjArray(2000);
46d29e70 151 fNtracks = 0;
5443e65e 152 fTracks = new TObjArray(1000);
a819a5f7 153
029cd327 154 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
155 Int_t trS = CookSectorIndex(geomS);
156 fTrSec[trS] = new AliTRDtrackingSector(fGeom, geomS, fPar);
3c625a9b 157 for (Int_t icham=0;icham<AliTRDgeometry::kNcham; icham++){
158 fHoles[icham][trS]=fGeom->IsHole(0,icham,geomS);
159 }
5443e65e 160 }
a819a5f7 161
029cd327 162 Float_t tiltAngle = TMath::Abs(fPar->GetTiltingAngle());
163 if(tiltAngle < 0.1) {
b8dc2353 164 fNoTilt = kTRUE;
165 }
166
167 fSY2corr = 0.2;
168 fSZ2corr = 120.;
169
029cd327 170 if(fNoTilt && (tiltAngle > 0.1)) fSY2corr = fSY2corr + tiltAngle * 0.05;
b8dc2353 171
bbf92647 172
5443e65e 173 // calculate max gap on track
a819a5f7 174
5443e65e 175 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
176 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
a819a5f7 177
5443e65e 178 Double_t dx = (Double_t) fPar->GetTimeBinSize();
179 Int_t tbAmp = fPar->GetTimeBefore();
180 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
b3a5a838 181 if(kTRUE) maxAmp = 0; // intentional until we change the parameter class
5443e65e 182 Int_t tbDrift = fPar->GetTimeMax();
183 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
a819a5f7 184
5443e65e 185 tbDrift = TMath::Min(tbDrift,maxDrift);
186 tbAmp = TMath::Min(tbAmp,maxAmp);
46d29e70 187
5443e65e 188 fTimeBinsPerPlane = tbAmp + tbDrift;
029cd327 189 fMaxGap = (Int_t) (fTimeBinsPerPlane * fGeom->Nplan() * fgkSkipDepth);
46d29e70 190
5443e65e 191 fVocal = kFALSE;
0a29d0f1 192
9c9d2487 193
194 // Barrel Tracks [SR, 03.04.2003]
195
196 fBarrelFile = 0;
197 fBarrelTree = 0;
198 fBarrelArray = 0;
199 fBarrelTrack = 0;
200
201 savedir->cd();
5443e65e 202}
46d29e70 203
5443e65e 204//___________________________________________________________________
205AliTRDtracker::~AliTRDtracker()
46d29e70 206{
029cd327 207 //
208 // Destructor of AliTRDtracker
209 //
210
89f05372 211 if (fClusters) {
212 fClusters->Delete();
213 delete fClusters;
214 }
215 if (fTracks) {
216 fTracks->Delete();
217 delete fTracks;
218 }
219 if (fSeeds) {
220 fSeeds->Delete();
221 delete fSeeds;
222 }
5443e65e 223 delete fGeom;
224 delete fPar;
0a29d0f1 225
029cd327 226 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
227 delete fTrSec[geomS];
5443e65e 228 }
229}
46d29e70 230
9c9d2487 231//_____________________________________________________________________
232
233void AliTRDtracker::SetBarrelTree(const char *mode) {
234 //
235 //
236 //
237
238 if (!IsStoringBarrel()) return;
239
240 TDirectory *sav = gDirectory;
241 if (!fBarrelFile) fBarrelFile = new TFile("AliBarrelTracks.root", "UPDATE");
242
243 char buff[40];
244 sprintf(buff, "BarrelTRD_%d_%s", GetEventNumber(), mode);
245
246 fBarrelFile->cd();
247 fBarrelTree = new TTree(buff, "Barrel TPC tracks");
248
029cd327 249 Int_t nRefs = fgkLastPlane - fgkFirstPlane + 1;
9c9d2487 250
251 if (!fBarrelArray) fBarrelArray = new TClonesArray("AliBarrelTrack", nRefs);
252 for(Int_t i=0; i<nRefs; i++) new((*fBarrelArray)[i]) AliBarrelTrack();
253
254 fBarrelTree->Branch("tracks", &fBarrelArray);
255 sav->cd();
256}
257
258//_____________________________________________________________________
259
260void AliTRDtracker::StoreBarrelTrack(AliTRDtrack *ps, Int_t refPlane, Int_t isIn) {
261 //
262 //
263 //
264
265 if (!IsStoringBarrel()) return;
266
267 static Int_t nClusters;
268 static Int_t nWrong;
269 static Double_t chi2;
270 static Int_t index;
271 static Bool_t wasLast = kTRUE;
272
273 Int_t newClusters, newWrong;
274 Double_t newChi2;
275
276 if (wasLast) {
277
278 fBarrelArray->Clear();
279 nClusters = nWrong = 0;
280 chi2 = 0.0;
281 index = 0;
282 wasLast = kFALSE;
283 }
284
285 fBarrelTrack = (AliBarrelTrack*)(*fBarrelArray)[index++];
286 ps->GetBarrelTrack(fBarrelTrack);
287
288 newClusters = ps->GetNumberOfClusters() - nClusters;
289 newWrong = ps->GetNWrong() - nWrong;
290 newChi2 = ps->GetChi2() - chi2;
291
292 nClusters = ps->GetNumberOfClusters();
293 nWrong = ps->GetNWrong();
294 chi2 = ps->GetChi2();
295
029cd327 296 if (refPlane != fgkLastPlane) {
9c9d2487 297 fBarrelTrack->SetNClusters(newClusters, newChi2);
298 fBarrelTrack->SetNWrongClusters(newWrong);
299 } else {
300 wasLast = kTRUE;
301 }
302
303 fBarrelTrack->SetRefPlane(refPlane, isIn);
304}
305
306//_____________________________________________________________________
307
308Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track) {
309 //
310 // Rotates the track when necessary
311 //
312
313 Double_t alpha = AliTRDgeometry::GetAlpha();
314 Double_t y = track->GetY();
315 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
316
c630aafd 317 //Int_t ns = AliTRDgeometry::kNsect;
9c9d2487 318 //Int_t s=Int_t(track->GetAlpha()/alpha)%ns;
319
320 if (y > ymax) {
321 //s = (s+1) % ns;
322 if (!track->Rotate(alpha)) return kFALSE;
323 } else if (y <-ymax) {
324 //s = (s-1+ns) % ns;
325 if (!track->Rotate(-alpha)) return kFALSE;
326 }
327
328 return kTRUE;
329}
330
46d29e70 331//_____________________________________________________________________
5443e65e 332inline Double_t f1trd(Double_t x1,Double_t y1,
a9814c09 333 Double_t x2,Double_t y2,
334 Double_t x3,Double_t y3)
46d29e70 335{
0a29d0f1 336 //
46d29e70 337 // Initial approximation of the track curvature
0a29d0f1 338 //
46d29e70 339 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
340 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
341 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
342 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
343 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
344
345 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
346
347 return -xr*yr/sqrt(xr*xr+yr*yr);
348}
349
350//_____________________________________________________________________
5443e65e 351inline Double_t f2trd(Double_t x1,Double_t y1,
a9814c09 352 Double_t x2,Double_t y2,
353 Double_t x3,Double_t y3)
46d29e70 354{
0a29d0f1 355 //
5443e65e 356 // Initial approximation of the track curvature times X coordinate
357 // of the center of curvature
0a29d0f1 358 //
46d29e70 359
360 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
361 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
362 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
363 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
364 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
365
366 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
367
368 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
369}
370
371//_____________________________________________________________________
5443e65e 372inline Double_t f3trd(Double_t x1,Double_t y1,
a9814c09 373 Double_t x2,Double_t y2,
374 Double_t z1,Double_t z2)
46d29e70 375{
0a29d0f1 376 //
46d29e70 377 // Initial approximation of the tangent of the track dip angle
0a29d0f1 378 //
46d29e70 379
380 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
381}
382
46e2d86c 383
384AliTRDcluster * AliTRDtracker::GetCluster(AliTRDtrack * track, Int_t plane, Int_t timebin){
385 //
386 //try to find cluster in the backup list
387 //
388 AliTRDcluster * cl =0;
389 UInt_t *indexes = track->GetBackupIndexes();
390 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
391 if (indexes[i]==0) break;
392 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
393 if (!cli) break;
394 if (cli->GetLocalTimeBin()!=timebin) continue;
395 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
396 if (iplane==plane) {
397 cl = cli;
398 break;
399 }
400 }
401 return cl;
402}
403
3c625a9b 404
405Int_t AliTRDtracker::GetLastPlane(AliTRDtrack * track){
406 //
407 //return last updated plane
408 Int_t lastplane=0;
409 UInt_t *indexes = track->GetBackupIndexes();
410 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
411 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
412 if (!cli) break;
413 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
414 if (iplane>lastplane) {
415 lastplane = iplane;
416 }
417 }
418 return lastplane;
419}
c630aafd 420//___________________________________________________________________
421Int_t AliTRDtracker::Clusters2Tracks(AliESD* event)
422{
423 //
424 // Finds tracks within the TRD. The ESD event is expected to contain seeds
425 // at the outer part of the TRD. The seeds
426 // are found within the TRD if fAddTRDseeds is TRUE.
427 // The tracks are propagated to the innermost time bin
428 // of the TRD and the ESD event is updated
429 //
430
431 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
029cd327 432 Float_t foundMin = fgkMinClustersInTrack * timeBins;
c630aafd 433 Int_t nseed = 0;
434 Int_t found = 0;
435 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
436
437 Int_t n = event->GetNumberOfTracks();
438 for (Int_t i=0; i<n; i++) {
439 AliESDtrack* seed=event->GetTrack(i);
440 ULong_t status=seed->GetStatus();
441 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
442 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
443 nseed++;
444
445 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
46e2d86c 446 //seed2->ResetCovariance();
c630aafd 447 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
448 AliTRDtrack &t=*pt;
449 FollowProlongation(t, innerTB);
450 if (t.GetNumberOfClusters() >= foundMin) {
451 UseClusters(&t);
029cd327 452 CookLabel(pt, 1-fgkLabelFraction);
c630aafd 453 // t.CookdEdx();
454 }
455 found++;
456// cout<<found<<'\r';
457
458 if(PropagateToTPC(t)) {
459 seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
460 }
461 delete seed2;
462 delete pt;
463 }
464
465 cout<<"Number of loaded seeds: "<<nseed<<endl;
466 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
467
468 // after tracks from loaded seeds are found and the corresponding
469 // clusters are used, look for additional seeds from TRD
470
471 if(fAddTRDseeds) {
472 // Find tracks for the seeds in the TRD
473 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
474
029cd327 475 Int_t nSteps = (Int_t) (fgkSeedDepth / fgkSeedStep);
476 Int_t gap = (Int_t) (timeBins * fgkSeedGap);
477 Int_t step = (Int_t) (timeBins * fgkSeedStep);
c630aafd 478
479 // make a first turn with tight cut on initial curvature
480 for(Int_t turn = 1; turn <= 2; turn++) {
481 if(turn == 2) {
029cd327 482 nSteps = (Int_t) (fgkSeedDepth / (3*fgkSeedStep));
483 step = (Int_t) (timeBins * (3*fgkSeedStep));
c630aafd 484 }
485 for(Int_t i=0; i<nSteps; i++) {
486 Int_t outer=timeBins-1-i*step;
487 Int_t inner=outer-gap;
488
489 nseed=fSeeds->GetEntriesFast();
490
491 MakeSeeds(inner, outer, turn);
492
493 nseed=fSeeds->GetEntriesFast();
7bed16a7 494 // printf("\n turn %d, step %d: number of seeds for TRD inward %d\n",
495 // turn, i, nseed);
c630aafd 496
497 for (Int_t i=0; i<nseed; i++) {
498 AliTRDtrack *pt=(AliTRDtrack*)fSeeds->UncheckedAt(i), &t=*pt;
499 FollowProlongation(t,innerTB);
500 if (t.GetNumberOfClusters() >= foundMin) {
501 UseClusters(&t);
029cd327 502 CookLabel(pt, 1-fgkLabelFraction);
c630aafd 503 t.CookdEdx();
504 found++;
505// cout<<found<<'\r';
506 if(PropagateToTPC(t)) {
507 AliESDtrack track;
508 track.UpdateTrackParams(pt,AliESDtrack::kTRDin);
509 event->AddTrack(&track);
c5a8e3df 510 // track.SetTRDtrack(new AliTRDtrack(*pt));
c630aafd 511 }
512 }
513 delete fSeeds->RemoveAt(i);
514 fNseeds--;
515 }
516 }
517 }
518 }
519
520 cout<<"Total number of found tracks: "<<found<<endl;
521
522 return 0;
523}
5443e65e 524
c630aafd 525
5443e65e 526
c630aafd 527//_____________________________________________________________________________
528Int_t AliTRDtracker::PropagateBack(AliESD* event) {
529 //
530 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
531 // backpropagated by the TPC tracker. Each seed is first propagated
532 // to the TRD, and then its prolongation is searched in the TRD.
533 // If sufficiently long continuation of the track is found in the TRD
534 // the track is updated, otherwise it's stored as originaly defined
535 // by the TPC tracker.
536 //
537
538 Int_t found=0;
c5a8e3df 539 Float_t foundMin = 20;
c630aafd 540
541 Int_t n = event->GetNumberOfTracks();
542 for (Int_t i=0; i<n; i++) {
543 AliESDtrack* seed=event->GetTrack(i);
544 ULong_t status=seed->GetStatus();
545 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
546 if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
547
548 Int_t lbl = seed->GetLabel();
549 AliTRDtrack *track = new AliTRDtrack(*seed);
550 track->SetSeedLabel(lbl);
f4e9508c 551 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup); //make backup
c630aafd 552 fNseeds++;
7bed16a7 553 Float_t p4 = track->GetC();
554 //
3c625a9b 555 Int_t expectedClr = FollowBackProlongation(*track);
f4e9508c 556 /*
557 // only debug purpose
3c625a9b 558 if (track->GetNumberOfClusters()<expectedClr/3){
559 AliTRDtrack *track1 = new AliTRDtrack(*seed);
560 track1->SetSeedLabel(lbl);
561 FollowBackProlongation(*track1);
562 AliTRDtrack *track2= new AliTRDtrack(*seed);
563 track->SetSeedLabel(lbl);
564 FollowBackProlongation(*track2);
565 delete track1;
566 delete track2;
567 }
f4e9508c 568 */
569 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)>0.2) {
7bed16a7 570 delete track;
571 continue; //too big change of curvature - to be checked
572 }
f4e9508c 573
c630aafd 574 Int_t foundClr = track->GetNumberOfClusters();
575 if (foundClr >= foundMin) {
eab5961e 576 track->CookdEdx(0.,1.);
577 CookdEdxTimBin(*track);
578
f4e9508c 579 CookLabel(track, 1-fgkLabelFraction);
580 if(track->GetChi2()/track->GetNumberOfClusters()<6) { // sign only gold tracks
c630aafd 581 UseClusters(track);
582 }
f4e9508c 583 Bool_t isGold = kFALSE;
584
585 if (track->GetChi2()/track->GetNumberOfClusters()<5) { //full gold track
586 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
587 isGold = kTRUE;
588 }
589 if (!isGold && track->GetNCross()==0&&track->GetChi2()/track->GetNumberOfClusters()<7){ //almost gold track
590 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
591 isGold = kTRUE;
592 }
593 if (!isGold && track->GetBackupTrack()){
594 if (track->GetBackupTrack()->GetNumberOfClusters()>foundMin&&
595 (track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1))<7){
596 seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
597 isGold = kTRUE;
598 }
16d9fbba 599 }
c630aafd 600 }
f4e9508c 601 else{
602 continue;
603 }
c630aafd 604
3c625a9b 605
606 if (track->GetStop()==kFALSE){
c5a8e3df 607
b94f0a96 608 Double_t xtof=371.;
3c625a9b 609 Double_t c2=track->GetC()*xtof - track->GetEta();
c5a8e3df 610 if (TMath::Abs(c2)>=0.85) {
611 delete track;
612 continue;
613 }
b94f0a96 614 Double_t xTOF0 = 371. ;
16d9fbba 615 PropagateToOuterPlane(*track,xTOF0);
7bed16a7 616 //
3c625a9b 617 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
618 Double_t y=track->GetYat(xtof);
619 if (y > ymax) {
7ac6fa52 620 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
621 delete track;
7bed16a7 622 continue;
7ac6fa52 623 }
3c625a9b 624 } else if (y <-ymax) {
7ac6fa52 625 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
626 delete track;
7bed16a7 627 continue;
7ac6fa52 628 }
3c625a9b 629 }
630
631 if (track->PropagateTo(xtof)) {
eab5961e 632 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
633 for (Int_t i=0;i<kNPlane;i++) {
634 seed->SetTRDsignals(track->GetPIDsignals(i),i);
635 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
636 }
c5a8e3df 637 seed->SetTRDtrack(new AliTRDtrack(*track));
3c625a9b 638 if (track->GetNumberOfClusters()>foundMin) found++;
639 }
640 }else{
641 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
642 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
16d9fbba 643 //seed->SetStatus(AliESDtrack::kTRDStop);
eab5961e 644 for (Int_t i=0;i<kNPlane;i++) {
645 seed->SetTRDsignals(track->GetPIDsignals(i),i);
646 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
647 }
c5a8e3df 648 seed->SetTRDtrack(new AliTRDtrack(*track));
3c625a9b 649 found++;
650 }
1e9bb598 651 }
652
d9b8978b 653 delete track;
3c625a9b 654
1e9bb598 655 //End of propagation to the TOF
3c625a9b 656 //if (foundClr>foundMin)
657 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
658
c630aafd 659
660 }
661
662 cerr<<"Number of seeds: "<<fNseeds<<endl;
663 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
664
1e9bb598 665 fSeeds->Clear(); fNseeds=0;
666
667 return 0;
668
669}
670
671//_____________________________________________________________________________
672Int_t AliTRDtracker::RefitInward(AliESD* event)
673{
674 //
675 // Refits tracks within the TRD. The ESD event is expected to contain seeds
676 // at the outer part of the TRD.
677 // The tracks are propagated to the innermost time bin
678 // of the TRD and the ESD event is updated
679 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
680 //
681
682 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
683 Float_t foundMin = fgkMinClustersInTrack * timeBins;
684 Int_t nseed = 0;
685 Int_t found = 0;
686 Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
687
688 Int_t n = event->GetNumberOfTracks();
689 for (Int_t i=0; i<n; i++) {
690 AliESDtrack* seed=event->GetTrack(i);
f4e9508c 691 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
692 if (seed2->GetX()<270){
693 seed->UpdateTrackParams(seed2, AliESDtrack::kTRDbackup); // backup TPC track - only update
694 continue;
695 }
696
1e9bb598 697 ULong_t status=seed->GetStatus();
698 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
46e2d86c 699 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
f4e9508c 700 nseed++;
3c625a9b 701 seed2->ResetCovariance(5.);
1e9bb598 702 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
eab5961e 703 for (Int_t i=0;i<kNPlane;i++) {
704 pt->SetPIDsignals(seed2->GetPIDsignals(i),i);
705 pt->SetPIDTimBin(seed2->GetPIDTimBin(i),i);
706 }
707
46e2d86c 708 UInt_t * indexes2 = seed2->GetIndexes();
709 UInt_t * indexes3 = pt->GetBackupIndexes();
710 for (Int_t i=0;i<200;i++) {
711 if (indexes2[i]==0) break;
712 indexes3[i] = indexes2[i];
713 }
714 //AliTRDtrack *pt = seed2;
1e9bb598 715 AliTRDtrack &t=*pt;
716 FollowProlongation(t, innerTB);
7bed16a7 717 /*
46e2d86c 718 if (t.GetNumberOfClusters()<seed->GetTRDclusters(indexes3)*0.5){
719 // debug - why we dont go back?
720 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
721 UInt_t * indexes2 = seed2->GetIndexes();
722 UInt_t * indexes3 = pt2->GetBackupIndexes();
723 for (Int_t i=0;i<200;i++) {
724 if (indexes2[i]==0) break;
725 indexes3[i] = indexes2[i];
726 }
3c625a9b 727 FollowProlongation(*pt2, innerTB);
728 delete pt2;
46e2d86c 729 }
7bed16a7 730 */
1e9bb598 731 if (t.GetNumberOfClusters() >= foundMin) {
46e2d86c 732 // UseClusters(&t);
733 //CookLabel(pt, 1-fgkLabelFraction);
1e9bb598 734 // t.CookdEdx();
735 }
736 found++;
737// cout<<found<<'\r';
738
739 if(PropagateToTPC(t)) {
0fa7dfa7 740 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
eab5961e 741 for (Int_t i=0;i<kNPlane;i++) {
742 seed->SetTRDsignals(pt->GetPIDsignals(i),i);
743 seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
744 }
7bed16a7 745 }else{
746 //if not prolongation to TPC - propagate without update
747 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
748 seed2->ResetCovariance(5.);
749 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
750 delete seed2;
eab5961e 751 if (PropagateToTPC(*pt2)) {
752 pt2->CookdEdx(0.,1.);
753 CookdEdxTimBin(*pt2);
7bed16a7 754 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
eab5961e 755 for (Int_t i=0;i<kNPlane;i++) {
756 seed->SetTRDsignals(seed2->GetPIDsignals(i),i);
757 seed->SetTRDTimBin(seed2->GetPIDTimBin(i),i);
758 }
759 }
7bed16a7 760 delete pt2;
1e9bb598 761 }
7bed16a7 762
1e9bb598 763 delete seed2;
764 delete pt;
eab5961e 765 }
1e9bb598 766
767 cout<<"Number of loaded seeds: "<<nseed<<endl;
768 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
769
c630aafd 770 return 0;
771
772}
773
bbf92647 774
5443e65e 775//---------------------------------------------------------------------------
776Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t, Int_t rf)
777{
778 // Starting from current position on track=t this function tries
779 // to extrapolate the track up to timeBin=0 and to confirm prolongation
780 // if a close cluster is found. Returns the number of clusters
781 // expected to be found in sensitive layers
bbf92647 782
5443e65e 783 Float_t wIndex, wTB, wChi2;
784 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
785 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
786 Float_t wPx, wPy, wPz, wC;
029cd327 787 Double_t px, py, pz;
5443e65e 788 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
3c625a9b 789 Int_t lastplane = GetLastPlane(&t);
46d29e70 790
5443e65e 791 Int_t trackIndex = t.GetLabel();
46d29e70 792
5443e65e 793 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
46d29e70 794
029cd327 795 Int_t tryAgain=fMaxGap;
46d29e70 796
5443e65e 797 Double_t alpha=t.GetAlpha();
c630aafd 798 alpha = TVector2::Phi_0_2pi(alpha);
46d29e70 799
5443e65e 800 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
029cd327 801 Double_t radLength, rho, x, dx, y, ymax, z;
46d29e70 802
5443e65e 803 Int_t expectedNumberOfClusters = 0;
804 Bool_t lookForCluster;
46d29e70 805
5443e65e 806 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
46d29e70 807
5443e65e 808
809 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
810
811 y = t.GetY(); z = t.GetZ();
812
813 // first propagate to the inner surface of the current time bin
029cd327 814 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
5443e65e 815 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
029cd327 816 if(!t.PropagateTo(x,radLength,rho)) break;
5443e65e 817 y = t.GetY();
818 ymax = x*TMath::Tan(0.5*alpha);
819 if (y > ymax) {
820 s = (s+1) % ns;
821 if (!t.Rotate(alpha)) break;
029cd327 822 if(!t.PropagateTo(x,radLength,rho)) break;
5443e65e 823 } else if (y <-ymax) {
824 s = (s-1+ns) % ns;
825 if (!t.Rotate(-alpha)) break;
029cd327 826 if(!t.PropagateTo(x,radLength,rho)) break;
5443e65e 827 }
46d29e70 828
5443e65e 829 y = t.GetY(); z = t.GetZ();
830
831 // now propagate to the middle plane of the next time bin
029cd327 832 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
5443e65e 833 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
029cd327 834 if(!t.PropagateTo(x,radLength,rho)) break;
5443e65e 835 y = t.GetY();
836 ymax = x*TMath::Tan(0.5*alpha);
837 if (y > ymax) {
838 s = (s+1) % ns;
839 if (!t.Rotate(alpha)) break;
029cd327 840 if(!t.PropagateTo(x,radLength,rho)) break;
5443e65e 841 } else if (y <-ymax) {
842 s = (s-1+ns) % ns;
843 if (!t.Rotate(-alpha)) break;
029cd327 844 if(!t.PropagateTo(x,radLength,rho)) break;
5443e65e 845 }
46d29e70 846
46d29e70 847
5443e65e 848 if(lookForCluster) {
a819a5f7 849
5443e65e 850 expectedNumberOfClusters++;
5443e65e 851 wIndex = (Float_t) t.GetLabel();
852 wTB = nr;
46d29e70 853
029cd327 854 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr-1));
46d29e70 855
5443e65e 856 Double_t sy2=ExpectedSigmaY2(x,t.GetTgl(),t.GetPt());
5443e65e 857 Double_t sz2=ExpectedSigmaZ2(x,t.GetTgl());
46d29e70 858
5443e65e 859 Double_t road;
860 if((t.GetSigmaY2() + sy2) > 0) road=10.*sqrt(t.GetSigmaY2() + sy2);
861 else return expectedNumberOfClusters;
862
863 wYrt = (Float_t) y;
864 wZrt = (Float_t) z;
865 wYwindow = (Float_t) road;
029cd327 866 t.GetPxPyPz(px,py,pz);
867 wPx = (Float_t) px;
868 wPy = (Float_t) py;
869 wPz = (Float_t) pz;
5443e65e 870 wC = (Float_t) t.GetC();
871 wSigmaC2 = (Float_t) t.GetSigmaC2();
872 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
873 wSigmaY2 = (Float_t) t.GetSigmaY2();
874 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
875 wChi2 = -1;
876
5443e65e 877
878 AliTRDcluster *cl=0;
879 UInt_t index=0;
880
029cd327 881 Double_t maxChi2=fgkMaxChi2;
5443e65e 882
883 wYclosest = 12345678;
884 wYcorrect = 12345678;
885 wZclosest = 12345678;
886 wZcorrect = 12345678;
887 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
888
889 // Find the closest correct cluster for debugging purposes
029cd327 890 if (timeBin) {
a9814c09 891 Float_t minDY = 1000000;
029cd327 892 for (Int_t i=0; i<timeBin; i++) {
893 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
a9814c09 894 if((c->GetLabel(0) != trackIndex) &&
895 (c->GetLabel(1) != trackIndex) &&
896 (c->GetLabel(2) != trackIndex)) continue;
897 if(TMath::Abs(c->GetY() - y) > minDY) continue;
898 minDY = TMath::Abs(c->GetY() - y);
899 wYcorrect = c->GetY();
900 wZcorrect = c->GetZ();
901
902 Double_t h01 = GetTiltFactor(c);
903 wChi2 = t.GetPredictedChi2(c, h01);
904 }
5443e65e 905 }
46d29e70 906
5443e65e 907 // Now go for the real cluster search
a819a5f7 908
029cd327 909 if (timeBin) {
46e2d86c 910 //
911 //find cluster in history
912 cl =0;
913
914 AliTRDcluster * cl0 = timeBin[0];
915 if (!cl0) {
916 continue;
917 }
918 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
3c625a9b 919 if (plane>lastplane) continue;
46e2d86c 920 Int_t timebin = cl0->GetLocalTimeBin();
921 AliTRDcluster * cl2= GetCluster(&t,plane, timebin);
3c625a9b 922 if (cl2) {
923 cl =cl2;
924 Double_t h01 = GetTiltFactor(cl);
925 maxChi2=t.GetPredictedChi2(cl,h01);
926 }
46e2d86c 927 if ((!cl) && road>fgkWideRoad) {
928 //if (t.GetNumberOfClusters()>4)
929 // cerr<<t.GetNumberOfClusters()
930 // <<"FindProlongation warning: Too broad road !\n";
931 continue;
932 }
933
934
935 if(!cl){
936
937 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
938 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
939 if (c->GetY() > y+road) break;
940 if (c->IsUsed() > 0) continue;
941 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
942
943 Double_t h01 = GetTiltFactor(c);
944 Double_t chi2=t.GetPredictedChi2(c,h01);
945
946 if (chi2 > maxChi2) continue;
947 maxChi2=chi2;
948 cl=c;
949 index=timeBin.GetIndex(i);
950 }
951 }
a9814c09 952
953 if(!cl) {
954
029cd327 955 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
956 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
a9814c09 957
958 if (c->GetY() > y+road) break;
959 if (c->IsUsed() > 0) continue;
960 if((c->GetZ()-z)*(c->GetZ()-z) > 12 * sz2) continue;
961
962 Double_t h01 = GetTiltFactor(c);
963 Double_t chi2=t.GetPredictedChi2(c, h01);
964
029cd327 965 if (chi2 > maxChi2) continue;
966 maxChi2=chi2;
a9814c09 967 cl=c;
029cd327 968 index=timeBin.GetIndex(i);
a9814c09 969 }
970 }
a9814c09 971 if (cl) {
972 wYclosest = cl->GetY();
973 wZclosest = cl->GetZ();
974 Double_t h01 = GetTiltFactor(cl);
975
976 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
46e2d86c 977 //printf("Track position\t%f\t%f\t%f\n",t.GetX(),t.GetY(),t.GetZ());
978 //printf("Cluster position\t%d\t%f\t%f\n",cl->GetLocalTimeBin(),cl->GetY(),cl->GetZ());
3c625a9b 979 Int_t det = cl->GetDetector();
980 Int_t plane = fGeom->GetPlane(det);
46e2d86c 981
3c625a9b 982 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
983 //if(!t.Update(cl,maxChi2,index,h01)) {
984 //if(!tryAgain--) return 0;
a9814c09 985 }
029cd327 986 else tryAgain=fMaxGap;
a9814c09 987 }
988 else {
3c625a9b 989 //if (tryAgain==0) break;
029cd327 990 tryAgain--;
a9814c09 991 }
992
993 /*
994 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
995
996 printf(" %f", wIndex); //1
997 printf(" %f", wTB); //2
998 printf(" %f", wYrt); //3
999 printf(" %f", wYclosest); //4
1000 printf(" %f", wYcorrect); //5
1001 printf(" %f", wYwindow); //6
1002 printf(" %f", wZrt); //7
1003 printf(" %f", wZclosest); //8
1004 printf(" %f", wZcorrect); //9
1005 printf(" %f", wZwindow); //10
1006 printf(" %f", wPx); //11
1007 printf(" %f", wPy); //12
1008 printf(" %f", wPz); //13
1009 printf(" %f", wSigmaC2*1000000); //14
1010 printf(" %f", wSigmaTgl2*1000); //15
1011 printf(" %f", wSigmaY2); //16
1012 // printf(" %f", wSigmaZ2); //17
1013 printf(" %f", wChi2); //17
1014 printf(" %f", wC); //18
1015 printf("\n");
1016 }
1017 */
5443e65e 1018 }
1019 }
1020 }
1021 return expectedNumberOfClusters;
1022
1023
1024}
a819a5f7 1025
5443e65e 1026//___________________________________________________________________
46d29e70 1027
5443e65e 1028Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
1029{
1030 // Starting from current radial position of track <t> this function
1031 // extrapolates the track up to outer timebin and in the sensitive
1032 // layers confirms prolongation if a close cluster is found.
1033 // Returns the number of clusters expected to be found in sensitive layers
46d29e70 1034
16d9fbba 1035
5443e65e 1036 Float_t wIndex, wTB, wChi2;
1037 Float_t wYrt, wYclosest, wYcorrect, wYwindow;
1038 Float_t wZrt, wZclosest, wZcorrect, wZwindow;
1039 Float_t wPx, wPy, wPz, wC;
029cd327 1040 Double_t px, py, pz;
5443e65e 1041 Float_t wSigmaC2, wSigmaTgl2, wSigmaY2, wSigmaZ2;
46d29e70 1042
5443e65e 1043 Int_t trackIndex = t.GetLabel();
029cd327 1044 Int_t tryAgain=fMaxGap;
46d29e70 1045
5443e65e 1046 Double_t alpha=t.GetAlpha();
9c9d2487 1047 TVector2::Phi_0_2pi(alpha);
46d29e70 1048
9c9d2487 1049 Int_t s;
46d29e70 1050
5443e65e 1051 Int_t outerTB = fTrSec[0]->GetOuterTimeBin();
029cd327 1052 Double_t radLength, rho, x, dx, y, ymax = 0, z;
5443e65e 1053 Bool_t lookForCluster;
a819a5f7 1054
5443e65e 1055 Int_t expectedNumberOfClusters = 0;
1056 x = t.GetX();
46d29e70 1057
5443e65e 1058 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
46d29e70 1059
029cd327 1060 Int_t nRefPlane = fgkFirstPlane;
9c9d2487 1061 Bool_t isNewLayer = kFALSE;
46d29e70 1062
9c9d2487 1063 Double_t chi2;
1064 Double_t minDY;
3c625a9b 1065 Int_t zone =-10;
1066 Int_t nr;
1067 for (nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr<outerTB+1; nr++) {
9c9d2487 1068
1069 y = t.GetY();
1070 z = t.GetZ();
46d29e70 1071
5443e65e 1072 // first propagate to the outer surface of the current time bin
46d29e70 1073
9c9d2487 1074 s = t.GetSector();
029cd327 1075 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
9c9d2487 1076 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2;
1077 y = t.GetY();
1078 z = t.GetZ();
46d29e70 1079
029cd327 1080 if(!t.PropagateTo(x,radLength,rho)) break;
3c625a9b 1081 // if (!AdjustSector(&t)) break;
1082 //
1083 // MI -fix untill correct material desription will be implemented
1084 //
16d9fbba 1085 Float_t angle = t.GetAlpha(); // MI - if rotation - we go through the material
9c9d2487 1086 if (!AdjustSector(&t)) break;
16d9fbba 1087 Int_t cross = kFALSE;
1088
1089 if (TMath::Abs(angle - t.GetAlpha())>0.000001) cross = kTRUE; //better to stop track
3c625a9b 1090 Int_t currentzone = fTrSec[s]->GetLayer(nr)->GetZone(z);
16d9fbba 1091 if (currentzone==-10) cross = kTRUE; // we are in the frame
3c625a9b 1092 if (currentzone>-10){ // layer knows where we are
1093 if (zone==-10) zone = currentzone;
16d9fbba 1094 if (zone!=currentzone) cross=kTRUE;
1095 }
1096 if (cross) {
1097 t.IncCross();
1098 if (t.GetNCross()==1) t.MakeBackupTrack();
1099 if (t.GetNCross()>2) break;
3c625a9b 1100 }
16d9fbba 1101
3c625a9b 1102 //
1103 //
9c9d2487 1104 s = t.GetSector();
029cd327 1105 if (!t.PropagateTo(x,radLength,rho)) break;
46d29e70 1106
5443e65e 1107 y = t.GetY();
9c9d2487 1108 z = t.GetZ();
a819a5f7 1109
9c9d2487 1110 // Barrel Tracks [SR, 04.04.2003]
46d29e70 1111
9c9d2487 1112 s = t.GetSector();
1113 if (fTrSec[s]->GetLayer(nr)->IsSensitive() !=
1114 fTrSec[s]->GetLayer(nr+1)->IsSensitive() ) {
46d29e70 1115
c630aafd 1116// if (IsStoringBarrel()) StoreBarrelTrack(&t, nRefPlane++, kTrackBack);
9c9d2487 1117 }
46d29e70 1118
9c9d2487 1119 if (fTrSec[s]->GetLayer(nr-1)->IsSensitive() &&
1120 ! fTrSec[s]->GetLayer(nr)->IsSensitive()) {
1121 isNewLayer = kTRUE;
1122 } else {isNewLayer = kFALSE;}
46d29e70 1123
5443e65e 1124 y = t.GetY();
9c9d2487 1125 z = t.GetZ();
a819a5f7 1126
9c9d2487 1127 // now propagate to the middle plane of the next time bin
029cd327 1128 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
5443e65e 1129
9c9d2487 1130 x = fTrSec[s]->GetLayer(nr+1)->GetX();
029cd327 1131 if(!t.PropagateTo(x,radLength,rho)) break;
9c9d2487 1132 if (!AdjustSector(&t)) break;
1133 s = t.GetSector();
029cd327 1134 if(!t.PropagateTo(x,radLength,rho)) break;
46d29e70 1135
9c9d2487 1136 y = t.GetY();
1137 z = t.GetZ();
1138
1139 if(fVocal) printf("nr+1=%d, x %f, z %f, y %f, ymax %f\n",nr+1,x,z,y,ymax);
5443e65e 1140 // printf("label %d, pl %d, lookForCluster %d \n",
a9814c09 1141 // trackIndex, nr+1, lookForCluster);
5443e65e 1142
1143 if(lookForCluster) {
1144 expectedNumberOfClusters++;
1145
1146 wIndex = (Float_t) t.GetLabel();
1147 wTB = fTrSec[s]->GetLayer(nr+1)->GetTimeBinIndex();
1148
029cd327 1149 AliTRDpropagationLayer& timeBin=*(fTrSec[s]->GetLayer(nr+1));
5443e65e 1150 Double_t sy2=ExpectedSigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
1151 Double_t sz2=ExpectedSigmaZ2(t.GetX(),t.GetTgl());
1152 if((t.GetSigmaY2() + sy2) < 0) break;
1153 Double_t road = 10.*sqrt(t.GetSigmaY2() + sy2);
1154 Double_t y=t.GetY(), z=t.GetZ();
1155
1156 wYrt = (Float_t) y;
1157 wZrt = (Float_t) z;
1158 wYwindow = (Float_t) road;
029cd327 1159 t.GetPxPyPz(px,py,pz);
1160 wPx = (Float_t) px;
1161 wPy = (Float_t) py;
1162 wPz = (Float_t) pz;
5443e65e 1163 wC = (Float_t) t.GetC();
1164 wSigmaC2 = (Float_t) t.GetSigmaC2();
1165 wSigmaTgl2 = (Float_t) t.GetSigmaTgl2();
1166 wSigmaY2 = (Float_t) t.GetSigmaY2();
1167 wSigmaZ2 = (Float_t) t.GetSigmaZ2();
1168 wChi2 = -1;
1169
029cd327 1170 if (road>fgkWideRoad) {
a9814c09 1171 if (t.GetNumberOfClusters()>4)
1172 cerr<<t.GetNumberOfClusters()
1173 <<"FindProlongation warning: Too broad road !\n";
1174 return 0;
5443e65e 1175 }
1176
1177 AliTRDcluster *cl=0;
1178 UInt_t index=0;
1179
029cd327 1180 Double_t maxChi2=fgkMaxChi2;
5443e65e 1181
9c9d2487 1182 if (isNewLayer) {
1183 road = 3 * road;
1184 //sz2 = 3 * sz2;
029cd327 1185 maxChi2 = 10 * fgkMaxChi2;
9c9d2487 1186 }
1187
029cd327 1188 if (nRefPlane == fgkFirstPlane) maxChi2 = 20 * fgkMaxChi2;
1189 if (nRefPlane == fgkFirstPlane+2) maxChi2 = 15 * fgkMaxChi2;
1190 if (t.GetNRotate() > 0) maxChi2 = 3 * maxChi2;
9c9d2487 1191
1192
5443e65e 1193 wYclosest = 12345678;
1194 wYcorrect = 12345678;
1195 wZclosest = 12345678;
1196 wZcorrect = 12345678;
1197 wZwindow = TMath::Sqrt(2.25 * 12 * sz2);
1198
1199 // Find the closest correct cluster for debugging purposes
029cd327 1200 if (timeBin) {
9c9d2487 1201 minDY = 1000000;
029cd327 1202 for (Int_t i=0; i<timeBin; i++) {
1203 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
a9814c09 1204 if((c->GetLabel(0) != trackIndex) &&
1205 (c->GetLabel(1) != trackIndex) &&
1206 (c->GetLabel(2) != trackIndex)) continue;
1207 if(TMath::Abs(c->GetY() - y) > minDY) continue;
9c9d2487 1208 //minDY = TMath::Abs(c->GetY() - y);
1209 minDY = c->GetY() - y;
a9814c09 1210 wYcorrect = c->GetY();
1211 wZcorrect = c->GetZ();
1212
1213 Double_t h01 = GetTiltFactor(c);
1214 wChi2 = t.GetPredictedChi2(c, h01);
1215 }
5443e65e 1216 }
46d29e70 1217
5443e65e 1218 // Now go for the real cluster search
46d29e70 1219
029cd327 1220 if (timeBin) {
5443e65e 1221
029cd327 1222 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1223 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
a9814c09 1224 if (c->GetY() > y+road) break;
9c9d2487 1225 if (c->IsUsed() > 0) continue;
a9814c09 1226 if((c->GetZ()-z)*(c->GetZ()-z) > 3 * sz2) continue;
1227
1228 Double_t h01 = GetTiltFactor(c);
9c9d2487 1229 chi2=t.GetPredictedChi2(c,h01);
a9814c09 1230
029cd327 1231 if (chi2 > maxChi2) continue;
1232 maxChi2=chi2;
a9814c09 1233 cl=c;
029cd327 1234 index=timeBin.GetIndex(i);
9c9d2487 1235
1236 //check is correct
1237 if((c->GetLabel(0) != trackIndex) &&
1238 (c->GetLabel(1) != trackIndex) &&
1239 (c->GetLabel(2) != trackIndex)) t.AddNWrong();
a9814c09 1240 }
5443e65e 1241
a9814c09 1242 if(!cl) {
1243
029cd327 1244 for (Int_t i=timeBin.Find(y-road); i<timeBin; i++) {
1245 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
a9814c09 1246
1247 if (c->GetY() > y+road) break;
9c9d2487 1248 if (c->IsUsed() > 0) continue;
a9814c09 1249 if((c->GetZ()-z)*(c->GetZ()-z) > 2.25 * 12 * sz2) continue;
1250
1251 Double_t h01 = GetTiltFactor(c);
9c9d2487 1252 chi2=t.GetPredictedChi2(c,h01);
a9814c09 1253
029cd327 1254 if (chi2 > maxChi2) continue;
1255 maxChi2=chi2;
a9814c09 1256 cl=c;
029cd327 1257 index=timeBin.GetIndex(i);
a9814c09 1258 }
1259 }
1260
1261 if (cl) {
1262 wYclosest = cl->GetY();
1263 wZclosest = cl->GetZ();
1264
1265 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1266 Double_t h01 = GetTiltFactor(cl);
3c625a9b 1267 Int_t det = cl->GetDetector();
1268 Int_t plane = fGeom->GetPlane(det);
1269
1270 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1271 //if(!t.Update(cl,maxChi2,index,h01)) {
029cd327 1272 if(!tryAgain--) return 0;
a9814c09 1273 }
029cd327 1274 else tryAgain=fMaxGap;
a9814c09 1275 }
1276 else {
029cd327 1277 if (tryAgain==0) break;
1278 tryAgain--;
9c9d2487 1279
1280 //if (minDY < 1000000 && isNewLayer)
1281 //cout << "\t" << nRefPlane << "\t" << "\t" << t.GetNRotate() << "\t" <<
029cd327 1282 // road << "\t" << minDY << "\t" << chi2 << "\t" << wChi2 << "\t" << maxChi2 << endl;
9c9d2487 1283
a9814c09 1284 }
1285
9c9d2487 1286 isNewLayer = kFALSE;
1287
a9814c09 1288 /*
1289 if((((Int_t) wTB)%15 == 0) || (((Int_t) wTB)%15 == 14)) {
1290
1291 printf(" %f", wIndex); //1
1292 printf(" %f", wTB); //2
1293 printf(" %f", wYrt); //3
1294 printf(" %f", wYclosest); //4
1295 printf(" %f", wYcorrect); //5
1296 printf(" %f", wYwindow); //6
1297 printf(" %f", wZrt); //7
1298 printf(" %f", wZclosest); //8
1299 printf(" %f", wZcorrect); //9
1300 printf(" %f", wZwindow); //10
1301 printf(" %f", wPx); //11
1302 printf(" %f", wPy); //12
1303 printf(" %f", wPz); //13
1304 printf(" %f", wSigmaC2*1000000); //14
1305 printf(" %f", wSigmaTgl2*1000); //15
1306 printf(" %f", wSigmaY2); //16
1307 // printf(" %f", wSigmaZ2); //17
1308 printf(" %f", wChi2); //17
1309 printf(" %f", wC); //18
1310 printf("\n");
1311 }
1312 */
5443e65e 1313 }
1314 }
1315 }
3c625a9b 1316 if (nr<outerTB)
1317 t.SetStop(kTRUE);
1318 else
1319 t.SetStop(kFALSE);
5443e65e 1320 return expectedNumberOfClusters;
9c9d2487 1321
1322
5443e65e 1323}
1324
1e9bb598 1325//---------------------------------------------------------------------------
1326Int_t AliTRDtracker::Refit(AliTRDtrack& t, Int_t rf)
1327{
1328 // Starting from current position on track=t this function tries
1329 // to extrapolate the track up to timeBin=0 and to reuse already
1330 // assigned clusters. Returns the number of clusters
1331 // expected to be found in sensitive layers
1332 // get indices of assigned clusters for each layer
1333 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
1334
1335 Int_t iCluster[90];
1336 for (Int_t i = 0; i < 90; i++) iCluster[i] = 0;
1337 for (Int_t i = 0; i < t.GetNumberOfClusters(); i++) {
1338 Int_t index = t.GetClusterIndex(i);
1339 AliTRDcluster *cl=(AliTRDcluster*) GetCluster(index);
1340 if (!cl) continue;
1341 Int_t detector=cl->GetDetector();
1342 Int_t localTimeBin=cl->GetLocalTimeBin();
1343 Int_t sector=fGeom->GetSector(detector);
1344 Int_t plane=fGeom->GetPlane(detector);
1345
1346 Int_t trackingSector = CookSectorIndex(sector);
1347
1348 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
1349 if(gtb < 0) continue;
1350 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
1351 iCluster[layer] = index;
1352 }
1353 t.ResetClusters();
1354
1355 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1356
1357 Double_t alpha=t.GetAlpha();
1358 alpha = TVector2::Phi_0_2pi(alpha);
1359
1360 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1361 Double_t radLength, rho, x, dx, y, ymax, z;
1362
1363 Int_t expectedNumberOfClusters = 0;
1364 Bool_t lookForCluster;
1365
1366 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1367
1368
1369 for (Int_t nr=fTrSec[0]->GetLayerNumber(t.GetX()); nr>rf; nr--) {
1370
1371 y = t.GetY(); z = t.GetZ();
1372
1373 // first propagate to the inner surface of the current time bin
1374 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1375 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2; y = t.GetY(); z = t.GetZ();
1376 if(!t.PropagateTo(x,radLength,rho)) break;
1377 y = t.GetY();
1378 ymax = x*TMath::Tan(0.5*alpha);
1379 if (y > ymax) {
1380 s = (s+1) % ns;
1381 if (!t.Rotate(alpha)) break;
1382 if(!t.PropagateTo(x,radLength,rho)) break;
1383 } else if (y <-ymax) {
1384 s = (s-1+ns) % ns;
1385 if (!t.Rotate(-alpha)) break;
1386 if(!t.PropagateTo(x,radLength,rho)) break;
1387 }
1388
1389 y = t.GetY(); z = t.GetZ();
1390
1391 // now propagate to the middle plane of the next time bin
1392 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
1393 x = fTrSec[s]->GetLayer(nr-1)->GetX(); y = t.GetY(); z = t.GetZ();
1394 if(!t.PropagateTo(x,radLength,rho)) break;
1395 y = t.GetY();
1396 ymax = x*TMath::Tan(0.5*alpha);
1397 if (y > ymax) {
1398 s = (s+1) % ns;
1399 if (!t.Rotate(alpha)) break;
1400 if(!t.PropagateTo(x,radLength,rho)) break;
1401 } else if (y <-ymax) {
1402 s = (s-1+ns) % ns;
1403 if (!t.Rotate(-alpha)) break;
1404 if(!t.PropagateTo(x,radLength,rho)) break;
1405 }
1406
1407 if(lookForCluster) expectedNumberOfClusters++;
1408
1409 // use assigned cluster
1410 if (!iCluster[nr-1]) continue;
1411 AliTRDcluster *cl=(AliTRDcluster*)GetCluster(iCluster[nr-1]);
1412 Double_t h01 = GetTiltFactor(cl);
1413 Double_t chi2=t.GetPredictedChi2(cl, h01);
1414 t.SetSampledEdx(cl->GetQ()/dx,t.GetNumberOfClusters());
1415 t.Update(cl,chi2,iCluster[nr-1],h01);
1416 }
1417
1418 return expectedNumberOfClusters;
1419}
1420
5443e65e 1421//___________________________________________________________________
1422
1423Int_t AliTRDtracker::PropagateToOuterPlane(AliTRDtrack& t, Double_t xToGo)
1424{
1425 // Starting from current radial position of track <t> this function
1426 // extrapolates the track up to radial position <xToGo>.
1427 // Returns 1 if track reaches the plane, and 0 otherwise
1428
1429 Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
1430
1431 Double_t alpha=t.GetAlpha();
1432
1433 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1434 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1435
1436 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1437
1438 Bool_t lookForCluster;
029cd327 1439 Double_t radLength, rho, x, dx, y, ymax, z;
5443e65e 1440
1441 x = t.GetX();
1442
1443 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
1444
1445 Int_t plToGo = fTrSec[0]->GetLayerNumber(xToGo);
1446
1447 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr<plToGo; nr++) {
1448
1449 y = t.GetY(); z = t.GetZ();
1450
1451 // first propagate to the outer surface of the current time bin
029cd327 1452 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
5443e65e 1453 x = fTrSec[s]->GetLayer(nr)->GetX()+dx/2; y = t.GetY(); z = t.GetZ();
029cd327 1454 if(!t.PropagateTo(x,radLength,rho)) return 0;
5443e65e 1455 y = t.GetY();
1456 ymax = x*TMath::Tan(0.5*alpha);
1457 if (y > ymax) {
1458 s = (s+1) % ns;
1459 if (!t.Rotate(alpha)) return 0;
1460 } else if (y <-ymax) {
1461 s = (s-1+ns) % ns;
1462 if (!t.Rotate(-alpha)) return 0;
1463 }
029cd327 1464 if(!t.PropagateTo(x,radLength,rho)) return 0;
5443e65e 1465
1466 y = t.GetY(); z = t.GetZ();
1467
1468 // now propagate to the middle plane of the next time bin
029cd327 1469 fTrSec[s]->GetLayer(nr+1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
5443e65e 1470 x = fTrSec[s]->GetLayer(nr+1)->GetX(); y = t.GetY(); z = t.GetZ();
029cd327 1471 if(!t.PropagateTo(x,radLength,rho)) return 0;
5443e65e 1472 y = t.GetY();
1473 ymax = x*TMath::Tan(0.5*alpha);
1474 if (y > ymax) {
1475 s = (s+1) % ns;
1476 if (!t.Rotate(alpha)) return 0;
1477 } else if (y <-ymax) {
1478 s = (s-1+ns) % ns;
1479 if (!t.Rotate(-alpha)) return 0;
1480 }
029cd327 1481 if(!t.PropagateTo(x,radLength,rho)) return 0;
5443e65e 1482 }
1483 return 1;
1484}
1485
1486//___________________________________________________________________
1487
1488Int_t AliTRDtracker::PropagateToTPC(AliTRDtrack& t)
1489{
1490 // Starting from current radial position of track <t> this function
1491 // extrapolates the track up to radial position of the outermost
1492 // padrow of the TPC.
1493 // Returns 1 if track reaches the TPC, and 0 otherwise
1494
c630aafd 1495 //Int_t ns=Int_t(2*TMath::Pi()/AliTRDgeometry::GetAlpha()+0.5);
5443e65e 1496
1497 Double_t alpha=t.GetAlpha();
c630aafd 1498 alpha = TVector2::Phi_0_2pi(alpha);
5443e65e 1499
1500 Int_t s=Int_t(alpha/AliTRDgeometry::GetAlpha())%AliTRDgeometry::kNsect;
1501
1502 Bool_t lookForCluster;
029cd327 1503 Double_t radLength, rho, x, dx, y, /*ymax,*/ z;
5443e65e 1504
1505 x = t.GetX();
1506
1507 alpha=AliTRDgeometry::GetAlpha(); // note: change in meaning
5443e65e 1508 Int_t plTPC = fTrSec[0]->GetLayerNumber(246.055);
1509
1510 for (Int_t nr=fTrSec[0]->GetLayerNumber(x); nr>plTPC; nr--) {
1511
9c9d2487 1512 y = t.GetY();
1513 z = t.GetZ();
5443e65e 1514
1515 // first propagate to the outer surface of the current time bin
029cd327 1516 fTrSec[s]->GetLayer(nr)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
9c9d2487 1517 x = fTrSec[s]->GetLayer(nr)->GetX()-dx/2;
1518
029cd327 1519 if(!t.PropagateTo(x,radLength,rho)) return 0;
9c9d2487 1520 AdjustSector(&t);
029cd327 1521 if(!t.PropagateTo(x,radLength,rho)) return 0;
5443e65e 1522
9c9d2487 1523 y = t.GetY();
1524 z = t.GetZ();
5443e65e 1525
1526 // now propagate to the middle plane of the next time bin
029cd327 1527 fTrSec[s]->GetLayer(nr-1)->GetPropagationParameters(y,z,dx,rho,radLength,lookForCluster);
9c9d2487 1528 x = fTrSec[s]->GetLayer(nr-1)->GetX();
1529
029cd327 1530 if(!t.PropagateTo(x,radLength,rho)) return 0;
9c9d2487 1531 AdjustSector(&t);
029cd327 1532 if(!t.PropagateTo(x,radLength,rho)) return 0;
c5a8e3df 1533 }
5443e65e 1534 return 1;
1535}
1536
c630aafd 1537//_____________________________________________________________________________
1538Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1539{
1540 // Fills clusters into TRD tracking_sectors
1541 // Note that the numbering scheme for the TRD tracking_sectors
1542 // differs from that of TRD sectors
1543
1544 if (ReadClusters(fClusters,cTree)) {
1545 Error("LoadClusters","Problem with reading the clusters !");
1546 return 1;
1547 }
1548 Int_t ncl=fClusters->GetEntriesFast();
b7a0917f 1549 fNclusters=ncl;
c630aafd 1550 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1551
1552 UInt_t index;
3c625a9b 1553 for (Int_t ichamber=0;ichamber<5;ichamber++)
1554 for (Int_t isector=0;isector<18;isector++){
1555 fHoles[ichamber][isector]=kTRUE;
1556 }
1557
1558
c630aafd 1559 while (ncl--) {
1560// printf("\r %d left ",ncl);
1561 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
029cd327 1562 Int_t detector=c->GetDetector();
1563 Int_t localTimeBin=c->GetLocalTimeBin();
c630aafd 1564 Int_t sector=fGeom->GetSector(detector);
1565 Int_t plane=fGeom->GetPlane(detector);
3c625a9b 1566
029cd327 1567 Int_t trackingSector = CookSectorIndex(sector);
3c625a9b 1568 if (c->GetLabel(0)>0){
1569 Int_t chamber = fGeom->GetChamber(detector);
1570 fHoles[chamber][trackingSector]=kFALSE;
1571 }
c630aafd 1572
029cd327 1573 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
c630aafd 1574 if(gtb < 0) continue;
029cd327 1575 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
c630aafd 1576
1577 index=ncl;
029cd327 1578 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
c630aafd 1579 }
7bed16a7 1580 // printf("\r\n");
3c625a9b 1581 //
1582 //
1583 /*
1584 for (Int_t isector=0;isector<18;isector++){
1585 for (Int_t ichamber=0;ichamber<5;ichamber++)
1586 if (fHoles[ichamber][isector]!=fGeom->IsHole(0,ichamber,17-isector))
1587 printf("Problem \t%d\t%d\t%d\t%d\n",isector,ichamber,fHoles[ichamber][isector],
1588 fGeom->IsHole(0,ichamber,17-isector));
1589 }
1590 */
c630aafd 1591 return 0;
1592}
1593
5443e65e 1594//_____________________________________________________________________________
b7a0917f 1595void AliTRDtracker::UnloadClusters()
5443e65e 1596{
1597 //
1598 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1599 //
1600
1601 Int_t i, nentr;
1602
1603 nentr = fClusters->GetEntriesFast();
1604 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
b7a0917f 1605 fNclusters = 0;
5443e65e 1606
1607 nentr = fSeeds->GetEntriesFast();
1608 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1609
1610 nentr = fTracks->GetEntriesFast();
1611 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1612
1613 Int_t nsec = AliTRDgeometry::kNsect;
1614
1615 for (i = 0; i < nsec; i++) {
1616 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1617 fTrSec[i]->GetLayer(pl)->Clear();
1618 }
1619 }
1620
1621}
1622
1623//__________________________________________________________________________
1624void AliTRDtracker::MakeSeeds(Int_t inner, Int_t outer, Int_t turn)
1625{
1626 // Creates track seeds using clusters in timeBins=i1,i2
1627
1628 if(turn > 2) {
1629 cerr<<"MakeSeeds: turn "<<turn<<" exceeds the limit of 2"<<endl;
1630 return;
1631 }
1632
1633 Double_t x[5], c[15];
029cd327 1634 Int_t maxSec=AliTRDgeometry::kNsect;
5443e65e 1635
1636 Double_t alpha=AliTRDgeometry::GetAlpha();
1637 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
1638 Double_t cs=cos(alpha), sn=sin(alpha);
1639 Double_t cs2=cos(2.*alpha), sn2=sin(2.*alpha);
1640
1641
1642 Int_t i2 = fTrSec[0]->GetLayerNumber(inner);
1643 Int_t i1 = fTrSec[0]->GetLayerNumber(outer);
1644
1645 Double_t x1 =fTrSec[0]->GetX(i1);
1646 Double_t xx2=fTrSec[0]->GetX(i2);
1647
029cd327 1648 for (Int_t ns=0; ns<maxSec; ns++) {
5443e65e 1649
029cd327 1650 Int_t nl2 = *(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
1651 Int_t nl=(*fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
5443e65e 1652 Int_t nm=(*fTrSec[ns]->GetLayer(i2));
029cd327 1653 Int_t nu=(*fTrSec[(ns+1)%maxSec]->GetLayer(i2));
1654 Int_t nu2=(*fTrSec[(ns+2)%maxSec]->GetLayer(i2));
5443e65e 1655
1656 AliTRDpropagationLayer& r1=*(fTrSec[ns]->GetLayer(i1));
1657
1658 for (Int_t is=0; is < r1; is++) {
1659 Double_t y1=r1[is]->GetY(), z1=r1[is]->GetZ();
1660
1661 for (Int_t js=0; js < nl2+nl+nm+nu+nu2; js++) {
a9814c09 1662
1663 const AliTRDcluster *cl;
1664 Double_t x2, y2, z2;
1665 Double_t x3=0., y3=0.;
1666
1667 if (js<nl2) {
1668 if(turn != 2) continue;
029cd327 1669 AliTRDpropagationLayer& r2=*(fTrSec[(ns-2+maxSec)%maxSec]->GetLayer(i2));
a9814c09 1670 cl=r2[js];
1671 y2=cl->GetY(); z2=cl->GetZ();
1672
1673 x2= xx2*cs2+y2*sn2;
1674 y2=-xx2*sn2+y2*cs2;
1675 }
1676 else if (js<nl2+nl) {
1677 if(turn != 1) continue;
029cd327 1678 AliTRDpropagationLayer& r2=*(fTrSec[(ns-1+maxSec)%maxSec]->GetLayer(i2));
a9814c09 1679 cl=r2[js-nl2];
1680 y2=cl->GetY(); z2=cl->GetZ();
1681
1682 x2= xx2*cs+y2*sn;
1683 y2=-xx2*sn+y2*cs;
1684 }
1685 else if (js<nl2+nl+nm) {
1686 if(turn != 1) continue;
1687 AliTRDpropagationLayer& r2=*(fTrSec[ns]->GetLayer(i2));
1688 cl=r2[js-nl2-nl];
1689 x2=xx2; y2=cl->GetY(); z2=cl->GetZ();
1690 }
1691 else if (js<nl2+nl+nm+nu) {
1692 if(turn != 1) continue;
029cd327 1693 AliTRDpropagationLayer& r2=*(fTrSec[(ns+1)%maxSec]->GetLayer(i2));
a9814c09 1694 cl=r2[js-nl2-nl-nm];
1695 y2=cl->GetY(); z2=cl->GetZ();
1696
1697 x2=xx2*cs-y2*sn;
1698 y2=xx2*sn+y2*cs;
1699 }
1700 else {
1701 if(turn != 2) continue;
029cd327 1702 AliTRDpropagationLayer& r2=*(fTrSec[(ns+2)%maxSec]->GetLayer(i2));
a9814c09 1703 cl=r2[js-nl2-nl-nm-nu];
1704 y2=cl->GetY(); z2=cl->GetZ();
1705
1706 x2=xx2*cs2-y2*sn2;
1707 y2=xx2*sn2+y2*cs2;
1708 }
1709
029cd327 1710 if(TMath::Abs(z1-z2) > fgkMaxSeedDeltaZ12) continue;
a9814c09 1711
1712 Double_t zz=z1 - z1/x1*(x1-x2);
1713
029cd327 1714 if (TMath::Abs(zz-z2)>fgkMaxSeedDeltaZ) continue;
a9814c09 1715
1716 Double_t d=(x2-x1)*(0.-y2)-(0.-x2)*(y2-y1);
1717 if (d==0.) {cerr<<"TRD MakeSeeds: Straight seed !\n"; continue;}
1718
1719 x[0]=y1;
1720 x[1]=z1;
1721 x[4]=f1trd(x1,y1,x2,y2,x3,y3);
1722
029cd327 1723 if (TMath::Abs(x[4]) > fgkMaxSeedC) continue;
a9814c09 1724
1725 x[2]=f2trd(x1,y1,x2,y2,x3,y3);
1726
1727 if (TMath::Abs(x[4]*x1-x[2]) >= 0.99999) continue;
1728
1729 x[3]=f3trd(x1,y1,x2,y2,z1,z2);
1730
029cd327 1731 if (TMath::Abs(x[3]) > fgkMaxSeedTan) continue;
a9814c09 1732
1733 Double_t a=asin(x[2]);
1734 Double_t zv=z1 - x[3]/x[4]*(a+asin(x[4]*x1-x[2]));
1735
029cd327 1736 if (TMath::Abs(zv)>fgkMaxSeedVertexZ) continue;
a9814c09 1737
1738 Double_t sy1=r1[is]->GetSigmaY2(), sz1=r1[is]->GetSigmaZ2();
1739 Double_t sy2=cl->GetSigmaY2(), sz2=cl->GetSigmaZ2();
029cd327 1740 Double_t sy3=fgkSeedErrorSY3, sy=fgkSeedErrorSY, sz=fgkSeedErrorSZ;
a9814c09 1741
1742 // Tilt changes
1743 Double_t h01 = GetTiltFactor(r1[is]);
029cd327 1744 Double_t xuFactor = 100.;
b8dc2353 1745 if(fNoTilt) {
1746 h01 = 0;
029cd327 1747 xuFactor = 1;
b8dc2353 1748 }
1749
fd621f36 1750 sy1=sy1+sz1*h01*h01;
b3a5a838 1751 Double_t syz=sz1*(-h01);
a9814c09 1752 // end of tilt changes
1753
b3a5a838 1754 Double_t f40=(f1trd(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
1755 Double_t f42=(f1trd(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
1756 Double_t f43=(f1trd(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
1757 Double_t f20=(f2trd(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
1758 Double_t f22=(f2trd(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
1759 Double_t f23=(f2trd(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
1760 Double_t f30=(f3trd(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
1761 Double_t f31=(f3trd(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
1762 Double_t f32=(f3trd(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
1763 Double_t f34=(f3trd(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
1764
a9814c09 1765
b3a5a838 1766 c[0]=sy1;
a9814c09 1767 // c[1]=0.; c[2]=sz1;
029cd327 1768 c[1]=syz; c[2]=sz1*xuFactor;
b3a5a838 1769 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
1770 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
1771 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
1772 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
1773 c[13]=f30*sy1*f40+f32*sy2*f42;
1774 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
a9814c09 1775
1776 UInt_t index=r1.GetIndex(is);
1777
1778 AliTRDtrack *track=new AliTRDtrack(r1[is],index,x,c,x1,ns*alpha+shift);
1779
1780 Int_t rc=FollowProlongation(*track, i2);
1781
1782 if ((rc < 1) ||
1783 (track->GetNumberOfClusters() <
029cd327 1784 (outer-inner)*fgkMinClustersInSeed)) delete track;
a9814c09 1785 else {
1786 fSeeds->AddLast(track); fNseeds++;
e24ea474 1787// cerr<<"\r found seed "<<fNseeds;
a9814c09 1788 }
5443e65e 1789 }
1790 }
1791 }
1792}
1793
1794//_____________________________________________________________________________
b7a0917f 1795Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
5443e65e 1796{
1797 //
a819a5f7 1798 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
1799 // from the file. The names of the cluster tree and branches
1800 // should match the ones used in AliTRDclusterizer::WriteClusters()
1801 //
029cd327 1802 TObjArray *clusterArray = new TObjArray(400);
5443e65e 1803
c630aafd 1804 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
1805 if (!branch) {
1806 Error("ReadClusters","Can't get the branch !");
1807 return 1;
1808 }
029cd327 1809 branch->SetAddress(&clusterArray);
5443e65e 1810
1811 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
19dd5b2f 1812 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
a819a5f7 1813
a819a5f7 1814 // Loop through all entries in the tree
eb187bed 1815 Int_t nbytes = 0;
a819a5f7 1816 AliTRDcluster *c = 0;
7bed16a7 1817 // printf("\n");
a819a5f7 1818
1819 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1820
1821 // Import the tree
5443e65e 1822 nbytes += ClusterTree->GetEvent(iEntry);
1823
a819a5f7 1824 // Get the number of points in the detector
029cd327 1825 Int_t nCluster = clusterArray->GetEntriesFast();
e24ea474 1826// printf("\r Read %d clusters from entry %d", nCluster, iEntry);
5443e65e 1827
a819a5f7 1828 // Loop through all TRD digits
1829 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
029cd327 1830 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
a819a5f7 1831 AliTRDcluster *co = new AliTRDcluster(*c);
1832 co->SetSigmaY2(c->GetSigmaY2() * fSY2corr);
5443e65e 1833 Int_t ltb = co->GetLocalTimeBin();
b8dc2353 1834 if(ltb == 19) co->SetSigmaZ2(c->GetSigmaZ2());
1835 else if(fNoTilt) co->SetSigmaZ2(c->GetSigmaZ2() * fSZ2corr);
a819a5f7 1836 array->AddLast(co);
029cd327 1837 delete clusterArray->RemoveAt(iCluster);
a819a5f7 1838 }
1839 }
1840
029cd327 1841 delete clusterArray;
5443e65e 1842
c630aafd 1843 return 0;
a819a5f7 1844}
1845
46d29e70 1846//__________________________________________________________________
029cd327 1847void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
1848{
1849 //
1850 // This cooks a label. Mmmmh, smells good...
1851 //
46d29e70 1852
1853 Int_t label=123456789, index, i, j;
5443e65e 1854 Int_t ncl=pt->GetNumberOfClusters();
029cd327 1855 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
5443e65e 1856
029cd327 1857 Bool_t labelAdded;
46d29e70 1858
029cd327 1859 // Int_t s[kRange][2];
1860 Int_t **s = new Int_t* [kRange];
1861 for (i=0; i<kRange; i++) {
d1b06c24 1862 s[i] = new Int_t[2];
1863 }
029cd327 1864 for (i=0; i<kRange; i++) {
46d29e70 1865 s[i][0]=-1;
1866 s[i][1]=0;
1867 }
1868
1869 Int_t t0,t1,t2;
1870 for (i=0; i<ncl; i++) {
5443e65e 1871 index=pt->GetClusterIndex(i);
46d29e70 1872 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
5443e65e 1873 t0=c->GetLabel(0);
1874 t1=c->GetLabel(1);
1875 t2=c->GetLabel(2);
46d29e70 1876 }
1877
1878 for (i=0; i<ncl; i++) {
5443e65e 1879 index=pt->GetClusterIndex(i);
46d29e70 1880 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1881 for (Int_t k=0; k<3; k++) {
5443e65e 1882 label=c->GetLabel(k);
029cd327 1883 labelAdded=kFALSE; j=0;
46d29e70 1884 if (label >= 0) {
029cd327 1885 while ( (!labelAdded) && ( j < kRange ) ) {
a9814c09 1886 if (s[j][0]==label || s[j][1]==0) {
1887 s[j][0]=label;
1888 s[j][1]=s[j][1]+1;
029cd327 1889 labelAdded=kTRUE;
a9814c09 1890 }
1891 j++;
1892 }
46d29e70 1893 }
1894 }
1895 }
1896
46d29e70 1897 Int_t max=0;
1898 label = -123456789;
1899
029cd327 1900 for (i=0; i<kRange; i++) {
46d29e70 1901 if (s[i][1]>max) {
1902 max=s[i][1]; label=s[i][0];
1903 }
1904 }
5443e65e 1905
029cd327 1906 for (i=0; i<kRange; i++) {
5443e65e 1907 delete []s[i];
1908 }
1909
d1b06c24 1910 delete []s;
5443e65e 1911
1912 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
1913
1914 pt->SetLabel(label);
1915
46d29e70 1916}
1917
c630aafd 1918
5443e65e 1919//__________________________________________________________________
029cd327 1920void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
1921{
1922 //
1923 // Use clusters, but don't abuse them!
1924 //
1925
5443e65e 1926 Int_t ncl=t->GetNumberOfClusters();
1927 for (Int_t i=from; i<ncl; i++) {
1928 Int_t index = t->GetClusterIndex(i);
1929 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
1930 c->Use();
1931 }
1932}
1933
1934
1935//_____________________________________________________________________
029cd327 1936Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
5443e65e 1937{
1938 // Parametrised "expected" error of the cluster reconstruction in Y
1939
1940 Double_t s = 0.08 * 0.08;
1941 return s;
1942}
1943
1944//_____________________________________________________________________
029cd327 1945Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
0a29d0f1 1946{
5443e65e 1947 // Parametrised "expected" error of the cluster reconstruction in Z
1948
a9814c09 1949 Double_t s = 9 * 9 /12.;
5443e65e 1950 return s;
1951}
1952
5443e65e 1953//_____________________________________________________________________
029cd327 1954Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
5443e65e 1955{
1956 //
029cd327 1957 // Returns radial position which corresponds to time bin <localTB>
5443e65e 1958 // in tracking sector <sector> and plane <plane>
1959 //
1960
029cd327 1961 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
5443e65e 1962 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
1963 return fTrSec[sector]->GetLayer(pl)->GetX();
1964
1965}
1966
c630aafd 1967
5443e65e 1968//_______________________________________________________
1969AliTRDtracker::AliTRDpropagationLayer::AliTRDpropagationLayer(Double_t x,
029cd327 1970 Double_t dx, Double_t rho, Double_t radLength, Int_t tbIndex)
5443e65e 1971{
0a29d0f1 1972 //
5443e65e 1973 // AliTRDpropagationLayer constructor
0a29d0f1 1974 //
46d29e70 1975
029cd327 1976 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
1977 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
5443e65e 1978
46d29e70 1979
029cd327 1980 for(Int_t i=0; i < (Int_t) kZones; i++) {
5443e65e 1981 fZc[i]=0; fZmax[i] = 0;
a819a5f7 1982 }
5443e65e 1983
1984 fYmax = 0;
1985
1986 if(fTimeBinIndex >= 0) {
029cd327 1987 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
1988 fIndex = new UInt_t[kMaxClusterPerTimeBin];
a819a5f7 1989 }
46d29e70 1990
3c625a9b 1991 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
5443e65e 1992 fHole = kFALSE;
1993 fHoleZc = 0;
1994 fHoleZmax = 0;
1995 fHoleYc = 0;
1996 fHoleYmax = 0;
1997 fHoleRho = 0;
1998 fHoleX0 = 0;
1999
2000}
2001
2002//_______________________________________________________
2003void AliTRDtracker::AliTRDpropagationLayer::SetHole(
a9814c09 2004 Double_t Zmax, Double_t Ymax, Double_t rho,
029cd327 2005 Double_t radLength, Double_t Yc, Double_t Zc)
5443e65e 2006{
2007 //
2008 // Sets hole in the layer
2009 //
5443e65e 2010 fHole = kTRUE;
2011 fHoleZc = Zc;
2012 fHoleZmax = Zmax;
2013 fHoleYc = Yc;
2014 fHoleYmax = Ymax;
2015 fHoleRho = rho;
029cd327 2016 fHoleX0 = radLength;
5443e65e 2017}
2018
46d29e70 2019
5443e65e 2020//_______________________________________________________
2021AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs, AliTRDparameter* par)
2022{
2023 //
2024 // AliTRDtrackingSector Constructor
2025 //
2026
2027 fGeom = geo;
2028 fPar = par;
2029 fGeomSector = gs;
2030 fTzeroShift = 0.13;
2031 fN = 0;
3c625a9b 2032 //
2033 // get holes description from geometry
2034 Bool_t holes[AliTRDgeometry::kNcham];
2035 //printf("sector\t%d\t",gs);
2036 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
2037 holes[icham] = fGeom->IsHole(0,icham,gs);
2038 //printf("%d",holes[icham]);
2039 }
2040 //printf("\n");
2041
029cd327 2042 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
5443e65e 2043
2044
2045 AliTRDpropagationLayer* ppl;
2046
029cd327 2047 Double_t x, xin, xout, dx, rho, radLength;
5443e65e 2048 Int_t steps;
46d29e70 2049
5443e65e 2050 // set time bins in the gas of the TPC
46d29e70 2051
5443e65e 2052 xin = 246.055; xout = 254.055; steps = 20; dx = (xout-xin)/steps;
029cd327 2053 rho = 0.9e-3; radLength = 28.94;
5443e65e 2054
2055 for(Int_t i=0; i<steps; i++) {
2056 x = xin + i*dx + dx/2;
029cd327 2057 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2058 InsertLayer(ppl);
46d29e70 2059 }
2060
5443e65e 2061 // set time bins in the outer field cage vessel
46d29e70 2062
029cd327 2063 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2064 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2065 InsertLayer(ppl);
46d29e70 2066
029cd327 2067 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2068 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2069 InsertLayer(ppl);
2070
029cd327 2071 dx = 2.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
5443e65e 2072 steps = 5; dx = (xout - xin)/steps;
2073 for(Int_t i=0; i<steps; i++) {
2074 x = xin + i*dx + dx/2;
029cd327 2075 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2076 InsertLayer(ppl);
2077 }
2078
029cd327 2079 dx = 0.02; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2080 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2081 InsertLayer(ppl);
2082
029cd327 2083 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2084 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2085 InsertLayer(ppl);
0a29d0f1 2086
5443e65e 2087
2088 // set time bins in CO2
2089
2090 xin = xout; xout = 275.0;
2091 steps = 50; dx = (xout - xin)/steps;
029cd327 2092 rho = 1.977e-3; radLength = 36.2;
5443e65e 2093
2094 for(Int_t i=0; i<steps; i++) {
2095 x = xin + i*dx + dx/2;
029cd327 2096 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2097 InsertLayer(ppl);
2098 }
2099
2100 // set time bins in the outer containment vessel
2101
029cd327 2102 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2103 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2104 InsertLayer(ppl);
2105
029cd327 2106 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2107 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2108 InsertLayer(ppl);
2109
029cd327 2110 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2111 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2112 InsertLayer(ppl);
2113
029cd327 2114 dx = 3.; xin = xout; xout = xin + dx; rho = 1.45*0.02; radLength = 41.28; // Nomex
5443e65e 2115 steps = 10; dx = (xout - xin)/steps;
2116 for(Int_t i=0; i<steps; i++) {
2117 x = xin + i*dx + dx/2;
029cd327 2118 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2119 InsertLayer(ppl);
2120 }
2121
029cd327 2122 dx = 0.06; xin = xout; xout = xin + dx; rho = 1.45; radLength = 44.86; // prepreg
2123 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2124 InsertLayer(ppl);
2125
029cd327 2126 dx = 50e-4; xin = xout; xout = xin + dx; rho = 1.71; radLength = 44.77; // Tedlar
2127 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2128 InsertLayer(ppl);
2129
029cd327 2130 dx = 50e-4; xin = xout; xout = xin + dx; rho = 2.7; radLength = 24.01; // Al
2131 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
5443e65e 2132 InsertLayer(ppl);
2133
2134 Double_t xtrd = (Double_t) fGeom->Rmin();
2135
2136 // add layers between TPC and TRD (Air temporarily)
2137 xin = xout; xout = xtrd;
2138 steps = 50; dx = (xout - xin)/steps;
029cd327 2139 rho = 1.2e-3; radLength = 36.66;
5443e65e 2140
2141 for(Int_t i=0; i<steps; i++) {
2142 x = xin + i*dx + dx/2;
029cd327 2143 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2144 InsertLayer(ppl);
2145 }
2146
2147
3c625a9b 2148 // Double_t alpha=AliTRDgeometry::GetAlpha();
5443e65e 2149
2150 // add layers for each of the planes
2151
2152 Double_t dxRo = (Double_t) fGeom->CroHght(); // Rohacell
2153 Double_t dxSpace = (Double_t) fGeom->Cspace(); // Spacing between planes
2154 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2155 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2156 Double_t dxRad = (Double_t) fGeom->CraHght(); // Radiator
2157 Double_t dxTEC = dxRad + dxDrift + dxAmp + dxRo;
2158 Double_t dxPlane = dxTEC + dxSpace;
2159
029cd327 2160 Int_t tb, tbIndex;
2161 const Int_t kNchambers = AliTRDgeometry::Ncham();
3c625a9b 2162 Double_t ymax = 0;
2163 //, holeYmax = 0;
2164 Double_t ymaxsensitive=0;
029cd327 2165 Double_t *zc = new Double_t[kNchambers];
2166 Double_t *zmax = new Double_t[kNchambers];
3c625a9b 2167 Double_t *zmaxsensitive = new Double_t[kNchambers];
2168 // Double_t holeZmax = 1000.; // the whole sector is missing
5443e65e 2169
5443e65e 2170 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
3c625a9b 2171 //
5443e65e 2172 // Radiator
2173 xin = xtrd + plane * dxPlane; xout = xin + dxRad;
029cd327 2174 steps = 12; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
5443e65e 2175 for(Int_t i=0; i<steps; i++) {
2176 x = xin + i*dx + dx/2;
3c625a9b 2177 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2178 InsertLayer(ppl);
2179 }
2180
3c625a9b 2181 ymax = fGeom->GetChamberWidth(plane)/2.;
2182 ymaxsensitive = (fPar->GetColPadSize(plane)*fPar->GetColMax(plane)-4)/2.;
2183
029cd327 2184 for(Int_t ch = 0; ch < kNchambers; ch++) {
2185 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
5443e65e 2186 Float_t pad = fPar->GetRowPadSize(plane,ch,0);
2187 Float_t row0 = fPar->GetRow0(plane,ch,0);
2188 Int_t nPads = fPar->GetRowMax(plane,ch,0);
857b3eb0 2189 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
3c625a9b 2190 // zc[ch] = (pad * nPads)/2 + row0 - pad/2;
2191 zc[ch] = (pad * nPads)/2 + row0;
2192 //zc[ch] = row0+zmax[ch]-AliTRDgeometry::RpadW();
2193
5443e65e 2194 }
2195
2196 dx = fPar->GetTimeBinSize();
029cd327 2197 rho = 0.00295 * 0.85; radLength = 11.0;
5443e65e 2198
2199 Double_t x0 = (Double_t) fPar->GetTime0(plane);
2200 Double_t xbottom = x0 - dxDrift;
2201 Double_t xtop = x0 + dxAmp;
3c625a9b 2202 //
5443e65e 2203 // Amplification region
5443e65e 2204 steps = (Int_t) (dxAmp/dx);
2205
2206 for(tb = 0; tb < steps; tb++) {
2207 x = x0 + tb * dx + dx/2;
029cd327 2208 tbIndex = CookTimeBinIndex(plane, -tb-1);
2209 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
3c625a9b 2210 ppl->SetYmax(ymax,ymaxsensitive);
2211 ppl->SetZ(zc, zmax, zmaxsensitive);
2212 ppl->SetHoles(holes);
5443e65e 2213 InsertLayer(ppl);
2214 }
029cd327 2215 tbIndex = CookTimeBinIndex(plane, -steps);
5443e65e 2216 x = (x + dx/2 + xtop)/2;
2217 dx = 2*(xtop-x);
029cd327 2218 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
3c625a9b 2219 ppl->SetYmax(ymax,ymaxsensitive);
2220 ppl->SetZ(zc, zmax,zmaxsensitive);
2221 ppl->SetHoles(holes);
5443e65e 2222 InsertLayer(ppl);
2223
2224 // Drift region
2225 dx = fPar->GetTimeBinSize();
2226 steps = (Int_t) (dxDrift/dx);
2227
2228 for(tb = 0; tb < steps; tb++) {
2229 x = x0 - tb * dx - dx/2;
029cd327 2230 tbIndex = CookTimeBinIndex(plane, tb);
5443e65e 2231
029cd327 2232 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
3c625a9b 2233 ppl->SetYmax(ymax,ymaxsensitive);
2234 ppl->SetZ(zc, zmax, zmaxsensitive);
2235 ppl->SetHoles(holes);
5443e65e 2236 InsertLayer(ppl);
2237 }
029cd327 2238 tbIndex = CookTimeBinIndex(plane, steps);
5443e65e 2239 x = (x - dx/2 + xbottom)/2;
2240 dx = 2*(x-xbottom);
029cd327 2241 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex);
3c625a9b 2242 ppl->SetYmax(ymax,ymaxsensitive);
2243 ppl->SetZ(zc, zmax, zmaxsensitive);
2244 ppl->SetHoles(holes);
5443e65e 2245 InsertLayer(ppl);
2246
2247 // Pad Plane
029cd327 2248 xin = xtop; dx = 0.025; xout = xin + dx; rho = 1.7; radLength = 33.0;
2249 ppl = new AliTRDpropagationLayer(xin+dx/2,dx,rho,radLength,-1);
3c625a9b 2250 ppl->SetYmax(ymax,ymaxsensitive);
2251 ppl->SetZ(zc, zmax,zmax);
2252 ppl->SetHoles(holes);
5443e65e 2253 InsertLayer(ppl);
2254
2255 // Rohacell
2256 xin = xout; xout = xtrd + (plane + 1) * dxPlane - dxSpace;
029cd327 2257 steps = 5; dx = (xout - xin)/steps; rho = 0.074; radLength = 40.6;
5443e65e 2258 for(Int_t i=0; i<steps; i++) {
2259 x = xin + i*dx + dx/2;
029cd327 2260 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
3c625a9b 2261 ppl->SetYmax(ymax,ymaxsensitive);
2262 ppl->SetZ(zc, zmax,zmax);
2263 ppl->SetHoles(holes);
5443e65e 2264 InsertLayer(ppl);
2265 }
2266
2267 // Space between the chambers, air
2268 xin = xout; xout = xtrd + (plane + 1) * dxPlane;
029cd327 2269 steps = 5; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
5443e65e 2270 for(Int_t i=0; i<steps; i++) {
2271 x = xin + i*dx + dx/2;
029cd327 2272 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2273 InsertLayer(ppl);
2274 }
2275 }
2276
2277 // Space between the TRD and RICH
2278 Double_t xRICH = 500.;
2279 xin = xout; xout = xRICH;
029cd327 2280 steps = 200; dx = (xout - xin)/steps; rho = 1.29e-3; radLength = 36.66;
5443e65e 2281 for(Int_t i=0; i<steps; i++) {
2282 x = xin + i*dx + dx/2;
029cd327 2283 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,-1);
5443e65e 2284 InsertLayer(ppl);
2285 }
2286
2287 MapTimeBinLayers();
029cd327 2288 delete [] zc;
2289 delete [] zmax;
5443e65e 2290
2291}
2292
2293//______________________________________________________
2294
029cd327 2295Int_t AliTRDtracker::AliTRDtrackingSector::CookTimeBinIndex(Int_t plane, Int_t localTB) const
5443e65e 2296{
2297 //
2298 // depending on the digitization parameters calculates "global"
029cd327 2299 // time bin index for timebin <localTB> in plane <plane>
5443e65e 2300 //
2301
2302 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
2303 Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
2304 Double_t dx = (Double_t) fPar->GetTimeBinSize();
2305
2306 Int_t tbAmp = fPar->GetTimeBefore();
2307 Int_t maxAmp = (Int_t) ((dxAmp+0.000001)/dx);
b3a5a838 2308 if(kTRUE) maxAmp = 0; // intentional until we change parameter class
5443e65e 2309 Int_t tbDrift = fPar->GetTimeMax();
2310 Int_t maxDrift = (Int_t) ((dxDrift+0.000001)/dx);
2311
029cd327 2312 Int_t tbPerPlane = TMath::Min(tbAmp,maxAmp) + TMath::Min(tbDrift,maxDrift);
5443e65e 2313
029cd327 2314 Int_t gtb = (plane+1) * tbPerPlane - localTB - 1 - TMath::Min(tbAmp,maxAmp);
5443e65e 2315
029cd327 2316 if((localTB < 0) &&
2317 (TMath::Abs(localTB) > TMath::Min(tbAmp,maxAmp))) return -1;
2318 if(localTB >= TMath::Min(tbDrift,maxDrift)) return -1;
5443e65e 2319
2320 return gtb;
2321
2322
2323}
2324
2325//______________________________________________________
2326
2327void AliTRDtracker::AliTRDtrackingSector::MapTimeBinLayers()
2328{
2329 //
2330 // For all sensitive time bins sets corresponding layer index
2331 // in the array fTimeBins
2332 //
2333
2334 Int_t index;
2335
2336 for(Int_t i = 0; i < fN; i++) {
2337 index = fLayers[i]->GetTimeBinIndex();
2338
2339 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2340
2341 if(index < 0) continue;
029cd327 2342 if(index >= (Int_t) kMaxTimeBinIndex) {
5443e65e 2343 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2344 printf(" index %d exceeds allowed maximum of %d!\n",
029cd327 2345 index, kMaxTimeBinIndex-1);
5443e65e 2346 continue;
2347 }
2348 fTimeBinIndex[index] = i;
2349 }
2350
2351 Double_t x1, dx1, x2, dx2, gap;
2352
2353 for(Int_t i = 0; i < fN-1; i++) {
2354 x1 = fLayers[i]->GetX();
2355 dx1 = fLayers[i]->GetdX();
2356 x2 = fLayers[i+1]->GetX();
2357 dx2 = fLayers[i+1]->GetdX();
2358 gap = (x2 - dx2/2) - (x1 + dx1/2);
2359 if(gap < -0.01) {
2360 printf("*** warning: layers %d and %d are overlayed:\n",i,i+1);
2361 printf(" %f + %f + %f > %f\n", x1, dx1/2, dx2/2, x2);
2362 }
2363 if(gap > 0.01) {
2364 printf("*** warning: layers %d and %d have a large gap:\n",i,i+1);
2365 printf(" (%f - %f) - (%f + %f) = %f\n",
a9814c09 2366 x2, dx2/2, x1, dx1, gap);
5443e65e 2367 }
2368 }
2369}
2370
2371
2372//______________________________________________________
2373
2374
2375Int_t AliTRDtracker::AliTRDtrackingSector::GetLayerNumber(Double_t x) const
2376{
2377 //
2378 // Returns the number of time bin which in radial position is closest to <x>
2379 //
2380
2381 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2382 if(x <= fLayers[0]->GetX()) return 0;
2383
2384 Int_t b=0, e=fN-1, m=(b+e)/2;
2385 for (; b<e; m=(b+e)/2) {
2386 if (x > fLayers[m]->GetX()) b=m+1;
2387 else e=m;
2388 }
2389 if(TMath::Abs(x - fLayers[m]->GetX()) >
2390 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
2391 else return m;
2392
2393}
2394
2395//______________________________________________________
2396
2397Int_t AliTRDtracker::AliTRDtrackingSector::GetInnerTimeBin() const
2398{
2399 //
2400 // Returns number of the innermost SENSITIVE propagation layer
2401 //
2402
2403 return GetLayerNumber(0);
2404}
2405
2406//______________________________________________________
2407
2408Int_t AliTRDtracker::AliTRDtrackingSector::GetOuterTimeBin() const
2409{
2410 //
2411 // Returns number of the outermost SENSITIVE time bin
2412 //
2413
2414 return GetLayerNumber(GetNumberOfTimeBins() - 1);
46d29e70 2415}
2416
5443e65e 2417//______________________________________________________
2418
2419Int_t AliTRDtracker::AliTRDtrackingSector::GetNumberOfTimeBins() const
2420{
2421 //
2422 // Returns number of SENSITIVE time bins
2423 //
2424
2425 Int_t tb, layer;
029cd327 2426 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
5443e65e 2427 layer = GetLayerNumber(tb);
2428 if(layer>=0) break;
2429 }
2430 return tb+1;
2431}
2432
2433//______________________________________________________
2434
2435void AliTRDtracker::AliTRDtrackingSector::InsertLayer(AliTRDpropagationLayer* pl)
2436{
2437 //
2438 // Insert layer <pl> in fLayers array.
2439 // Layers are sorted according to X coordinate.
2440
029cd327 2441 if ( fN == ((Int_t) kMaxLayersPerSector)) {
5443e65e 2442 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2443 return;
2444 }
2445 if (fN==0) {fLayers[fN++] = pl; return;}
2446 Int_t i=Find(pl->GetX());
2447
2448 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2449 fLayers[i]=pl; fN++;
2450
2451}
2452
2453//______________________________________________________
2454
2455Int_t AliTRDtracker::AliTRDtrackingSector::Find(Double_t x) const
2456{
2457 //
2458 // Returns index of the propagation layer nearest to X
2459 //
2460
2461 if (x <= fLayers[0]->GetX()) return 0;
2462 if (x > fLayers[fN-1]->GetX()) return fN;
2463 Int_t b=0, e=fN-1, m=(b+e)/2;
2464 for (; b<e; m=(b+e)/2) {
2465 if (x > fLayers[m]->GetX()) b=m+1;
2466 else e=m;
2467 }
2468 return m;
2469}
2470
3c625a9b 2471//______________________________________________________
2472void AliTRDtracker::AliTRDpropagationLayer::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
2473{
2474 //
2475 // set centers and the width of sectors
2476 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2477 fZc[icham] = center[icham];
2478 fZmax[icham] = w[icham];
2479 fZmaxSensitive[icham] = wsensitive[icham];
2480 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
2481 }
2482}
5443e65e 2483//______________________________________________________
2484
3c625a9b 2485void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2486{
2487 //
2488 // set centers and the width of sectors
2489 fHole = kFALSE;
2490 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2491 fIsHole[icham] = holes[icham];
2492 if (holes[icham]) fHole = kTRUE;
2493 }
2494}
2495
2496
2497
5443e65e 2498void AliTRDtracker::AliTRDpropagationLayer::GetPropagationParameters(
029cd327 2499 Double_t y, Double_t z, Double_t &dx, Double_t &rho, Double_t &radLength,
a9814c09 2500 Bool_t &lookForCluster) const
5443e65e 2501{
2502 //
029cd327 2503 // Returns radial step <dx>, density <rho>, rad. length <radLength>,
5443e65e 2504 // and sensitivity <lookForCluster> in point <y,z>
2505 //
2506
2507 dx = fdX;
2508 rho = fRho;
029cd327 2509 radLength = fX0;
5443e65e 2510 lookForCluster = kFALSE;
3c625a9b 2511 //
2512 // check dead regions in sensitive volume
5443e65e 2513 if(fTimeBinIndex >= 0) {
3c625a9b 2514 //
2515 Int_t zone=-1;
029cd327 2516 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
3c625a9b 2517 if (TMath::Abs(z - fZc[ch]) < fZmaxSensitive[ch]){
2518 zone = ch;
2519 lookForCluster = !(fIsHole[zone]);
2520 if(TMath::Abs(y) > fYmaxSensitive){
2521 lookForCluster = kFALSE;
2522 }
2523 if (fIsHole[zone]) {
2524 //if hole
2525 rho = 1.29e-3;
2526 radLength = 36.66;
2527 }
2528 }
5443e65e 2529 }
3c625a9b 2530 return;
5443e65e 2531 }
3c625a9b 2532 //
2533 //
5443e65e 2534 // check hole
3c625a9b 2535 if (fHole==kFALSE) return;
2536 //
2537 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2538 if (TMath::Abs(z - fZc[ch]) < fZmax[ch]){
2539 if (fIsHole[ch]) {
2540 //if hole
2541 rho = 1.29e-3;
2542 radLength = 36.66;
2543 }
2544 }
2545 }
5443e65e 2546 return;
2547}
2548
3c625a9b 2549Int_t AliTRDtracker::AliTRDpropagationLayer::GetZone( Double_t z) const
2550{
2551 //
2552 //
2553 if (fTimeBinIndex < 0) return -20; //unknown
2554 Int_t zone=-10; // dead zone
2555 for(Int_t ch = 0; ch < (Int_t) kZones; ch++) {
2556 if(TMath::Abs(z - fZc[ch]) < fZmax[ch])
2557 zone = ch;
2558 }
2559 return zone;
2560}
2561
2562
5443e65e 2563//______________________________________________________
2564
2565void AliTRDtracker::AliTRDpropagationLayer::InsertCluster(AliTRDcluster* c,
a9814c09 2566 UInt_t index) {
5443e65e 2567
2568// Insert cluster in cluster array.
2569// Clusters are sorted according to Y coordinate.
2570
2571 if(fTimeBinIndex < 0) {
2572 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2573 return;
2574 }
2575
029cd327 2576 if (fN== (Int_t) kMaxClusterPerTimeBin) {
5443e65e 2577 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2578 return;
2579 }
2580 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2581 Int_t i=Find(c->GetY());
2582 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2583 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2584 fIndex[i]=index; fClusters[i]=c; fN++;
2585}
2586
2587//______________________________________________________
2588
2589Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Double_t y) const {
2590
2591// Returns index of the cluster nearest in Y
2592
2593 if (y <= fClusters[0]->GetY()) return 0;
2594 if (y > fClusters[fN-1]->GetY()) return fN;
2595 Int_t b=0, e=fN-1, m=(b+e)/2;
2596 for (; b<e; m=(b+e)/2) {
2597 if (y > fClusters[m]->GetY()) b=m+1;
2598 else e=m;
2599 }
2600 return m;
2601}
2602
fd621f36 2603//---------------------------------------------------------
5443e65e 2604
fd621f36 2605Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c) {
2606//
2607// Returns correction factor for tilted pads geometry
2608//
5443e65e 2609
fd621f36 2610 Double_t h01 = sin(TMath::Pi() / 180.0 * fPar->GetTiltingAngle());
2611 Int_t det = c->GetDetector();
2612 Int_t plane = fGeom->GetPlane(det);
b3a5a838 2613
c5a8e3df 2614 //if((plane == 1) || (plane == 3) || (plane == 5)) h01=-h01;
2615 if((plane == 0) || (plane == 2) || (plane == 4)) h01=-h01;
b8dc2353 2616
2617 if(fNoTilt) h01 = 0;
fd621f36 2618
2619 return h01;
2620}
5443e65e 2621
c630aafd 2622
eab5961e 2623void AliTRDtracker::CookdEdxTimBin(AliTRDtrack& TRDtrack)
2624{
2625 // *** ADDED TO GET MORE INFORMATION FOR TRD PID ---- PS
2626 // This is setting fdEdxPlane and fTimBinPlane
2627 // Sums up the charge in each plane for track TRDtrack and also get the
2628 // Time bin for Max. Cluster
2629 // Prashant Shukla (shukla@physi.uni-heidelberg.de)
2630
2631 // const Int_t kNPlane = AliTRDgeometry::Nplan();
2632 // const Int_t kNPlane = 6;
2633 Double_t clscharge[kNPlane], maxclscharge[kNPlane];
2634 Int_t nCluster[kNPlane], timebin[kNPlane];
2635
2636 //Initialization of cluster charge per plane.
2637 for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2638 clscharge[iPlane] = 0.0;
2639 nCluster[iPlane] = 0;
2640 timebin[iPlane] = -1;
2641 maxclscharge[iPlane] = 0.0;
2642 }
2643
2644 // Loop through all clusters associated to track TRDtrack
2645 Int_t nClus = TRDtrack.GetNumberOfClusters(); // from Kalmantrack
2646 for (Int_t iClus = 0; iClus < nClus; iClus++) {
2647 Double_t charge = TRDtrack.GetClusterdQdl(iClus);
2648 Int_t index = TRDtrack.GetClusterIndex(iClus);
2649 AliTRDcluster *TRDcluster = (AliTRDcluster *) GetCluster(index);
2650 if (!TRDcluster) continue;
2651 Int_t tb = TRDcluster->GetLocalTimeBin();
2652 if (!tb) continue;
2653 Int_t detector = TRDcluster->GetDetector();
2654 Int_t iPlane = fGeom->GetPlane(detector);
2655 clscharge[iPlane] = clscharge[iPlane]+charge;
2656 if(charge > maxclscharge[iPlane]) {
2657 maxclscharge[iPlane] = charge;
2658 timebin[iPlane] = tb;
2659 }
2660 nCluster[iPlane]++;
2661 } // end of loop over cluster
2662
2663 // Setting the fdEdxPlane and fTimBinPlane variabales
2664 Double_t Total_ch = 0;
2665 for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2666 if (nCluster[iPlane]) clscharge[iPlane] /= nCluster[iPlane];
2667 TRDtrack.SetPIDsignals(clscharge[iPlane], iPlane);
2668 TRDtrack.SetPIDTimBin(timebin[iPlane], iPlane);
2669 Total_ch= Total_ch+clscharge[iPlane];
2670 }
2671 // Int_t i;
2672 // Int_t nc=TRDtrack.GetNumberOfClusters();
2673 // Float_t dedx=0;
2674 // for (i=0; i<nc; i++) dedx += TRDtrack.GetClusterdQdl(i);
2675 // dedx /= nc;
2676 // for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2677 // TRDtrack.SetPIDsignals(dedx, iPlane);
2678 // TRDtrack.SetPIDTimBin(timbin[iPlane], iPlane);
2679 // }
2680
2681} // end of function
2682
c630aafd 2683
2684
2685
2686