]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDtracker.cxx
Updated TRD PID for electrons. Additional information from START
[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// //
c6f438c0 20// The standard TRD tracker //
21// Based on Kalman filltering approach //
029cd327 22// //
23///////////////////////////////////////////////////////////////////////////////
24
a2cb5b3d 25#include <Riostream.h>
46d29e70 26#include <TFile.h>
46d29e70 27#include <TBranch.h>
5443e65e 28#include <TTree.h>
c630aafd 29#include <TObjArray.h>
46d29e70 30
46d29e70 31#include "AliTRDgeometry.h"
a5cadd36 32#include "AliTRDpadPlane.h"
dde59437 33#include "AliTRDgeometryFull.h"
46d29e70 34#include "AliTRDcluster.h"
35#include "AliTRDtrack.h"
75bd7f81 36#include "AliTRDseed.h"
b7a0917f 37#include "AliESD.h"
46d29e70 38
3551db50 39#include "AliTRDcalibDB.h"
40#include "AliTRDCommonParam.h"
41
7ad19338 42#include "TTreeStream.h"
43#include "TGraph.h"
46d29e70 44#include "AliTRDtracker.h"
69b55c55 45#include "TLinearFitter.h"
46#include "AliRieman.h"
3551db50 47#include "AliTrackPointArray.h"
48#include "AliAlignObj.h"
7b580082 49#include "AliTRDReconstructor.h"
46d29e70 50
51ClassImp(AliTRDtracker)
bbf92647 52
029cd327 53 const Float_t AliTRDtracker::fgkMinClustersInTrack = 0.5;
029cd327 54 const Float_t AliTRDtracker::fgkLabelFraction = 0.8;
029cd327 55 const Double_t AliTRDtracker::fgkMaxChi2 = 12.;
75bd7f81 56 const Double_t AliTRDtracker::fgkMaxSnp = 0.95; // correspond to tan = 3
57 const Double_t AliTRDtracker::fgkMaxStep = 2.; // maximal step size in propagation
7ad19338 58
75bd7f81 59//_____________________________________________________________________________
89f05372 60AliTRDtracker::AliTRDtracker():AliTracker(),
61 fGeom(0),
89f05372 62 fNclusters(0),
63 fClusters(0),
64 fNseeds(0),
65 fSeeds(0),
66 fNtracks(0),
67 fTracks(0),
89f05372 68 fTimeBinsPerPlane(0),
89f05372 69 fAddTRDseeds(kFALSE),
70 fNoTilt(kFALSE)
71{
75bd7f81 72 //
b7a0917f 73 // Default constructor
75bd7f81 74 //
b7a0917f 75
89f05372 76 for(Int_t i=0;i<kTrackingSectors;i++) fTrSec[i]=0;
77 for(Int_t j=0;j<5;j++)
75bd7f81 78 for(Int_t k=0;k<18;k++) fHoles[j][k]=kFALSE;
7ad19338 79 fDebugStreamer = 0;
75bd7f81 80
89f05372 81}
75bd7f81 82
83//_____________________________________________________________________________
c630aafd 84AliTRDtracker::AliTRDtracker(const TFile *geomfile):AliTracker()
46d29e70 85{
5443e65e 86 //
87 // Main constructor
88 //
b8dc2353 89
5443e65e 90 fAddTRDseeds = kFALSE;
5443e65e 91 fGeom = NULL;
b8dc2353 92 fNoTilt = kFALSE;
5443e65e 93
94 TDirectory *savedir=gDirectory;
95 TFile *in=(TFile*)geomfile;
96 if (!in->IsOpen()) {
97 printf("AliTRDtracker::AliTRDtracker(): geometry file is not open!\n");
dde59437 98 printf(" FULL TRD geometry and DEFAULT TRD parameter will be used\n");
5443e65e 99 }
100 else {
101 in->cd();
5443e65e 102 fGeom = (AliTRDgeometry*) in->Get("TRDgeometry");
5443e65e 103 }
46d29e70 104
5443e65e 105 if(fGeom) {
7c1698cb 106 // printf("Found geometry version %d on file \n", fGeom->IsVersion());
5443e65e 107 }
108 else {
c630aafd 109 printf("AliTRDtracker::AliTRDtracker(): can't find TRD geometry!\n");
6d45eaef 110
dde59437 111 fGeom = new AliTRDgeometryFull();
6d45eaef 112
c630aafd 113 }
c6f438c0 114 fGeom->ReadGeoMatrices();
c630aafd 115
5443e65e 116 savedir->cd();
46d29e70 117
46d29e70 118 fNclusters = 0;
119 fClusters = new TObjArray(2000);
120 fNseeds = 0;
5443e65e 121 fSeeds = new TObjArray(2000);
46d29e70 122 fNtracks = 0;
5443e65e 123 fTracks = new TObjArray(1000);
a819a5f7 124
029cd327 125 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
126 Int_t trS = CookSectorIndex(geomS);
59393e34 127 fTrSec[trS] = new AliTRDtrackingSector(fGeom, geomS);
3c625a9b 128 for (Int_t icham=0;icham<AliTRDgeometry::kNcham; icham++){
129 fHoles[icham][trS]=fGeom->IsHole(0,icham,geomS);
130 }
5443e65e 131 }
3551db50 132 AliTRDpadPlane *padPlane = AliTRDCommonParam::Instance()->GetPadPlane(0,0);
7ad19338 133 Float_t tiltAngle = TMath::Abs(padPlane->GetTiltingAngle());
029cd327 134 if(tiltAngle < 0.1) {
b8dc2353 135 fNoTilt = kTRUE;
136 }
137
59393e34 138 fTimeBinsPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
46d29e70 139
7ad19338 140 fDebugStreamer = new TTreeSRedirector("TRDdebug.root");
0a29d0f1 141
9c9d2487 142 savedir->cd();
75bd7f81 143
5443e65e 144}
46d29e70 145
75bd7f81 146//_____________________________________________________________________________
5443e65e 147AliTRDtracker::~AliTRDtracker()
46d29e70 148{
029cd327 149 //
150 // Destructor of AliTRDtracker
151 //
152
89f05372 153 if (fClusters) {
154 fClusters->Delete();
155 delete fClusters;
156 }
157 if (fTracks) {
158 fTracks->Delete();
159 delete fTracks;
160 }
161 if (fSeeds) {
162 fSeeds->Delete();
163 delete fSeeds;
164 }
5443e65e 165 delete fGeom;
0a29d0f1 166
029cd327 167 for(Int_t geomS = 0; geomS < kTrackingSectors; geomS++) {
168 delete fTrSec[geomS];
5443e65e 169 }
7ad19338 170 if (fDebugStreamer) {
171 //fDebugStreamer->Close();
172 delete fDebugStreamer;
173 }
9c9d2487 174
75bd7f81 175}
59393e34 176
75bd7f81 177//_____________________________________________________________________________
178Int_t AliTRDtracker::LocalToGlobalID(Int_t lid)
179{
59393e34 180 //
75bd7f81 181 // Transform internal TRD ID to global detector ID
59393e34 182 //
75bd7f81 183
59393e34 184 Int_t isector = fGeom->GetSector(lid);
185 Int_t ichamber= fGeom->GetChamber(lid);
186 Int_t iplan = fGeom->GetPlane(lid);
187 //
188 AliAlignObj::ELayerID iLayer = AliAlignObj::kTRD1;
189 switch (iplan) {
190 case 0:
191 iLayer = AliAlignObj::kTRD1;
192 break;
193 case 1:
194 iLayer = AliAlignObj::kTRD2;
195 break;
196 case 2:
197 iLayer = AliAlignObj::kTRD3;
198 break;
199 case 3:
200 iLayer = AliAlignObj::kTRD4;
201 break;
202 case 4:
203 iLayer = AliAlignObj::kTRD5;
204 break;
205 case 5:
206 iLayer = AliAlignObj::kTRD6;
207 break;
208 };
209 Int_t modId = isector*fGeom->Ncham()+ichamber;
210 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer,modId);
75bd7f81 211
59393e34 212 return volid;
75bd7f81 213
59393e34 214}
215
75bd7f81 216//_____________________________________________________________________________
217Int_t AliTRDtracker::GlobalToLocalID(Int_t gid)
218{
59393e34 219 //
75bd7f81 220 // Transform global detector ID to local detector ID
59393e34 221 //
75bd7f81 222
59393e34 223 Int_t modId=0;
224 AliAlignObj::ELayerID layerId = AliAlignObj::VolUIDToLayer(gid, modId);
225 Int_t isector = modId/fGeom->Ncham();
226 Int_t ichamber = modId%fGeom->Ncham();
227 Int_t iLayer = -1;
228 switch (layerId) {
229 case AliAlignObj::kTRD1:
230 iLayer = 0;
231 break;
232 case AliAlignObj::kTRD2:
233 iLayer = 1;
234 break;
235 case AliAlignObj::kTRD3:
236 iLayer = 2;
237 break;
238 case AliAlignObj::kTRD4:
239 iLayer = 3;
240 break;
241 case AliAlignObj::kTRD5:
242 iLayer = 4;
243 break;
244 case AliAlignObj::kTRD6:
245 iLayer = 5;
246 break;
247 default:
248 iLayer =-1;
249 }
250 if (iLayer<0) return -1;
251 Int_t lid = fGeom->GetDetector(iLayer,ichamber,isector);
75bd7f81 252
59393e34 253 return lid;
59393e34 254
75bd7f81 255}
59393e34 256
75bd7f81 257//_____________________________________________________________________________
258Bool_t AliTRDtracker::Transform(AliTRDcluster * cluster)
259{
59393e34 260 //
75bd7f81 261 // Transform something ... whatever ...
c6f438c0 262 //
75bd7f81 263
c6f438c0 264 const Double_t kX0shift = 2.52; // magic constants for geo manager transformation
265 const Double_t kX0shift5 = 3.05; //
266 //
59393e34 267 //
268 // apply alignment and calibration to transform cluster
269 //
270 //
c6f438c0 271 Int_t detector = cluster->GetDetector();
272 Int_t plane = fGeom->GetPlane(cluster->GetDetector());
273 Int_t chamber = fGeom->GetChamber(cluster->GetDetector());
274 Int_t sector = fGeom->GetSector(cluster->GetDetector());
275
59393e34 276 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
277 Double_t driftX = TMath::Max(cluster->GetX()-dxAmp*0.5,0.); // drift distance
278 //
59393e34 279 // ExB correction
280 //
281 Double_t vdrift = AliTRDcalibDB::Instance()->GetVdrift(cluster->GetDetector(),0,0);
282 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(vdrift);
283 //
c6f438c0 284 AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
285 AliTRDpadPlane * padPlane = commonParam->GetPadPlane(plane,chamber);
286 Double_t zshiftIdeal = 0.5*(padPlane->GetRow0()+padPlane->GetRowEnd());
287 Double_t localPos[3], localPosTracker[3];
288 localPos[0] = -cluster->GetX();
289 localPos[1] = cluster->GetY() - driftX*exB;
290 localPos[2] = cluster->GetZ() -zshiftIdeal;
291 //
292 cluster->SetY(cluster->GetY() - driftX*exB);
293 Double_t xplane = (Double_t) AliTRDgeometry::GetTime0(plane);
294 cluster->SetX(xplane- cluster->GetX());
295 //
296 TGeoHMatrix * matrix = fGeom->GetCorrectionMatrix(cluster->GetDetector());
297 if (!matrix){
298 // no matrix found - if somebody used geometry with holes
299 AliError("Invalid Geometry - Default Geometry used\n");
300 return kTRUE;
301 }
302 matrix->LocalToMaster(localPos, localPosTracker);
303 //
304 //
305 //
d337ef8d 306 if (AliTRDReconstructor::StreamLevel()>1){
c6f438c0 307 (*fDebugStreamer)<<"Transform"<<
308 "Cl.="<<cluster<<
309 "matrix.="<<matrix<<
310 "Detector="<<detector<<
311 "Sector="<<sector<<
312 "Plane="<<plane<<
313 "Chamber="<<chamber<<
314 "lx0="<<localPosTracker[0]<<
315 "ly0="<<localPosTracker[1]<<
316 "lz0="<<localPosTracker[2]<<
317 "\n";
318 }
319 //
320 if (plane==5)
321 cluster->SetX(localPosTracker[0]+kX0shift5);
322 else
323 cluster->SetX(localPosTracker[0]+kX0shift);
324
325 cluster->SetY(localPosTracker[1]);
326 cluster->SetZ(localPosTracker[2]);
75bd7f81 327
59393e34 328 return kTRUE;
75bd7f81 329
59393e34 330}
331
75bd7f81 332//_____________________________________________________________________________
333// Bool_t AliTRDtracker::Transform(AliTRDcluster * cluster)
334//{
c6f438c0 335// //
336// //
337// const Double_t kDriftCorrection = 1.01; // drift coeficient correction
338// const Double_t kTime0Cor = 0.32; // time0 correction
339// //
340// const Double_t kX0shift = 2.52;
341// const Double_t kX0shift5 = 3.05;
342
343// //
344// // apply alignment and calibration to transform cluster
345// //
346// //
347// Int_t detector = cluster->GetDetector();
348// Int_t plane = fGeom->GetPlane(cluster->GetDetector());
349// Int_t chamber = fGeom->GetChamber(cluster->GetDetector());
350// Int_t sector = fGeom->GetSector(cluster->GetDetector());
351
352// Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
353// Double_t driftX = TMath::Max(cluster->GetX()-dxAmp*0.5,0.); // drift distance
354// //
355// // ExB correction
356// //
357// Double_t vdrift = AliTRDcalibDB::Instance()->GetVdrift(cluster->GetDetector(),0,0);
358// Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(vdrift);
359// //
360
361// AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
362// AliTRDpadPlane * padPlane = commonParam->GetPadPlane(plane,chamber);
363// Double_t zshiftIdeal = 0.5*(padPlane->GetRow0()+padPlane->GetRowEnd());
364// Double_t localPos[3], globalPos[3], localPosTracker[3], localPosTracker2[3];
365// localPos[2] = -cluster->GetX();
366// localPos[0] = cluster->GetY() - driftX*exB;
367// localPos[1] = cluster->GetZ() -zshiftIdeal;
368// TGeoHMatrix * matrix = fGeom->GetGeoMatrix(cluster->GetDetector());
369// matrix->LocalToMaster(localPos, globalPos);
370
371// Double_t sectorAngle = 20.*(sector%18)+10;
372// TGeoHMatrix rotSector;
373// rotSector.RotateZ(sectorAngle);
374// rotSector.LocalToMaster(globalPos, localPosTracker);
375// //
376// //
377// TGeoHMatrix matrix2(*matrix);
378// matrix2.MultiplyLeft(&rotSector);
379// matrix2.LocalToMaster(localPos,localPosTracker2);
380// //
381// //
382// //
383// cluster->SetY(cluster->GetY() - driftX*exB);
384// Double_t xplane = (Double_t) AliTRDgeometry::GetTime0(plane);
385// cluster->SetX(xplane- kDriftCorrection*(cluster->GetX()-kTime0Cor));
386// (*fDebugStreamer)<<"Transform"<<
387// "Cl.="<<cluster<<
388// "matrix.="<<matrix<<
389// "matrix2.="<<&matrix2<<
390// "Detector="<<detector<<
391// "Sector="<<sector<<
392// "Plane="<<plane<<
393// "Chamber="<<chamber<<
394// "lx0="<<localPosTracker[0]<<
395// "ly0="<<localPosTracker[1]<<
396// "lz0="<<localPosTracker[2]<<
397// "lx2="<<localPosTracker2[0]<<
398// "ly2="<<localPosTracker2[1]<<
399// "lz2="<<localPosTracker2[2]<<
400// "\n";
401// //
402// if (plane==5)
403// cluster->SetX(localPosTracker[0]+kX0shift5);
404// else
405// cluster->SetX(localPosTracker[0]+kX0shift);
406
407// cluster->SetY(localPosTracker[1]);
408// cluster->SetZ(localPosTracker[2]);
409// return kTRUE;
410// }
411
75bd7f81 412//_____________________________________________________________________________
413Bool_t AliTRDtracker::AdjustSector(AliTRDtrack *track)
414{
9c9d2487 415 //
416 // Rotates the track when necessary
417 //
418
419 Double_t alpha = AliTRDgeometry::GetAlpha();
420 Double_t y = track->GetY();
421 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
422
c630aafd 423 //Int_t ns = AliTRDgeometry::kNsect;
9c9d2487 424 //Int_t s=Int_t(track->GetAlpha()/alpha)%ns;
425
426 if (y > ymax) {
427 //s = (s+1) % ns;
428 if (!track->Rotate(alpha)) return kFALSE;
429 } else if (y <-ymax) {
430 //s = (s-1+ns) % ns;
431 if (!track->Rotate(-alpha)) return kFALSE;
432 }
433
434 return kTRUE;
9c9d2487 435
75bd7f81 436}
46e2d86c 437
75bd7f81 438//_____________________________________________________________________________
439AliTRDcluster *AliTRDtracker::GetCluster(AliTRDtrack *track, Int_t plane
440 , Int_t timebin, UInt_t &index)
441{
46e2d86c 442 //
75bd7f81 443 // Try to find cluster in the backup list
46e2d86c 444 //
75bd7f81 445
46e2d86c 446 AliTRDcluster * cl =0;
ef7253ac 447 Int_t *indexes = track->GetBackupIndexes();
46e2d86c 448 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
449 if (indexes[i]==0) break;
450 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
451 if (!cli) break;
452 if (cli->GetLocalTimeBin()!=timebin) continue;
453 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
454 if (iplane==plane) {
455 cl = cli;
7ad19338 456 index = indexes[i];
46e2d86c 457 break;
458 }
459 }
75bd7f81 460
46e2d86c 461 return cl;
46e2d86c 462
75bd7f81 463}
3c625a9b 464
75bd7f81 465//_____________________________________________________________________________
466Int_t AliTRDtracker::GetLastPlane(AliTRDtrack * track)
467{
3c625a9b 468 //
75bd7f81 469 // Return last updated plane
470 //
471
3c625a9b 472 Int_t lastplane=0;
ef7253ac 473 Int_t *indexes = track->GetBackupIndexes();
3c625a9b 474 for (UInt_t i=0;i<kMaxTimeBinIndex;i++){
475 AliTRDcluster * cli = (AliTRDcluster*)fClusters->UncheckedAt(indexes[i]);
476 if (!cli) break;
477 Int_t iplane = fGeom->GetPlane(cli->GetDetector());
478 if (iplane>lastplane) {
479 lastplane = iplane;
480 }
481 }
75bd7f81 482
3c625a9b 483 return lastplane;
75bd7f81 484
3c625a9b 485}
75bd7f81 486
487//_____________________________________________________________________________
c630aafd 488Int_t AliTRDtracker::Clusters2Tracks(AliESD* event)
489{
490 //
491 // Finds tracks within the TRD. The ESD event is expected to contain seeds
492 // at the outer part of the TRD. The seeds
493 // are found within the TRD if fAddTRDseeds is TRUE.
494 // The tracks are propagated to the innermost time bin
495 // of the TRD and the ESD event is updated
496 //
497
498 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
029cd327 499 Float_t foundMin = fgkMinClustersInTrack * timeBins;
c630aafd 500 Int_t nseed = 0;
501 Int_t found = 0;
7b580082 502 // Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
c630aafd 503
504 Int_t n = event->GetNumberOfTracks();
505 for (Int_t i=0; i<n; i++) {
506 AliESDtrack* seed=event->GetTrack(i);
507 ULong_t status=seed->GetStatus();
508 if ( (status & AliESDtrack::kTRDout ) == 0 ) continue;
509 if ( (status & AliESDtrack::kTRDin) != 0 ) continue;
510 nseed++;
7ad19338 511
c630aafd 512 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
46e2d86c 513 //seed2->ResetCovariance();
c630aafd 514 AliTRDtrack *pt = new AliTRDtrack(*seed2,seed2->GetAlpha());
515 AliTRDtrack &t=*pt;
7b580082 516 FollowProlongation(t);
c630aafd 517 if (t.GetNumberOfClusters() >= foundMin) {
518 UseClusters(&t);
029cd327 519 CookLabel(pt, 1-fgkLabelFraction);
c630aafd 520 // t.CookdEdx();
521 }
522 found++;
523// cout<<found<<'\r';
524
59393e34 525 Double_t xTPC = 250;
526 if (PropagateToX(t,xTPC,fgkMaxStep)) {
c630aafd 527 seed->UpdateTrackParams(pt, AliESDtrack::kTRDin);
528 }
529 delete seed2;
530 delete pt;
531 }
532
533 cout<<"Number of loaded seeds: "<<nseed<<endl;
534 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
535
536 // after tracks from loaded seeds are found and the corresponding
537 // clusters are used, look for additional seeds from TRD
c630aafd 538
c630aafd 539
540 cout<<"Total number of found tracks: "<<found<<endl;
541
542 return 0;
75bd7f81 543
c630aafd 544}
5443e65e 545
c630aafd 546//_____________________________________________________________________________
75bd7f81 547Int_t AliTRDtracker::PropagateBack(AliESD* event)
548{
c630aafd 549 //
550 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
551 // backpropagated by the TPC tracker. Each seed is first propagated
552 // to the TRD, and then its prolongation is searched in the TRD.
553 // If sufficiently long continuation of the track is found in the TRD
554 // the track is updated, otherwise it's stored as originaly defined
555 // by the TPC tracker.
556 //
557
558 Int_t found=0;
c5a8e3df 559 Float_t foundMin = 20;
c630aafd 560 Int_t n = event->GetNumberOfTracks();
4f1c04d3 561 //
562 //Sort tracks
563 Float_t *quality =new Float_t[n];
564 Int_t *index =new Int_t[n];
c630aafd 565 for (Int_t i=0; i<n; i++) {
566 AliESDtrack* seed=event->GetTrack(i);
4f1c04d3 567 Double_t covariance[15];
568 seed->GetExternalCovariance(covariance);
569 quality[i] = covariance[0]+covariance[2];
570 }
571 TMath::Sort(n,quality,index,kFALSE);
572 //
573 for (Int_t i=0; i<n; i++) {
574 // AliESDtrack* seed=event->GetTrack(i);
575 AliESDtrack* seed=event->GetTrack(index[i]);
576
c630aafd 577 ULong_t status=seed->GetStatus();
578 if ( (status & AliESDtrack::kTPCout ) == 0 ) continue;
579 if ( (status & AliESDtrack::kTRDout) != 0 ) continue;
580
581 Int_t lbl = seed->GetLabel();
582 AliTRDtrack *track = new AliTRDtrack(*seed);
583 track->SetSeedLabel(lbl);
f4e9508c 584 seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup); //make backup
c630aafd 585 fNseeds++;
4f1c04d3 586 Float_t p4 = track->GetC();
7bed16a7 587 //
f6625211 588 Int_t expectedClr = FollowBackProlongation(*track);
4f1c04d3 589 if (TMath::Abs(track->GetC()-p4)/TMath::Abs(p4)<0.2 || TMath::Abs(track->GetPt())>0.8 ) {
590 //
591 //make backup for back propagation
592 //
593 Int_t foundClr = track->GetNumberOfClusters();
594 if (foundClr >= foundMin) {
595 track->CookdEdx();
8979685e 596 CookdEdxTimBin(*track);
4f1c04d3 597 CookLabel(track, 1-fgkLabelFraction);
69b55c55 598 if (track->GetBackupTrack()) UseClusters(track->GetBackupTrack());
4f1c04d3 599 if(track->GetChi2()/track->GetNumberOfClusters()<4) { // sign only gold tracks
600 if (seed->GetKinkIndex(0)==0&&TMath::Abs(track->GetPt())<1.5 ) UseClusters(track);
601 }
602 Bool_t isGold = kFALSE;
603
604 if (track->GetChi2()/track->GetNumberOfClusters()<5) { //full gold track
7ad19338 605 // seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
606 if (track->GetBackupTrack()) seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
4f1c04d3 607 isGold = kTRUE;
608 }
609 if (!isGold && track->GetNCross()==0&&track->GetChi2()/track->GetNumberOfClusters()<7){ //almost gold track
7ad19338 610 // seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
611 if (track->GetBackupTrack()) seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
f4e9508c 612 isGold = kTRUE;
613 }
4f1c04d3 614 if (!isGold && track->GetBackupTrack()){
615 if (track->GetBackupTrack()->GetNumberOfClusters()>foundMin&&
616 (track->GetBackupTrack()->GetChi2()/(track->GetBackupTrack()->GetNumberOfClusters()+1))<7){
617 seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
618 isGold = kTRUE;
619 }
620 }
7ad19338 621 if (track->StatusForTOF()>0 &&track->fNCross==0 && Float_t(track->fN)/Float_t(track->fNExpected)>0.4){
59393e34 622 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
7ad19338 623 }
16d9fbba 624 }
c630aafd 625 }
8979685e 626 // Debug part of tracking
627 TTreeSRedirector& cstream = *fDebugStreamer;
628 Int_t eventNr = event->GetEventNumber();
d337ef8d 629 if (AliTRDReconstructor::StreamLevel()>0){
630 if (track->GetBackupTrack()){
631 cstream<<"Tracks"<<
632 "EventNr="<<eventNr<<
633 "ESD.="<<seed<<
634 "trd.="<<track<<
635 "trdback.="<<track->GetBackupTrack()<<
636 "\n";
637 }else{
638 cstream<<"Tracks"<<
639 "EventNr="<<eventNr<<
640 "ESD.="<<seed<<
641 "trd.="<<track<<
642 "trdback.="<<track<<
643 "\n";
644 }
8979685e 645 }
646 //
647 //Propagation to the TOF (I.Belikov)
3c625a9b 648 if (track->GetStop()==kFALSE){
4f1c04d3 649
b94f0a96 650 Double_t xtof=371.;
3c625a9b 651 Double_t c2=track->GetC()*xtof - track->GetEta();
4f1c04d3 652 if (TMath::Abs(c2)>=0.99) {
c5a8e3df 653 delete track;
654 continue;
655 }
59393e34 656 Double_t xTOF0 = 370. ;
657 PropagateToX(*track,xTOF0,fgkMaxStep);
4f1c04d3 658 //
659 //energy losses taken to the account - check one more time
660 c2=track->GetC()*xtof - track->GetEta();
661 if (TMath::Abs(c2)>=0.99) {
662 delete track;
663 continue;
664 }
665
7bed16a7 666 //
3c625a9b 667 Double_t ymax=xtof*TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
668 Double_t y=track->GetYat(xtof);
669 if (y > ymax) {
7ac6fa52 670 if (!track->Rotate(AliTRDgeometry::GetAlpha())) {
671 delete track;
7bed16a7 672 continue;
7ac6fa52 673 }
3c625a9b 674 } else if (y <-ymax) {
7ac6fa52 675 if (!track->Rotate(-AliTRDgeometry::GetAlpha())) {
676 delete track;
7bed16a7 677 continue;
7ac6fa52 678 }
3c625a9b 679 }
680
681 if (track->PropagateTo(xtof)) {
eab5961e 682 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
f122c485 683 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
6d45eaef 684 for (Int_t j=0;j<AliESDtrack::kNSlice;j++) {
685 seed->SetTRDsignals(track->GetPIDsignals(i,j),i,j);
686 }
687 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
eab5961e 688 }
7ad19338 689 // seed->SetTRDtrack(new AliTRDtrack(*track));
3c625a9b 690 if (track->GetNumberOfClusters()>foundMin) found++;
691 }
692 }else{
693 if (track->GetNumberOfClusters()>15&&track->GetNumberOfClusters()>0.5*expectedClr){
694 seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
16d9fbba 695 //seed->SetStatus(AliESDtrack::kTRDStop);
f122c485 696 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
6d45eaef 697 for (Int_t j=0;j<AliESDtrack::kNSlice;j++) {
698 seed->SetTRDsignals(track->GetPIDsignals(i,j),i,j);
699 }
700 seed->SetTRDTimBin(track->GetPIDTimBin(i),i);
eab5961e 701 }
7ad19338 702 //seed->SetTRDtrack(new AliTRDtrack(*track));
3c625a9b 703 found++;
704 }
1e9bb598 705 }
7ad19338 706 seed->SetTRDQuality(track->StatusForTOF());
8979685e 707 seed->SetTRDBudget(track->fBudget[0]);
708
d9b8978b 709 delete track;
7ad19338 710 //
1e9bb598 711 //End of propagation to the TOF
3c625a9b 712 //if (foundClr>foundMin)
713 // seed->UpdateTrackParams(track, AliESDtrack::kTRDout);
714
c630aafd 715
716 }
717
718 cerr<<"Number of seeds: "<<fNseeds<<endl;
719 cerr<<"Number of back propagated TRD tracks: "<<found<<endl;
69b55c55 720
7b580082 721 if (AliTRDReconstructor::SeedingOn()) MakeSeedsMI(3,5,event); //new seeding
7ad19338 722
1e9bb598 723 fSeeds->Clear(); fNseeds=0;
4f1c04d3 724 delete [] index;
725 delete [] quality;
726
1e9bb598 727 return 0;
728
729}
730
731//_____________________________________________________________________________
732Int_t AliTRDtracker::RefitInward(AliESD* event)
733{
734 //
735 // Refits tracks within the TRD. The ESD event is expected to contain seeds
736 // at the outer part of the TRD.
737 // The tracks are propagated to the innermost time bin
738 // of the TRD and the ESD event is updated
739 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
740 //
741
742 Int_t timeBins = fTrSec[0]->GetNumberOfTimeBins();
743 Float_t foundMin = fgkMinClustersInTrack * timeBins;
744 Int_t nseed = 0;
745 Int_t found = 0;
7b580082 746 // Int_t innerTB = fTrSec[0]->GetInnerTimeBin();
4f1c04d3 747 AliTRDtrack seed2;
1e9bb598 748
749 Int_t n = event->GetNumberOfTracks();
750 for (Int_t i=0; i<n; i++) {
751 AliESDtrack* seed=event->GetTrack(i);
4f1c04d3 752 new(&seed2) AliTRDtrack(*seed);
753 if (seed2.GetX()<270){
754 seed->UpdateTrackParams(&seed2, AliESDtrack::kTRDbackup); // backup TPC track - only update
f4e9508c 755 continue;
756 }
757
1e9bb598 758 ULong_t status=seed->GetStatus();
0dd7d129 759 if ( (status & AliESDtrack::kTRDout ) == 0 ) {
0dd7d129 760 continue;
761 }
762 if ( (status & AliESDtrack::kTRDin) != 0 ) {
0dd7d129 763 continue;
764 }
f4e9508c 765 nseed++;
7ad19338 766// if (1/seed2.Get1Pt()>1.5&& seed2.GetX()>260.) {
767// Double_t oldx = seed2.GetX();
768// seed2.PropagateTo(500.);
769// seed2.ResetCovariance(1.);
770// seed2.PropagateTo(oldx);
771// }
772// else{
773// seed2.ResetCovariance(5.);
774// }
4f1c04d3 775
776 AliTRDtrack *pt = new AliTRDtrack(seed2,seed2.GetAlpha());
ef7253ac 777 Int_t * indexes2 = seed2.GetIndexes();
f122c485 778 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
6d45eaef 779 for (Int_t j=0;j<AliESDtrack::kNSlice;j++) {
780 pt->SetPIDsignals(seed2.GetPIDsignals(i,j),i,j);
781 }
7ad19338 782 pt->SetPIDTimBin(seed2.GetPIDTimBin(i),i);
783 }
eab5961e 784
ef7253ac 785 Int_t * indexes3 = pt->GetBackupIndexes();
46e2d86c 786 for (Int_t i=0;i<200;i++) {
787 if (indexes2[i]==0) break;
788 indexes3[i] = indexes2[i];
789 }
790 //AliTRDtrack *pt = seed2;
1e9bb598 791 AliTRDtrack &t=*pt;
7b580082 792 FollowProlongation(t);
1e9bb598 793 if (t.GetNumberOfClusters() >= foundMin) {
46e2d86c 794 // UseClusters(&t);
795 //CookLabel(pt, 1-fgkLabelFraction);
7ad19338 796 t.CookdEdx();
797 CookdEdxTimBin(t);
1e9bb598 798 }
799 found++;
800// cout<<found<<'\r';
59393e34 801 Double_t xTPC = 250;
802 if(PropagateToX(t,xTPC,fgkMaxStep)) {
0fa7dfa7 803 seed->UpdateTrackParams(pt, AliESDtrack::kTRDrefit);
f122c485 804 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
6d45eaef 805 for (Int_t j=0;j<AliESDtrack::kNSlice;j++) {
806 seed->SetTRDsignals(pt->GetPIDsignals(i,j),i,j);
807 }
7ad19338 808 seed->SetTRDTimBin(pt->GetPIDTimBin(i),i);
809 }
7bed16a7 810 }else{
811 //if not prolongation to TPC - propagate without update
812 AliTRDtrack* seed2 = new AliTRDtrack(*seed);
813 seed2->ResetCovariance(5.);
814 AliTRDtrack *pt2 = new AliTRDtrack(*seed2,seed2->GetAlpha());
815 delete seed2;
59393e34 816 if (PropagateToX(*pt2,xTPC,fgkMaxStep)) {
3551db50 817 //pt2->CookdEdx(0.,1.);
818 pt2->CookdEdx( ); // Modification by PS
eab5961e 819 CookdEdxTimBin(*pt2);
7bed16a7 820 seed->UpdateTrackParams(pt2, AliESDtrack::kTRDrefit);
f122c485 821 for (Int_t i=0;i<AliESDtrack::kNPlane;i++) {
6d45eaef 822 for (Int_t j=0;j<AliESDtrack::kNSlice;j++) {
823 seed->SetTRDsignals(pt2->GetPIDsignals(i,j),i,j);
824 }
7ad19338 825 seed->SetTRDTimBin(pt2->GetPIDTimBin(i),i);
826 }
eab5961e 827 }
7bed16a7 828 delete pt2;
1e9bb598 829 }
1e9bb598 830 delete pt;
eab5961e 831 }
1e9bb598 832
833 cout<<"Number of loaded seeds: "<<nseed<<endl;
834 cout<<"Number of found tracks from loaded seeds: "<<found<<endl;
835
c630aafd 836 return 0;
837
838}
839
75bd7f81 840//_____________________________________________________________________________
7b580082 841Int_t AliTRDtracker::FollowProlongation(AliTRDtrack& t)
8979685e 842{
75bd7f81 843 //
8979685e 844 // Starting from current position on track=t this function tries
845 // to extrapolate the track up to timeBin=0 and to confirm prolongation
846 // if a close cluster is found. Returns the number of clusters
847 // expected to be found in sensitive layers
848 // GeoManager used to estimate mean density
75bd7f81 849 //
850
8979685e 851 Int_t sector;
852 Int_t lastplane = GetLastPlane(&t);
3551db50 853 Double_t radLength = 0.0;
854 Double_t rho = 0.0;
8979685e 855 Int_t expectedNumberOfClusters = 0;
8979685e 856 //
857 //
7b580082 858 //
859 for (Int_t iplane = lastplane; iplane>=0; iplane--){
8979685e 860 //
7b580082 861 Int_t row0 = GetGlobalTimeBin(0, iplane,GetTimeBinsPerPlane()-1);
862 Int_t rowlast = GetGlobalTimeBin(0, iplane,0);
8979685e 863 //
7b580082 864 // propagate track close to the plane if neccessary
865 //
866 Double_t currentx = fTrSec[0]->GetLayer(rowlast)->GetX();
59393e34 867 if (currentx < -fgkMaxStep +t.GetX()){
868 //propagate closer to chamber - safety space fgkMaxStep
869 if (!PropagateToX(t, currentx+fgkMaxStep, fgkMaxStep)) break;
8979685e 870 }
8979685e 871 if (!AdjustSector(&t)) break;
59393e34 872 //
7b580082 873 // get material budget
8979685e 874 //
7b580082 875 Double_t xyz0[3],xyz1[3],param[7],x,y,z;
876 t.GetGlobalXYZ(xyz0[0],xyz0[1],xyz0[2]); //starting global position
877 // end global position
878 x = fTrSec[0]->GetLayer(row0)->GetX();
879 if (!t.GetProlongation(x,y,z)) break;
880 xyz1[0] = x*TMath::Cos(t.GetAlpha())-y*TMath::Sin(t.GetAlpha());
881 xyz1[1] = +x*TMath::Sin(t.GetAlpha())+y*TMath::Cos(t.GetAlpha());
882 xyz1[2] = z;
883 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
884 rho = param[0];
885 radLength = param[1]; // get mean propagation parameters
8979685e 886 //
7b580082 887 // propagate nad update
8979685e 888 //
8979685e 889 sector = t.GetSector();
7b580082 890 // for (Int_t itime=GetTimeBinsPerPlane()-1;itime>=0;itime--) {
891 for (Int_t itime=0 ;itime<GetTimeBinsPerPlane();itime++) {
892 Int_t ilayer = GetGlobalTimeBin(0, iplane,itime);
8979685e 893 expectedNumberOfClusters++;
894 t.fNExpected++;
895 if (t.fX>345) t.fNExpectedLast++;
896 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(ilayer));
897 AliTRDcluster *cl=0;
898 UInt_t index=0;
899 Double_t maxChi2=fgkMaxChi2;
8979685e 900 x = timeBin.GetX();
8979685e 901 if (timeBin) {
902 AliTRDcluster * cl0 = timeBin[0];
903 if (!cl0) continue; // no clusters in given time bin
904 Int_t plane = fGeom->GetPlane(cl0->GetDetector());
905 if (plane>lastplane) continue;
906 Int_t timebin = cl0->GetLocalTimeBin();
907 AliTRDcluster * cl2= GetCluster(&t,plane, timebin,index);
908 //
909 if (cl2) {
910 cl =cl2;
911 Double_t h01 = GetTiltFactor(cl);
912 maxChi2=t.GetPredictedChi2(cl,h01);
7b580082 913 }
8979685e 914 if (cl) {
915 // if (cl->GetNPads()<5)
59393e34 916 Double_t dxsample = timeBin.GetdX();
917 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
8979685e 918 Double_t h01 = GetTiltFactor(cl);
919 Int_t det = cl->GetDetector();
920 Int_t plane = fGeom->GetPlane(det);
921 if (t.fX>345){
922 t.fNLast++;
923 t.fChi2Last+=maxChi2;
924 }
59393e34 925 Double_t xcluster = cl->GetX();
926 t.PropagateTo(xcluster,radLength,rho);
8979685e 927 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
7b580082 928 }
8979685e 929 }
930 }
7b580082 931 }
8979685e 932 }
8979685e 933
75bd7f81 934 return expectedNumberOfClusters;
69b55c55 935
75bd7f81 936}
7b580082 937
75bd7f81 938//_____________________________________________________________________________
f6625211 939Int_t AliTRDtracker::FollowBackProlongation(AliTRDtrack& t)
8979685e 940{
75bd7f81 941 //
8979685e 942 // Starting from current radial position of track <t> this function
943 // extrapolates the track up to outer timebin and in the sensitive
944 // layers confirms prolongation if a close cluster is found.
945 // Returns the number of clusters expected to be found in sensitive layers
946 // Use GEO manager for material Description
75bd7f81 947 //
59393e34 948
8979685e 949 Int_t sector;
950 Int_t clusters[1000];
951 for (Int_t i=0;i<1000;i++) clusters[i]=-1;
3551db50 952 Double_t radLength = 0.0;
953 Double_t rho = 0.0;
8979685e 954 Int_t expectedNumberOfClusters = 0;
8979685e 955 Float_t ratio0=0;
956 AliTRDtracklet tracklet;
957 //
958 //
f122c485 959 for (Int_t iplane = 0; iplane<AliESDtrack::kNPlane; iplane++){
7b580082 960 Int_t row0 = GetGlobalTimeBin(0, iplane,GetTimeBinsPerPlane()-1);
961 Int_t rowlast = GetGlobalTimeBin(0, iplane,0);
8979685e 962 //
7b580082 963 Double_t currentx = fTrSec[0]->GetLayer(row0)->GetX();
964 if (currentx<t.GetX()) continue;
965 //
966 // propagate closer to chamber if neccessary
8979685e 967 //
59393e34 968 if (currentx > fgkMaxStep +t.GetX()){
59393e34 969 if (!PropagateToX(t, currentx-fgkMaxStep, fgkMaxStep)) break;
8979685e 970 }
8979685e 971 if (!AdjustSector(&t)) break;
59393e34 972 if (TMath::Abs(t.GetSnp())>fgkMaxSnp) break;
8979685e 973 //
7b580082 974 // get material budget inside of chamber
59393e34 975 //
7b580082 976 Double_t xyz0[3],xyz1[3],param[7],x,y,z;
977 t.GetGlobalXYZ(xyz0[0],xyz0[1],xyz0[2]); //starting global position
978 // end global position
979 x = fTrSec[0]->GetLayer(rowlast)->GetX();
980 if (!t.GetProlongation(x,y,z)) break;
981 xyz1[0] = x*TMath::Cos(t.GetAlpha())-y*TMath::Sin(t.GetAlpha());
982 xyz1[1] = +x*TMath::Sin(t.GetAlpha())+y*TMath::Cos(t.GetAlpha());
983 xyz1[2] = z;
984 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
985 rho = param[0];
986 radLength = param[1]; // get mean propagation parameters
8979685e 987 //
7b580082 988 // Find clusters
8979685e 989 //
990 sector = t.GetSector();
7b580082 991 Float_t ncl = FindClusters(sector,row0,rowlast,&t,clusters,tracklet);
992 if (tracklet.GetN()<GetTimeBinsPerPlane()/3) continue;
8979685e 993 //
7b580082 994 // Propagate and update track
8979685e 995 //
7b580082 996 for (Int_t itime= GetTimeBinsPerPlane()-1;itime>=0;itime--) {
997 Int_t ilayer = GetGlobalTimeBin(0, iplane,itime);
8979685e 998 expectedNumberOfClusters++;
999 t.fNExpected++;
1000 if (t.fX>345) t.fNExpectedLast++;
1001 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(ilayer));
1002 AliTRDcluster *cl=0;
1003 UInt_t index=0;
1004 Double_t maxChi2=fgkMaxChi2;
8979685e 1005 x = timeBin.GetX();
59393e34 1006 //
8979685e 1007 if (timeBin) {
1008 if (clusters[ilayer]>0) {
1009 index = clusters[ilayer];
1010 cl = (AliTRDcluster*)GetCluster(index);
1011 Double_t h01 = GetTiltFactor(cl);
1012 maxChi2=t.GetPredictedChi2(cl,h01);
1013 }
1014
1015 if (cl) {
1016 // if (cl->GetNPads()<5)
59393e34 1017 Double_t dxsample = timeBin.GetdX();
1018 t.SetSampledEdx(TMath::Abs(cl->GetQ()/dxsample));
8979685e 1019 Double_t h01 = GetTiltFactor(cl);
1020 Int_t det = cl->GetDetector();
1021 Int_t plane = fGeom->GetPlane(det);
1022 if (t.fX>345){
1023 t.fNLast++;
1024 t.fChi2Last+=maxChi2;
1025 }
59393e34 1026 Double_t xcluster = cl->GetX();
1027 t.PropagateTo(xcluster,radLength,rho);
8979685e 1028 if(!t.UpdateMI(cl,maxChi2,index,h01,plane)) {
1029 if(!t.Update(cl,maxChi2,index,h01)) {
8979685e 1030 }
1031 }
7b580082 1032 //
8979685e 1033 // reset material budget if 2 consecutive gold
1034 if (plane>0)
1035 if (t.fTracklets[plane].GetN()+t.fTracklets[plane-1].GetN()>20){
1036 t.fBudget[2] = 0;
1037 }
1038 }
1039 }
59393e34 1040 }
1041 ratio0 = ncl/Float_t(fTimeBinsPerPlane);
1042 Float_t ratio1 = Float_t(t.fN+1)/Float_t(t.fNExpected+1.);
1043 if (tracklet.GetChi2()<18.&&ratio0>0.8 && ratio1>0.6 && ratio0+ratio1>1.5 && t.GetNCross()==0 && TMath::Abs(t.GetSnp())<0.85&&t.fN>20){
1044 t.MakeBackupTrack(); // make backup of the track until is gold
1045 }
7b580082 1046
8979685e 1047 }
5443e65e 1048
75bd7f81 1049 return expectedNumberOfClusters;
1e9bb598 1050
75bd7f81 1051}
1e9bb598 1052
75bd7f81 1053//_____________________________________________________________________________
59393e34 1054Int_t AliTRDtracker::PropagateToX(AliTRDtrack& t, Double_t xToGo, Double_t maxStep)
5443e65e 1055{
75bd7f81 1056 //
5443e65e 1057 // Starting from current radial position of track <t> this function
1058 // extrapolates the track up to radial position <xToGo>.
1059 // Returns 1 if track reaches the plane, and 0 otherwise
75bd7f81 1060 //
1061
59393e34 1062 const Double_t kEpsilon = 0.00001;
1063 // Double_t tanmax = TMath::Tan(0.5*AliTRDgeometry::GetAlpha());
1064 Double_t xpos = t.GetX();
1065 Double_t dir = (xpos<xToGo) ? 1.:-1.;
1066 //
1067 while ( (xToGo-xpos)*dir > kEpsilon){
1068 Double_t step = dir*TMath::Min(TMath::Abs(xToGo-xpos), maxStep);
1069 //
1070 Double_t xyz0[3],xyz1[3],param[7],x,y,z;
1071 t.GetGlobalXYZ(xyz0[0],xyz0[1],xyz0[2]); //starting global position
1072 x = xpos+step;
1073 //
1074 if (!t.GetProlongation(x,y,z)) return 0; // no prolongation
1075 //
1076 xyz1[0] = x*TMath::Cos(t.GetAlpha())-y*TMath::Sin(t.GetAlpha());
1077 xyz1[1] = +x*TMath::Sin(t.GetAlpha())+y*TMath::Cos(t.GetAlpha());
1078 xyz1[2] = z;
1079 //
1080 AliKalmanTrack::MeanMaterialBudget(xyz0,xyz1,param);
1081 if (!t.PropagateTo(x,param[1],param[0])) return 0;
1082 AdjustSector(&t);
1083 xpos = t.GetX();
5443e65e 1084 }
75bd7f81 1085
5443e65e 1086 return 1;
5443e65e 1087
59393e34 1088}
5443e65e 1089
c630aafd 1090//_____________________________________________________________________________
1091Int_t AliTRDtracker::LoadClusters(TTree *cTree)
1092{
75bd7f81 1093 //
c630aafd 1094 // Fills clusters into TRD tracking_sectors
1095 // Note that the numbering scheme for the TRD tracking_sectors
1096 // differs from that of TRD sectors
75bd7f81 1097 //
1098
4f1c04d3 1099 cout<<"\n Read Sectors clusters"<<endl;
c630aafd 1100 if (ReadClusters(fClusters,cTree)) {
1101 Error("LoadClusters","Problem with reading the clusters !");
1102 return 1;
1103 }
1104 Int_t ncl=fClusters->GetEntriesFast();
b7a0917f 1105 fNclusters=ncl;
c630aafd 1106 cout<<"\n LoadSectors: sorting "<<ncl<<" clusters"<<endl;
1107
1108 UInt_t index;
3c625a9b 1109 for (Int_t ichamber=0;ichamber<5;ichamber++)
1110 for (Int_t isector=0;isector<18;isector++){
1111 fHoles[ichamber][isector]=kTRUE;
1112 }
1113
1114
c630aafd 1115 while (ncl--) {
1116// printf("\r %d left ",ncl);
1117 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(ncl);
029cd327 1118 Int_t detector=c->GetDetector();
1119 Int_t localTimeBin=c->GetLocalTimeBin();
c630aafd 1120 Int_t sector=fGeom->GetSector(detector);
1121 Int_t plane=fGeom->GetPlane(detector);
3c625a9b 1122
029cd327 1123 Int_t trackingSector = CookSectorIndex(sector);
3c625a9b 1124 if (c->GetLabel(0)>0){
1125 Int_t chamber = fGeom->GetChamber(detector);
1126 fHoles[chamber][trackingSector]=kFALSE;
1127 }
c630aafd 1128
029cd327 1129 Int_t gtb = fTrSec[trackingSector]->CookTimeBinIndex(plane,localTimeBin);
c630aafd 1130 if(gtb < 0) continue;
029cd327 1131 Int_t layer = fTrSec[trackingSector]->GetLayerNumber(gtb);
c630aafd 1132
1133 index=ncl;
69b55c55 1134 //
1135 // apply pos correction
59393e34 1136 Transform(c);
029cd327 1137 fTrSec[trackingSector]->GetLayer(layer)->InsertCluster(c,index);
c630aafd 1138 }
75bd7f81 1139
c630aafd 1140 return 0;
75bd7f81 1141
c630aafd 1142}
1143
5443e65e 1144//_____________________________________________________________________________
b7a0917f 1145void AliTRDtracker::UnloadClusters()
5443e65e 1146{
1147 //
1148 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1149 //
1150
1151 Int_t i, nentr;
1152
1153 nentr = fClusters->GetEntriesFast();
1154 for (i = 0; i < nentr; i++) delete fClusters->RemoveAt(i);
b7a0917f 1155 fNclusters = 0;
5443e65e 1156
1157 nentr = fSeeds->GetEntriesFast();
1158 for (i = 0; i < nentr; i++) delete fSeeds->RemoveAt(i);
1159
1160 nentr = fTracks->GetEntriesFast();
1161 for (i = 0; i < nentr; i++) delete fTracks->RemoveAt(i);
1162
1163 Int_t nsec = AliTRDgeometry::kNsect;
1164
1165 for (i = 0; i < nsec; i++) {
1166 for(Int_t pl = 0; pl < fTrSec[i]->GetNumberOfLayers(); pl++) {
1167 fTrSec[i]->GetLayer(pl)->Clear();
1168 }
1169 }
1170
1171}
1172
75bd7f81 1173//_____________________________________________________________________________
69b55c55 1174void AliTRDtracker::MakeSeedsMI(Int_t /*inner*/, Int_t /*outer*/, AliESD * esd)
7ad19338 1175{
1176 //
1177 // Creates seeds using clusters between position inner plane and outer plane
1178 //
75bd7f81 1179
c6f438c0 1180 const Double_t kMaxTheta = 1;
1181 const Double_t kMaxPhi = 2.0;
69b55c55 1182 //
1183 const Double_t kRoad0y = 6; // road for middle cluster
1184 const Double_t kRoad0z = 8.5; // road for middle cluster
1185 //
1186 const Double_t kRoad1y = 2; // road in y for seeded cluster
1187 const Double_t kRoad1z = 20; // road in z for seeded cluster
1188 //
1189 const Double_t kRoad2y = 3; // road in y for extrapolated cluster
1190 const Double_t kRoad2z = 20; // road in z for extrapolated cluster
c6f438c0 1191 const Int_t kMaxSeed = 3000;
69b55c55 1192 Int_t maxSec=AliTRDgeometry::kNsect;
7ad19338 1193
69b55c55 1194 //
1195 // linear fitters in planes
1196 TLinearFitter fitterTC(2,"hyp2"); // fitting with tilting pads - kz fixed - kz= Z/x, + vertex const
1197 TLinearFitter fitterT2(4,"hyp4"); // fitting with tilting pads - kz not fixed
1198 fitterTC.StoreData(kTRUE);
1199 fitterT2.StoreData(kTRUE);
1200 AliRieman rieman(1000); // rieman fitter
1201 AliRieman rieman2(1000); // rieman fitter
7ad19338 1202 //
1203 // find the maximal and minimal layer for the planes
7ad19338 1204 //
1205 Int_t layers[6][2];
69b55c55 1206 AliTRDpropagationLayer* reflayers[6];
7ad19338 1207 for (Int_t i=0;i<6;i++){layers[i][0]=10000; layers[i][1]=0;}
7ad19338 1208 for (Int_t ns=0;ns<maxSec;ns++){
1209 for (Int_t ilayer=0;ilayer<fTrSec[ns]->GetNumberOfLayers();ilayer++){
1210 AliTRDpropagationLayer& layer=*(fTrSec[ns]->GetLayer(ilayer));
1211 if (layer==0) continue;
1212 Int_t det = layer[0]->GetDetector();
1213 Int_t plane = fGeom->GetPlane(det);
1214 if (ilayer<layers[plane][0]) layers[plane][0] = ilayer;
1215 if (ilayer>layers[plane][1]) layers[plane][1] = ilayer;
1216 }
1217 }
1218 //
3551db50 1219 AliTRDpadPlane *padPlane = AliTRDCommonParam::Instance()->GetPadPlane(0,0);
69b55c55 1220 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
1221 Double_t hL[6]; // tilting angle
1222 Double_t xcl[6]; // x - position of reference cluster
1223 Double_t ycl[6]; // y - position of reference cluster
1224 Double_t zcl[6]; // z - position of reference cluster
1225 AliTRDcluster *cl[6]={0,0,0,0,0,0}; // seeding clusters
1226 Float_t padlength[6]={10,10,10,10,10,10}; //current pad-length
1227 Double_t chi2R =0, chi2Z=0;
1228 Double_t chi2RF =0, chi2ZF=0;
1229 //
1230 Int_t nclusters; // total number of clusters
1231 for (Int_t i=0;i<6;i++) {hL[i]=h01; if (i%2==1) hL[i]*=-1.;}
1232 //
1233 //
1234 // registered seed
c6f438c0 1235 AliTRDseed *pseed = new AliTRDseed[kMaxSeed*6];
1236 AliTRDseed *seed[kMaxSeed];
1237 for (Int_t iseed=0;iseed<kMaxSeed;iseed++) seed[iseed]= &pseed[iseed*6];
69b55c55 1238 AliTRDseed *cseed = seed[0];
1239 //
c6f438c0 1240 Double_t seedquality[kMaxSeed];
1241 Double_t seedquality2[kMaxSeed];
1242 Double_t seedparams[kMaxSeed][7];
1243 Int_t seedlayer[kMaxSeed];
69b55c55 1244 Int_t registered =0;
c6f438c0 1245 Int_t sort[kMaxSeed];
69b55c55 1246 //
1247 // seeding part
1248 //
1249 for (Int_t ns = 0; ns<maxSec; ns++){ //loop over sectors
1250 //for (Int_t ns = 0; ns<5; ns++){ //loop over sectors
1251 registered = 0; // reset registerd seed counter
1252 cseed = seed[registered];
1253 Float_t iter=0;
1254 for (Int_t sLayer=2; sLayer>=0;sLayer--){
1255 //for (Int_t dseed=5;dseed<15; dseed+=3){ //loop over central seeding time bins
1256 iter+=1.;
1257 Int_t dseed = 5+Int_t(iter)*3;
1258 // Initialize seeding layers
1259 for (Int_t ilayer=0;ilayer<6;ilayer++){
1260 reflayers[ilayer] = fTrSec[ns]->GetLayer(layers[ilayer][1]-dseed);
1261 xcl[ilayer] = reflayers[ilayer]->GetX();
1262 }
7ad19338 1263 //
69b55c55 1264 Double_t xref = (xcl[sLayer+1] + xcl[sLayer+2])*0.5;
1265 AliTRDpropagationLayer& layer0=*reflayers[sLayer+0];
1266 AliTRDpropagationLayer& layer1=*reflayers[sLayer+1];
1267 AliTRDpropagationLayer& layer2=*reflayers[sLayer+2];
1268 AliTRDpropagationLayer& layer3=*reflayers[sLayer+3];
1269 //
1270 Int_t maxn3 = layer3;
1271 for (Int_t icl3=0;icl3<maxn3;icl3++){
1272 AliTRDcluster *cl3 = layer3[icl3];
1273 if (!cl3) continue;
1274 padlength[sLayer+3] = TMath::Sqrt(cl3->GetSigmaZ2()*12.);
1275 ycl[sLayer+3] = cl3->GetY();
1276 zcl[sLayer+3] = cl3->GetZ();
c6f438c0 1277 Float_t yymin0 = ycl[sLayer+3] - 1- kMaxPhi *(xcl[sLayer+3]-xcl[sLayer+0]);
1278 Float_t yymax0 = ycl[sLayer+3] + 1+ kMaxPhi *(xcl[sLayer+3]-xcl[sLayer+0]);
69b55c55 1279 Int_t maxn0 = layer0; //
1280 for (Int_t icl0=layer0.Find(yymin0);icl0<maxn0;icl0++){
1281 AliTRDcluster *cl0 = layer0[icl0];
1282 if (!cl0) continue;
1283 if (cl3->IsUsed()&&cl0->IsUsed()) continue;
1284 ycl[sLayer+0] = cl0->GetY();
1285 zcl[sLayer+0] = cl0->GetZ();
1286 if ( ycl[sLayer+0]>yymax0) break;
1287 Double_t tanphi = (ycl[sLayer+3]-ycl[sLayer+0])/(xcl[sLayer+3]-xcl[sLayer+0]);
c6f438c0 1288 if (TMath::Abs(tanphi)>kMaxPhi) continue;
69b55c55 1289 Double_t tantheta = (zcl[sLayer+3]-zcl[sLayer+0])/(xcl[sLayer+3]-xcl[sLayer+0]);
c6f438c0 1290 if (TMath::Abs(tantheta)>kMaxTheta) continue;
69b55c55 1291 padlength[sLayer+0] = TMath::Sqrt(cl0->GetSigmaZ2()*12.);
1292 //
1293 // expected position in 1 layer
1294 Double_t y1exp = ycl[sLayer+0]+(tanphi) *(xcl[sLayer+1]-xcl[sLayer+0]);
1295 Double_t z1exp = zcl[sLayer+0]+(tantheta)*(xcl[sLayer+1]-xcl[sLayer+0]);
1296 Float_t yymin1 = y1exp - kRoad0y-tanphi;
1297 Float_t yymax1 = y1exp + kRoad0y+tanphi;
1298 Int_t maxn1 = layer1; //
1299 //
1300 for (Int_t icl1=layer1.Find(yymin1);icl1<maxn1;icl1++){
1301 AliTRDcluster *cl1 = layer1[icl1];
1302 if (!cl1) continue;
1303 Int_t nusedCl = 0;
1304 if (cl3->IsUsed()) nusedCl++;
1305 if (cl0->IsUsed()) nusedCl++;
1306 if (cl1->IsUsed()) nusedCl++;
1307 if (nusedCl>1) continue;
1308 ycl[sLayer+1] = cl1->GetY();
1309 zcl[sLayer+1] = cl1->GetZ();
1310 if ( ycl[sLayer+1]>yymax1) break;
1311 if (TMath::Abs(ycl[sLayer+1]-y1exp)>kRoad0y+tanphi) continue;
1312 if (TMath::Abs(zcl[sLayer+1]-z1exp)>kRoad0z) continue;
1313 padlength[sLayer+1] = TMath::Sqrt(cl1->GetSigmaZ2()*12.);
1314 //
1315 Double_t y2exp = ycl[sLayer+0]+(tanphi) *(xcl[sLayer+2]-xcl[sLayer+0])+(ycl[sLayer+1]-y1exp);
1316 Double_t z2exp = zcl[sLayer+0]+(tantheta)*(xcl[sLayer+2]-xcl[sLayer+0]);
1317 Int_t index2 = layer2.FindNearestCluster(y2exp,z2exp,kRoad1y, kRoad1z);
1318 if (index2<=0) continue;
1319 AliTRDcluster *cl2 = (AliTRDcluster*)GetCluster(index2);
1320 padlength[sLayer+2] = TMath::Sqrt(cl2->GetSigmaZ2()*12.);
1321 ycl[sLayer+2] = cl2->GetY();
1322 zcl[sLayer+2] = cl2->GetZ();
1323 if (TMath::Abs(cl2->GetZ()-z2exp)>kRoad0z) continue;
1324 //
1325 rieman.Reset();
1326 rieman.AddPoint(xcl[sLayer+0],ycl[sLayer+0],zcl[sLayer+0],1,10);
1327 rieman.AddPoint(xcl[sLayer+1],ycl[sLayer+1],zcl[sLayer+1],1,10);
1328 rieman.AddPoint(xcl[sLayer+3],ycl[sLayer+3],zcl[sLayer+3],1,10);
1329 rieman.AddPoint(xcl[sLayer+2],ycl[sLayer+2],zcl[sLayer+2],1,10);
1330 rieman.Update();
1331 //
1332 // reset fitter
1333 for (Int_t iLayer=0;iLayer<6;iLayer++){
1334 cseed[iLayer].Reset();
1335 }
1336 chi2Z =0.; chi2R=0.;
1337 for (Int_t iLayer=0;iLayer<4;iLayer++){
1338 cseed[sLayer+iLayer].fZref[0] = rieman.GetZat(xcl[sLayer+iLayer]);
1339 chi2Z += (cseed[sLayer+iLayer].fZref[0]- zcl[sLayer+iLayer])*
1340 (cseed[sLayer+iLayer].fZref[0]- zcl[sLayer+iLayer]);
1341 cseed[sLayer+iLayer].fZref[1] = rieman.GetDZat(xcl[sLayer+iLayer]);
1342 cseed[sLayer+iLayer].fYref[0] = rieman.GetYat(xcl[sLayer+iLayer]);
1343 chi2R += (cseed[sLayer+iLayer].fYref[0]- ycl[sLayer+iLayer])*
1344 (cseed[sLayer+iLayer].fYref[0]- ycl[sLayer+iLayer]);
1345 cseed[sLayer+iLayer].fYref[1] = rieman.GetDYat(xcl[sLayer+iLayer]);
1346 }
1347 if (TMath::Sqrt(chi2R)>1./iter) continue;
1348 if (TMath::Sqrt(chi2Z)>7./iter) continue;
1349 //
1350 //
1351 //
1352 Float_t minmax[2]={-100,100};
1353 for (Int_t iLayer=0;iLayer<4;iLayer++){
1354 Float_t max = zcl[sLayer+iLayer]+padlength[sLayer+iLayer]*0.5+1 -cseed[sLayer+iLayer].fZref[0];
1355 if (max<minmax[1]) minmax[1]=max;
1356 Float_t min = zcl[sLayer+iLayer]-padlength[sLayer+iLayer]*0.5-1 -cseed[sLayer+iLayer].fZref[0];
1357 if (min>minmax[0]) minmax[0]=min;
1358 }
1359 Bool_t isFake = kFALSE;
1360 if (cl0->GetLabel(0)!=cl3->GetLabel(0)) isFake = kTRUE;
1361 if (cl1->GetLabel(0)!=cl3->GetLabel(0)) isFake = kTRUE;
1362 if (cl2->GetLabel(0)!=cl3->GetLabel(0)) isFake = kTRUE;
d337ef8d 1363 if (AliTRDReconstructor::StreamLevel()>0){
1364 if ((!isFake) || (icl3%10)==0 ){ //debugging print
1365 TTreeSRedirector& cstream = *fDebugStreamer;
1366 cstream<<"Seeds0"<<
1367 "isFake="<<isFake<<
1368 "Cl0.="<<cl0<<
1369 "Cl1.="<<cl1<<
1370 "Cl2.="<<cl2<<
1371 "Cl3.="<<cl3<<
1372 "Xref="<<xref<<
1373 "X0="<<xcl[sLayer+0]<<
1374 "X1="<<xcl[sLayer+1]<<
1375 "X2="<<xcl[sLayer+2]<<
1376 "X3="<<xcl[sLayer+3]<<
1377 "Y2exp="<<y2exp<<
1378 "Z2exp="<<z2exp<<
1379 "Chi2R="<<chi2R<<
1380 "Chi2Z="<<chi2Z<<
1381 "Seed0.="<<&cseed[sLayer+0]<<
1382 "Seed1.="<<&cseed[sLayer+1]<<
1383 "Seed2.="<<&cseed[sLayer+2]<<
1384 "Seed3.="<<&cseed[sLayer+3]<<
1385 "Zmin="<<minmax[0]<<
1386 "Zmax="<<minmax[1]<<
1387 "\n";
1388 }
69b55c55 1389 }
1390
1391 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1392 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1393 //<<<<<<<<<<<<<<<<<< FIT SEEDING PART <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1394 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1395 cl[sLayer+0] = cl0;
1396 cl[sLayer+1] = cl1;
1397 cl[sLayer+2] = cl2;
1398 cl[sLayer+3] = cl3;
1399 Bool_t isOK=kTRUE;
1400 for (Int_t jLayer=0;jLayer<4;jLayer++){
1401 cseed[sLayer+jLayer].fTilt = hL[sLayer+jLayer];
1402 cseed[sLayer+jLayer].fPadLength = padlength[sLayer+jLayer];
1403 cseed[sLayer+jLayer].fX0 = xcl[sLayer+jLayer];
1404 for (Int_t iter=0; iter<2; iter++){
1405 //
1406 // in iteration 0 we try only one pad-row
1407 // if quality not sufficient we try 2 pad-rows - about 5% of tracks cross 2 pad-rows
1408 //
1409 AliTRDseed tseed = cseed[sLayer+jLayer];
1410 Float_t roadz = padlength[sLayer+jLayer]*0.5;
1411 if (iter>0) roadz = padlength[sLayer+jLayer];
1412 //
1413 Float_t quality =10000;
1414 for (Int_t iTime=2;iTime<20;iTime++){
1415 AliTRDpropagationLayer& layer = *(fTrSec[ns]->GetLayer(layers[sLayer+jLayer][1]-iTime));
1416 Double_t dxlayer= layer.GetX()-xcl[sLayer+jLayer];
1417 Double_t zexp = cl[sLayer+jLayer]->GetZ() ;
1418 if (iter>0){
1419 // try 2 pad-rows in second iteration
1420 zexp = tseed.fZref[0]+ tseed.fZref[1]*dxlayer;
1421 if (zexp>cl[sLayer+jLayer]->GetZ()) zexp = cl[sLayer+jLayer]->GetZ()+padlength[sLayer+jLayer]*0.5;
1422 if (zexp<cl[sLayer+jLayer]->GetZ()) zexp = cl[sLayer+jLayer]->GetZ()-padlength[sLayer+jLayer]*0.5;
1423 }
1424 //
1425 Double_t yexp = tseed.fYref[0]+
1426 tseed.fYref[1]*dxlayer;
1427 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y, roadz);
1428 if (index<=0) continue;
1429 AliTRDcluster *cl = (AliTRDcluster*)GetCluster(index);
1430 //
1431 tseed.fIndexes[iTime] = index;
1432 tseed.fClusters[iTime] = cl; // register cluster
1433 tseed.fX[iTime] = dxlayer; // register cluster
1434 tseed.fY[iTime] = cl->GetY(); // register cluster
1435 tseed.fZ[iTime] = cl->GetZ(); // register cluster
1436 }
1437 tseed.Update();
1438 //count the number of clusters and distortions into quality
1439 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1440 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1441 TMath::Abs(tseed.fYfit[0]-tseed.fYref[0])/0.2+
1442 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
c6f438c0 1443 if (iter==0 && tseed.IsOK()) {
69b55c55 1444 cseed[sLayer+jLayer] = tseed;
1445 quality = tquality;
1446 if (tquality<5) break;
1447 }
c6f438c0 1448 if (tseed.IsOK() && tquality<quality)
69b55c55 1449 cseed[sLayer+jLayer] = tseed;
1450 }
c6f438c0 1451 if (!cseed[sLayer+jLayer].IsOK()){
69b55c55 1452 isOK = kFALSE;
1453 break;
1454 }
1455 cseed[sLayer+jLayer].CookLabels();
1456 cseed[sLayer+jLayer].UpdateUsed();
1457 nusedCl+= cseed[sLayer+jLayer].fNUsed;
1458 if (nusedCl>25){
1459 isOK = kFALSE;
1460 break;
1461 }
1462 }
1463 //
1464 if (!isOK) continue;
1465 nclusters=0;
1466 for (Int_t iLayer=0;iLayer<4;iLayer++){
c6f438c0 1467 if (cseed[sLayer+iLayer].IsOK()){
69b55c55 1468 nclusters+=cseed[sLayer+iLayer].fN2;
1469 }
1470 }
1471 //
1472 // iteration 0
1473 rieman.Reset();
1474 for (Int_t iLayer=0;iLayer<4;iLayer++){
1475 rieman.AddPoint(xcl[sLayer+iLayer],cseed[sLayer+iLayer].fYfitR[0],
1476 cseed[sLayer+iLayer].fZProb,1,10);
1477 }
1478 rieman.Update();
1479 //
1480 //
1481 chi2R =0; chi2Z=0;
1482 for (Int_t iLayer=0;iLayer<4;iLayer++){
1483 cseed[sLayer+iLayer].fYref[0] = rieman.GetYat(xcl[sLayer+iLayer]);
1484 chi2R += (cseed[sLayer+iLayer].fYref[0]-cseed[sLayer+iLayer].fYfitR[0])*
1485 (cseed[sLayer+iLayer].fYref[0]-cseed[sLayer+iLayer].fYfitR[0]);
1486 cseed[sLayer+iLayer].fYref[1] = rieman.GetDYat(xcl[sLayer+iLayer]);
1487 cseed[sLayer+iLayer].fZref[0] = rieman.GetZat(xcl[sLayer+iLayer]);
1488 chi2Z += (cseed[sLayer+iLayer].fZref[0]- cseed[sLayer+iLayer].fMeanz)*
1489 (cseed[sLayer+iLayer].fZref[0]- cseed[sLayer+iLayer].fMeanz);
1490 cseed[sLayer+iLayer].fZref[1] = rieman.GetDZat(xcl[sLayer+iLayer]);
1491 }
1492 Double_t curv = rieman.GetC();
1493 //
1494 // likelihoods
1495 //
1496 Double_t sumda =
1497 TMath::Abs(cseed[sLayer+0].fYfitR[1]- cseed[sLayer+0].fYref[1])+
1498 TMath::Abs(cseed[sLayer+1].fYfitR[1]- cseed[sLayer+1].fYref[1])+
1499 TMath::Abs(cseed[sLayer+2].fYfitR[1]- cseed[sLayer+2].fYref[1])+
1500 TMath::Abs(cseed[sLayer+3].fYfitR[1]- cseed[sLayer+3].fYref[1]);
1501 Double_t likea = TMath::Exp(-sumda*10.6);
1502 Double_t likechi2 = 0.0000000001;
1503 if (chi2R<0.5) likechi2+=TMath::Exp(-TMath::Sqrt(chi2R)*7.73);
1504 Double_t likechi2z = TMath::Exp(-chi2Z*0.088)/TMath::Exp(-chi2Z*0.019);
1505 Double_t likeN = TMath::Exp(-(72-nclusters)*0.19);
1506 Double_t like = likea*likechi2*likechi2z*likeN;
1507 //
1508 Double_t likePrimY = TMath::Exp(-TMath::Abs(cseed[sLayer+0].fYref[1]-130*curv)*1.9);
1509 Double_t likePrimZ = TMath::Exp(-TMath::Abs(cseed[sLayer+0].fZref[1]-
1510 cseed[sLayer+0].fZref[0]/xcl[sLayer+0])*5.9);
1511 Double_t likePrim = TMath::Max(likePrimY*likePrimZ,0.0005);
1512
1513 seedquality[registered] = like;
1514 seedlayer[registered] = sLayer;
1515 if (TMath::Log(0.000000000000001+like)<-15) continue;
1516 AliTRDseed seedb[6];
1517 for (Int_t iLayer=0;iLayer<6;iLayer++){
1518 seedb[iLayer] = cseed[iLayer];
1519 }
1520 //
1521 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1522 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1523 //<<<<<<<<<<<<<<< FULL TRACK FIT PART <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1524 //<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
1525 //
1526 Int_t nlayers = 0;
1527 Int_t nusedf = 0;
1528 Int_t findable = 0;
1529 //
1530 // add new layers - avoid long extrapolation
1531 //
1532 Int_t tLayer[2]={0,0};
1533 if (sLayer==2) {tLayer[0]=1; tLayer[1]=0;}
1534 if (sLayer==1) {tLayer[0]=5; tLayer[1]=0;}
1535 if (sLayer==0) {tLayer[0]=4; tLayer[1]=5;}
1536 //
1537 for (Int_t iLayer=0;iLayer<2;iLayer++){
1538 Int_t jLayer = tLayer[iLayer]; // set tracking layer
1539 cseed[jLayer].Reset();
1540 cseed[jLayer].fTilt = hL[jLayer];
1541 cseed[jLayer].fPadLength = padlength[jLayer];
1542 cseed[jLayer].fX0 = xcl[jLayer];
1543 // get pad length and rough cluster
1544 Int_t indexdummy = reflayers[jLayer]->FindNearestCluster(cseed[jLayer].fYref[0],
1545 cseed[jLayer].fZref[0],kRoad2y,kRoad2z);
1546 if (indexdummy<=0) continue;
1547 AliTRDcluster *cldummy = (AliTRDcluster*)GetCluster(indexdummy);
1548 padlength[jLayer] = TMath::Sqrt(cldummy->GetSigmaZ2()*12.);
1549 }
1550 AliTRDseed::FitRiemanTilt(cseed, kTRUE);
1551 //
1552 for (Int_t iLayer=0;iLayer<2;iLayer++){
1553 Int_t jLayer = tLayer[iLayer]; // set tracking layer
c6f438c0 1554 if ( (jLayer==0) && !(cseed[1].IsOK())) continue; // break not allowed
1555 if ( (jLayer==5) && !(cseed[4].IsOK())) continue; // break not allowed
69b55c55 1556 Float_t zexp = cseed[jLayer].fZref[0];
1557 Double_t zroad = padlength[jLayer]*0.5+1.;
1558 //
1559 //
1560 for (Int_t iter=0;iter<2;iter++){
1561 AliTRDseed tseed = cseed[jLayer];
1562 Float_t quality = 10000;
1563 for (Int_t iTime=2;iTime<20;iTime++){
1564 AliTRDpropagationLayer& layer = *(fTrSec[ns]->GetLayer(layers[jLayer][1]-iTime));
1565 Double_t dxlayer = layer.GetX()-xcl[jLayer];
1566 Double_t yexp = tseed.fYref[0]+tseed.fYref[1]*dxlayer;
1567 Float_t yroad = kRoad1y;
1568 Int_t index = layer.FindNearestCluster(yexp,zexp, yroad, zroad);
1569 if (index<=0) continue;
1570 AliTRDcluster *cl = (AliTRDcluster*)GetCluster(index);
1571 //
1572 tseed.fIndexes[iTime] = index;
1573 tseed.fClusters[iTime] = cl; // register cluster
1574 tseed.fX[iTime] = dxlayer; // register cluster
1575 tseed.fY[iTime] = cl->GetY(); // register cluster
1576 tseed.fZ[iTime] = cl->GetZ(); // register cluster
1577 }
1578 tseed.Update();
c6f438c0 1579 if (tseed.IsOK()){
69b55c55 1580 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1581 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1582 TMath::Abs(tseed.fYfit[0]-tseed.fYref[0])/0.2+
1583 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1584 //
1585 if (tquality<quality){
1586 cseed[jLayer]=tseed;
1587 quality = tquality;
1588 }
1589 }
1590 zroad*=2.;
1591 }
c6f438c0 1592 if ( cseed[jLayer].IsOK()){
69b55c55 1593 cseed[jLayer].CookLabels();
1594 cseed[jLayer].UpdateUsed();
1595 nusedf+= cseed[jLayer].fNUsed;
1596 AliTRDseed::FitRiemanTilt(cseed, kTRUE);
1597 }
1598 }
1599 //
1600 //
1601 // make copy
1602 AliTRDseed bseed[6];
1603 for (Int_t jLayer=0;jLayer<6;jLayer++){
1604 bseed[jLayer] = cseed[jLayer];
1605 }
1606 Float_t lastquality = 10000;
1607 Float_t lastchi2 = 10000;
1608 Float_t chi2 = 1000;
1609
1610 //
1611 for (Int_t iter =0; iter<4;iter++){
1612 //
1613 // sort tracklets according "quality", try to "improve" 4 worst
1614 //
1615 Float_t sumquality = 0;
1616 Float_t squality[6];
1617 Int_t sortindexes[6];
1618 for (Int_t jLayer=0;jLayer<6;jLayer++){
c6f438c0 1619 if (bseed[jLayer].IsOK()){
69b55c55 1620 AliTRDseed &tseed = bseed[jLayer];
1621 Double_t zcor = tseed.fTilt*(tseed.fZProb-tseed.fZref[0]);
1622 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1623 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1624 TMath::Abs(tseed.fYfit[0]-(tseed.fYref[0]-zcor))/0.2+
1625 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1626 squality[jLayer] = tquality;
1627 }
1628 else squality[jLayer]=-1;
1629 sumquality +=squality[jLayer];
1630 }
1631
1632 if (sumquality>=lastquality || chi2>lastchi2) break;
1633 lastquality = sumquality;
1634 lastchi2 = chi2;
1635 if (iter>0){
1636 for (Int_t jLayer=0;jLayer<6;jLayer++){
1637 cseed[jLayer] = bseed[jLayer];
1638 }
1639 }
1640 TMath::Sort(6,squality,sortindexes,kFALSE);
1641 //
1642 //
1643 for (Int_t jLayer=5;jLayer>1;jLayer--){
1644 Int_t bLayer = sortindexes[jLayer];
1645 AliTRDseed tseed = bseed[bLayer];
1646 for (Int_t iTime=2;iTime<20;iTime++){
1647 AliTRDpropagationLayer& layer = *(fTrSec[ns]->GetLayer(layers[bLayer][1]-iTime));
1648 Double_t dxlayer= layer.GetX()-xcl[bLayer];
1649 //
1650 Double_t zexp = tseed.fZref[0];
1651 Double_t zcor = tseed.fTilt*(tseed.fZProb-tseed.fZref[0]);
1652 //
1653 Float_t roadz = padlength[bLayer]+1;
1654 if (TMath::Abs(tseed.fZProb-zexp)> padlength[bLayer]*0.5) {roadz = padlength[bLayer]*0.5;}
1655 if (tseed.fZfit[1]*tseed.fZref[1]<0) {roadz = padlength[bLayer]*0.5;}
1656 if (TMath::Abs(tseed.fZProb-zexp)<0.1*padlength[bLayer]) {
1657 zexp = tseed.fZProb;
1658 roadz = padlength[bLayer]*0.5;
1659 }
1660 //
1661 Double_t yexp = tseed.fYref[0]+
1662 tseed.fYref[1]*dxlayer-zcor;
1663 Int_t index = layer.FindNearestCluster(yexp,zexp,kRoad1y, roadz);
1664 if (index<=0) continue;
1665 AliTRDcluster *cl = (AliTRDcluster*)GetCluster(index);
1666 //
1667 tseed.fIndexes[iTime] = index;
1668 tseed.fClusters[iTime] = cl; // register cluster
1669 tseed.fX[iTime] = dxlayer; // register cluster
1670 tseed.fY[iTime] = cl->GetY(); // register cluster
1671 tseed.fZ[iTime] = cl->GetZ(); // register cluster
1672 }
1673 tseed.Update();
c6f438c0 1674 if (tseed.IsOK()) {
69b55c55 1675 Float_t dangle = tseed.fYfit[1]-tseed.fYref[1];
1676 Double_t zcor = tseed.fTilt*(tseed.fZProb-tseed.fZref[0]);
1677 //
1678 Float_t tquality = (18-tseed.fN2)/2. + TMath::Abs(dangle)/0.1+
1679 TMath::Abs(tseed.fYfit[0]-(tseed.fYref[0]-zcor))/0.2+
1680 2.*TMath::Abs(tseed.fMeanz-tseed.fZref[0])/padlength[jLayer];
1681 //
1682 if (tquality<squality[bLayer])
1683 bseed[bLayer] = tseed;
1684 }
1685 }
1686 chi2 = AliTRDseed::FitRiemanTilt(bseed, kTRUE);
1687 }
1688 //
1689 //
1690 //
1691 nclusters = 0;
1692 nlayers = 0;
1693 findable = 0;
1694 for (Int_t iLayer=0;iLayer<6;iLayer++) {
1695 if (TMath::Abs(cseed[iLayer].fYref[0]/cseed[iLayer].fX0)<0.15)
1696 findable++;
c6f438c0 1697 if (cseed[iLayer].IsOK()){
69b55c55 1698 nclusters+=cseed[iLayer].fN2;
1699 nlayers++;
1700 }
1701 }
1702 if (nlayers<3) continue;
1703 rieman.Reset();
1704 for (Int_t iLayer=0;iLayer<6;iLayer++){
c6f438c0 1705 if (cseed[iLayer].IsOK()) rieman.AddPoint(xcl[iLayer],cseed[iLayer].fYfitR[0],
69b55c55 1706 cseed[iLayer].fZProb,1,10);
1707 }
1708 rieman.Update();
1709 //
1710 chi2RF =0;
1711 chi2ZF =0;
1712 for (Int_t iLayer=0;iLayer<6;iLayer++){
c6f438c0 1713 if (cseed[iLayer].IsOK()){
69b55c55 1714 cseed[iLayer].fYref[0] = rieman.GetYat(xcl[iLayer]);
1715 chi2RF += (cseed[iLayer].fYref[0]-cseed[iLayer].fYfitR[0])*
1716 (cseed[iLayer].fYref[0]-cseed[iLayer].fYfitR[0]);
1717 cseed[iLayer].fYref[1] = rieman.GetDYat(xcl[iLayer]);
1718 cseed[iLayer].fZref[0] = rieman.GetZat(xcl[iLayer]);
1719 chi2ZF += (cseed[iLayer].fZref[0]- cseed[iLayer].fMeanz)*
1720 (cseed[iLayer].fZref[0]- cseed[iLayer].fMeanz);
1721 cseed[iLayer].fZref[1] = rieman.GetDZat(xcl[iLayer]);
1722 }
1723 }
1724 chi2RF/=TMath::Max((nlayers-3.),1.);
1725 chi2ZF/=TMath::Max((nlayers-3.),1.);
1726 curv = rieman.GetC();
1727
1728 //
1729
1730 Double_t xref2 = (xcl[2]+xcl[3])*0.5; // middle of the chamber
1731 Double_t dzmf = rieman.GetDZat(xref2);
1732 Double_t zmf = rieman.GetZat(xref2);
1733 //
1734 // fit hyperplane
1735 //
1736 Int_t npointsT =0;
1737 fitterTC.ClearPoints();
1738 fitterT2.ClearPoints();
1739 rieman2.Reset();
1740 for (Int_t iLayer=0; iLayer<6;iLayer++){
c6f438c0 1741 if (!cseed[iLayer].IsOK()) continue;
69b55c55 1742 for (Int_t itime=0;itime<25;itime++){
1743 if (!cseed[iLayer].fUsable[itime]) continue;
1744 Double_t x = cseed[iLayer].fX[itime]+cseed[iLayer].fX0-xref2; // x relative to the midle chamber
1745 Double_t y = cseed[iLayer].fY[itime];
1746 Double_t z = cseed[iLayer].fZ[itime];
1747 // ExB correction to the correction
1748 // tilted rieman
1749 //
1750 Double_t uvt[6];
1751 Double_t x2 = cseed[iLayer].fX[itime]+cseed[iLayer].fX0; // global x
1752 //
1753 Double_t t = 1./(x2*x2+y*y);
1754 uvt[1] = t; // t
1755 uvt[0] = 2.*x2*uvt[1]; // u
1756 //
1757 uvt[2] = 2.0*hL[iLayer]*uvt[1];
1758 uvt[3] = 2.0*hL[iLayer]*x*uvt[1];
1759 uvt[4] = 2.0*(y+hL[iLayer]*z)*uvt[1];
1760 //
1761 Double_t error = 2*0.2*uvt[1];
1762 fitterT2.AddPoint(uvt,uvt[4],error);
1763 //
1764 // constrained rieman
1765 //
1766 z =cseed[iLayer].fZ[itime];
1767 uvt[0] = 2.*x2*t; // u
1768 uvt[1] = 2*hL[iLayer]*x2*uvt[1];
1769 uvt[2] = 2*(y+hL[iLayer]*(z-GetZ()))*t;
1770 fitterTC.AddPoint(uvt,uvt[2],error);
1771 //
1772 rieman2.AddPoint(x2,y,z,1,10);
1773 npointsT++;
1774 }
1775 }
1776 rieman2.Update();
1777 fitterTC.Eval();
1778 fitterT2.Eval();
1779 Double_t rpolz0 = fitterT2.GetParameter(3);
1780 Double_t rpolz1 = fitterT2.GetParameter(4);
1781 //
1782 // linear fitter - not possible to make boundaries
1783 // non accept non possible z and dzdx combination
1784 //
1785 Bool_t acceptablez =kTRUE;
1786 for (Int_t iLayer=0; iLayer<6;iLayer++){
c6f438c0 1787 if (cseed[iLayer].IsOK()){
69b55c55 1788 Double_t zT2 = rpolz0+rpolz1*(xcl[iLayer] - xref2);
1789 if (TMath::Abs(cseed[iLayer].fZProb-zT2)>padlength[iLayer]*0.5+1)
1790 acceptablez = kFALSE;
1791 }
1792 }
1793 if (!acceptablez){
1794 fitterT2.FixParameter(3,zmf);
1795 fitterT2.FixParameter(4,dzmf);
1796 fitterT2.Eval();
1797 fitterT2.ReleaseParameter(3);
1798 fitterT2.ReleaseParameter(4);
1799 rpolz0 = fitterT2.GetParameter(3);
1800 rpolz1 = fitterT2.GetParameter(4);
1801 }
1802 //
1803 Double_t chi2TR = fitterT2.GetChisquare()/Float_t(npointsT);
1804 Double_t chi2TC = fitterTC.GetChisquare()/Float_t(npointsT);
1805 //
1806 Double_t polz1c = fitterTC.GetParameter(2);
1807 Double_t polz0c = polz1c*xref2;
1808 //
1809 Double_t aC = fitterTC.GetParameter(0);
1810 Double_t bC = fitterTC.GetParameter(1);
c6f438c0 1811 Double_t cC = aC/TMath::Sqrt(bC*bC+1.); // curvature
69b55c55 1812 //
1813 Double_t aR = fitterT2.GetParameter(0);
1814 Double_t bR = fitterT2.GetParameter(1);
1815 Double_t dR = fitterT2.GetParameter(2);
c6f438c0 1816 Double_t cR = 1+bR*bR-dR*aR;
69b55c55 1817 Double_t dca = 0.;
c6f438c0 1818 if (cR>0){
69b55c55 1819 dca = -dR/(TMath::Sqrt(1+bR*bR-dR*aR)+TMath::Sqrt(1+bR*bR));
c6f438c0 1820 cR = aR/TMath::Sqrt(cR);
69b55c55 1821 }
1822 //
1823 Double_t chi2ZT2=0, chi2ZTC=0;
1824 for (Int_t iLayer=0; iLayer<6;iLayer++){
c6f438c0 1825 if (cseed[iLayer].IsOK()){
69b55c55 1826 Double_t zT2 = rpolz0+rpolz1*(xcl[iLayer] - xref2);
1827 Double_t zTC = polz0c+polz1c*(xcl[iLayer] - xref2);
1828 chi2ZT2 += TMath::Abs(cseed[iLayer].fMeanz-zT2);
1829 chi2ZTC += TMath::Abs(cseed[iLayer].fMeanz-zTC);
1830 }
1831 }
1832 chi2ZT2/=TMath::Max((nlayers-3.),1.);
1833 chi2ZTC/=TMath::Max((nlayers-3.),1.);
1834 //
1835 //
1836 //
1837 AliTRDseed::FitRiemanTilt(cseed, kTRUE);
1838 Float_t sumdaf = 0;
1839 for (Int_t iLayer=0;iLayer<6;iLayer++){
c6f438c0 1840 if (cseed[iLayer].IsOK())
69b55c55 1841 sumdaf += TMath::Abs((cseed[iLayer].fYfit[1]-cseed[iLayer].fYref[1])/cseed[iLayer].fSigmaY2);
1842 }
1843 sumdaf /= Float_t (nlayers-2.);
1844 //
1845 // likelihoods for full track
1846 //
1847 Double_t likezf = TMath::Exp(-chi2ZF*0.14);
1848 Double_t likechi2C = TMath::Exp(-chi2TC*0.677);
1849 Double_t likechi2TR = TMath::Exp(-chi2TR*0.78);
1850 Double_t likeaf = TMath::Exp(-sumdaf*3.23);
1851 seedquality2[registered] = likezf*likechi2TR*likeaf;
1852// Bool_t isGold = kFALSE;
1853//
1854// if (nlayers == 6 && TMath::Log(0.000000001+seedquality2[index])<-5.) isGold =kTRUE; // gold
1855// if (nlayers == findable && TMath::Log(0.000000001+seedquality2[index])<-4.) isGold =kTRUE; // gold
1856// if (isGold &&nusedf<10){
1857// for (Int_t jLayer=0;jLayer<6;jLayer++){
c6f438c0 1858// if ( seed[index][jLayer].IsOK()&&TMath::Abs(seed[index][jLayer].fYfit[1]-seed[index][jLayer].fYfit[1])<0.1)
69b55c55 1859// seed[index][jLayer].UseClusters(); //sign gold
1860// }
1861// }
1862 //
1863 //
1864 //
1865 Int_t index0=0;
c6f438c0 1866 if (!cseed[0].IsOK()){
69b55c55 1867 index0 = 1;
c6f438c0 1868 if (!cseed[1].IsOK()) index0 = 2;
69b55c55 1869 }
1870 seedparams[registered][0] = cseed[index0].fX0;
1871 seedparams[registered][1] = cseed[index0].fYref[0];
1872 seedparams[registered][2] = cseed[index0].fZref[0];
c6f438c0 1873 seedparams[registered][5] = cR;
1874 seedparams[registered][3] = cseed[index0].fX0*cR - TMath::Sin(TMath::ATan(cseed[0].fYref[1]));
69b55c55 1875 seedparams[registered][4] = cseed[index0].fZref[1]/
1876 TMath::Sqrt(1+cseed[index0].fYref[1]*cseed[index0].fYref[1]);
1877 seedparams[registered][6] = ns;
1878 //
1879 //
1880 Int_t labels[12], outlab[24];
1881 Int_t nlab=0;
1882 for (Int_t iLayer=0;iLayer<6;iLayer++){
c6f438c0 1883 if (!cseed[iLayer].IsOK()) continue;
69b55c55 1884 if (cseed[iLayer].fLabels[0]>=0) {
1885 labels[nlab] = cseed[iLayer].fLabels[0];
1886 nlab++;
1887 }
1888 if (cseed[iLayer].fLabels[1]>=0) {
1889 labels[nlab] = cseed[iLayer].fLabels[1];
1890 nlab++;
1891 }
1892 }
1893 Freq(nlab,labels,outlab,kFALSE);
1894 Int_t label = outlab[0];
1895 Int_t frequency = outlab[1];
1896 for (Int_t iLayer=0;iLayer<6;iLayer++){
1897 cseed[iLayer].fFreq = frequency;
c6f438c0 1898 cseed[iLayer].fC = cR;
1899 cseed[iLayer].fCC = cC;
69b55c55 1900 cseed[iLayer].fChi2 = chi2TR;
1901 cseed[iLayer].fChi2Z = chi2ZF;
1902 }
1903 //
1904 if (1||(!isFake)){ //debugging print
1905 Float_t zvertex = GetZ();
1906 TTreeSRedirector& cstream = *fDebugStreamer;
d337ef8d 1907 if (AliTRDReconstructor::StreamLevel()>0)
69b55c55 1908 cstream<<"Seeds1"<<
1909 "isFake="<<isFake<<
1910 "Vertex="<<zvertex<<
1911 "Rieman2.="<<&rieman2<<
1912 "Rieman.="<<&rieman<<
1913 "Xref="<<xref<<
1914 "X0="<<xcl[0]<<
1915 "X1="<<xcl[1]<<
1916 "X2="<<xcl[2]<<
1917 "X3="<<xcl[3]<<
1918 "X4="<<xcl[4]<<
1919 "X5="<<xcl[5]<<
1920 "Chi2R="<<chi2R<<
1921 "Chi2Z="<<chi2Z<<
1922 "Chi2RF="<<chi2RF<< //chi2 of trackletes on full track
1923 "Chi2ZF="<<chi2ZF<< //chi2 z on tracklets on full track
1924 "Chi2ZT2="<<chi2ZT2<< //chi2 z on tracklets on full track - rieman tilt
1925 "Chi2ZTC="<<chi2ZTC<< //chi2 z on tracklets on full track - rieman tilt const
1926 //
1927 "Chi2TR="<<chi2TR<< //chi2 without vertex constrain
1928 "Chi2TC="<<chi2TC<< //chi2 with vertex constrain
1929 "C="<<curv<< // non constrained - no tilt correction
1930 "DR="<<dR<< // DR parameter - tilt correction
1931 "DCA="<<dca<< // DCA - tilt correction
c6f438c0 1932 "CR="<<cR<< // non constrained curvature - tilt correction
1933 "CC="<<cC<< // constrained curvature
69b55c55 1934 "Polz0="<<polz0c<<
1935 "Polz1="<<polz1c<<
1936 "RPolz0="<<rpolz0<<
1937 "RPolz1="<<rpolz1<<
1938 "Ncl="<<nclusters<<
1939 "Nlayers="<<nlayers<<
1940 "NUsedS="<<nusedCl<<
1941 "NUsed="<<nusedf<<
1942 "Findable="<<findable<<
1943 "Like="<<like<<
1944 "LikePrim="<<likePrim<<
1945 "Likechi2C="<<likechi2C<<
1946 "Likechi2TR="<<likechi2TR<<
1947 "Likezf="<<likezf<<
1948 "LikeF="<<seedquality2[registered]<<
1949 "S0.="<<&cseed[0]<<
1950 "S1.="<<&cseed[1]<<
1951 "S2.="<<&cseed[2]<<
1952 "S3.="<<&cseed[3]<<
1953 "S4.="<<&cseed[4]<<
1954 "S5.="<<&cseed[5]<<
1955 "SB0.="<<&seedb[0]<<
1956 "SB1.="<<&seedb[1]<<
1957 "SB2.="<<&seedb[2]<<
1958 "SB3.="<<&seedb[3]<<
1959 "SB4.="<<&seedb[4]<<
1960 "SB5.="<<&seedb[5]<<
1961 "Label="<<label<<
1962 "Freq="<<frequency<<
1963 "sLayer="<<sLayer<<
1964 "\n";
1965 }
c6f438c0 1966 if (registered<kMaxSeed-1) {
69b55c55 1967 registered++;
1968 cseed = seed[registered];
1969 }
1970 }// end of loop over layer 1
1971 } // end of loop over layer 0
1972 } // end of loop over layer 3
1973 } // end of loop over seeding time bins
1974 //
1975 // choos best
1976 //
1977 TMath::Sort(registered,seedquality2,sort,kTRUE);
c6f438c0 1978 Bool_t signedseed[kMaxSeed];
69b55c55 1979 for (Int_t i=0;i<registered;i++){
1980 signedseed[i]= kFALSE;
1981 }
1982 for (Int_t iter=0; iter<5; iter++){
1983 for (Int_t iseed=0;iseed<registered;iseed++){
1984 Int_t index = sort[iseed];
1985 if (signedseed[index]) continue;
1986 Int_t labelsall[1000];
1987 Int_t nlabelsall=0;
1988 Int_t naccepted=0;;
1989 Int_t sLayer = seedlayer[index];
1990 Int_t ncl = 0;
1991 Int_t nused = 0;
1992 Int_t nlayers =0;
1993 Int_t findable = 0;
1994 for (Int_t jLayer=0;jLayer<6;jLayer++){
1995 if (TMath::Abs(seed[index][jLayer].fYref[0]/xcl[jLayer])<0.15)
1996 findable++;
c6f438c0 1997 if (seed[index][jLayer].IsOK()){
69b55c55 1998 seed[index][jLayer].UpdateUsed();
1999 ncl +=seed[index][jLayer].fN2;
2000 nused +=seed[index][jLayer].fNUsed;
2001 nlayers++;
2002 //cooking label
2003 for (Int_t itime=0;itime<25;itime++){
2004 if (seed[index][jLayer].fUsable[itime]){
2005 naccepted++;
2006 for (Int_t ilab=0;ilab<3;ilab++){
2007 Int_t tindex = seed[index][jLayer].fClusters[itime]->GetLabel(ilab);
2008 if (tindex>=0){
2009 labelsall[nlabelsall] = tindex;
2010 nlabelsall++;
2011 }
2012 }
2013 }
2014 }
2015 }
2016 }
7ad19338 2017 //
69b55c55 2018 if (nused>30) continue;
7ad19338 2019 //
69b55c55 2020 if (iter==0){
2021 if (nlayers<6) continue;
2022 if (TMath::Log(0.000000001+seedquality2[index])<-5.) continue; // gold
2023 }
2024 //
2025 if (iter==1){
2026 if (nlayers<findable) continue;
2027 if (TMath::Log(0.000000001+seedquality2[index])<-4.) continue; //
7ad19338 2028 }
7ad19338 2029 //
7ad19338 2030 //
69b55c55 2031 if (iter==2){
2032 if (nlayers==findable || nlayers==6) continue;
2033 if (TMath::Log(0.000000001+seedquality2[index])<-6.) continue;
2034 }
7ad19338 2035 //
69b55c55 2036 if (iter==3){
2037 if (TMath::Log(0.000000001+seedquality2[index])<-5.) continue;
2038 }
7ad19338 2039 //
69b55c55 2040 if (iter==4){
2041 if (TMath::Log(0.000000001+seedquality2[index])-nused/(nlayers-3.)<-15.) continue;
2042 }
7ad19338 2043 //
69b55c55 2044 signedseed[index] = kTRUE;
2045 //
2046 Int_t labels[1000], outlab[1000];
2047 Int_t nlab=0;
2048 for (Int_t iLayer=0;iLayer<6;iLayer++){
c6f438c0 2049 if (seed[index][iLayer].IsOK()){
69b55c55 2050 if (seed[index][iLayer].fLabels[0]>=0) {
2051 labels[nlab] = seed[index][iLayer].fLabels[0];
2052 nlab++;
2053 }
2054 if (seed[index][iLayer].fLabels[1]>=0) {
2055 labels[nlab] = seed[index][iLayer].fLabels[1];
2056 nlab++;
2057 }
2058 }
7ad19338 2059 }
69b55c55 2060 Freq(nlab,labels,outlab,kFALSE);
2061 Int_t label = outlab[0];
2062 Int_t frequency = outlab[1];
2063 Freq(nlabelsall,labelsall,outlab,kFALSE);
2064 Int_t label1 = outlab[0];
2065 Int_t label2 = outlab[2];
2066 Float_t fakeratio = (naccepted-outlab[1])/Float_t(naccepted);
2067 Float_t ratio = Float_t(nused)/Float_t(ncl);
2068 if (ratio<0.25){
2069 for (Int_t jLayer=0;jLayer<6;jLayer++){
c6f438c0 2070 if ( seed[index][jLayer].IsOK()&&TMath::Abs(seed[index][jLayer].fYfit[1]-seed[index][jLayer].fYfit[1])<0.2 )
69b55c55 2071 seed[index][jLayer].UseClusters(); //sign gold
2072 }
7ad19338 2073 }
2074 //
69b55c55 2075 Int_t eventNr = esd->GetEventNumber();
2076 TTreeSRedirector& cstream = *fDebugStreamer;
2077 //
2078 // register seed
2079 //
2080 AliTRDtrack * track = RegisterSeed(seed[index],seedparams[index]);
2081 AliTRDtrack dummy;
2082 if (!track) track=&dummy;
2083 else{
2084 AliESDtrack esdtrack;
2085 esdtrack.UpdateTrackParams(track, AliESDtrack::kTRDout);
2086 esdtrack.SetLabel(label);
2087 esd->AddTrack(&esdtrack);
7ad19338 2088 TTreeSRedirector& cstream = *fDebugStreamer;
d337ef8d 2089 if (AliTRDReconstructor::StreamLevel()>0)
2090 cstream<<"Tracks"<<
2091 "EventNr="<<eventNr<<
2092 "ESD.="<<&esdtrack<<
2093 "trd.="<<track<<
2094 "trdback.="<<track<<
2095 "\n";
7ad19338 2096 }
d337ef8d 2097 if (AliTRDReconstructor::StreamLevel()>0)
2098 cstream<<"Seeds2"<<
69b55c55 2099 "Iter="<<iter<<
2100 "Track.="<<track<<
2101 "Like="<<seedquality[index]<<
2102 "LikeF="<<seedquality2[index]<<
2103 "S0.="<<&seed[index][0]<<
2104 "S1.="<<&seed[index][1]<<
2105 "S2.="<<&seed[index][2]<<
2106 "S3.="<<&seed[index][3]<<
2107 "S4.="<<&seed[index][4]<<
2108 "S5.="<<&seed[index][5]<<
2109 "Label="<<label<<
2110 "Label1="<<label1<<
2111 "Label2="<<label2<<
2112 "FakeRatio="<<fakeratio<<
2113 "Freq="<<frequency<<
2114 "Ncl="<<ncl<<
2115 "Nlayers="<<nlayers<<
2116 "Findable="<<findable<<
2117 "NUsed="<<nused<<
2118 "sLayer="<<sLayer<<
2119 "EventNr="<<eventNr<<
2120 "\n";
7ad19338 2121 }
2122 }
69b55c55 2123 } // end of loop over sectors
75bd7f81 2124
69b55c55 2125 delete [] pseed;
75bd7f81 2126
69b55c55 2127}
2128
5443e65e 2129//_____________________________________________________________________________
b7a0917f 2130Int_t AliTRDtracker::ReadClusters(TObjArray *array, TTree *ClusterTree) const
5443e65e 2131{
2132 //
a819a5f7 2133 // Reads AliTRDclusters (option >= 0) or AliTRDrecPoints (option < 0)
2134 // from the file. The names of the cluster tree and branches
2135 // should match the ones used in AliTRDclusterizer::WriteClusters()
2136 //
75bd7f81 2137
4f1c04d3 2138 Int_t nsize = Int_t(ClusterTree->GetTotBytes()/(sizeof(AliTRDcluster)));
2139 TObjArray *clusterArray = new TObjArray(nsize+1000);
5443e65e 2140
c630aafd 2141 TBranch *branch=ClusterTree->GetBranch("TRDcluster");
2142 if (!branch) {
2143 Error("ReadClusters","Can't get the branch !");
2144 return 1;
2145 }
029cd327 2146 branch->SetAddress(&clusterArray);
5443e65e 2147
2148 Int_t nEntries = (Int_t) ClusterTree->GetEntries();
19dd5b2f 2149 // printf("found %d entries in %s.\n",nEntries,ClusterTree->GetName());
a819a5f7 2150
a819a5f7 2151 // Loop through all entries in the tree
eb187bed 2152 Int_t nbytes = 0;
a819a5f7 2153 AliTRDcluster *c = 0;
7bed16a7 2154 // printf("\n");
a819a5f7 2155 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
2156
2157 // Import the tree
5443e65e 2158 nbytes += ClusterTree->GetEvent(iEntry);
2159
a819a5f7 2160 // Get the number of points in the detector
029cd327 2161 Int_t nCluster = clusterArray->GetEntriesFast();
e24ea474 2162// printf("\r Read %d clusters from entry %d", nCluster, iEntry);
5443e65e 2163
a819a5f7 2164 // Loop through all TRD digits
2165 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
029cd327 2166 c = (AliTRDcluster*)clusterArray->UncheckedAt(iCluster);
4f1c04d3 2167 AliTRDcluster *co = c;
a819a5f7 2168 array->AddLast(co);
4f1c04d3 2169 // delete clusterArray->RemoveAt(iCluster);
2170 clusterArray->RemoveAt(iCluster);
a819a5f7 2171 }
2172 }
7c1698cb 2173// cout<<"Allocated"<<nsize<<"\tLoaded"<<array->GetEntriesFast()<<"\n";
a819a5f7 2174
029cd327 2175 delete clusterArray;
5443e65e 2176
c630aafd 2177 return 0;
75bd7f81 2178
a819a5f7 2179}
2180
75bd7f81 2181//_____________________________________________________________________________
3551db50 2182Bool_t AliTRDtracker::GetTrackPoint(Int_t index, AliTrackPoint& p) const
2183{
2184 //
2185 // Get track space point with index i
2186 // Origin: C.Cheshkov
2187 //
2188
2189 AliTRDcluster *cl = (AliTRDcluster*)fClusters->UncheckedAt(index);
2190 Int_t idet = cl->GetDetector();
2191 Int_t isector = fGeom->GetSector(idet);
2192 Int_t ichamber= fGeom->GetChamber(idet);
2193 Int_t iplan = fGeom->GetPlane(idet);
2194 Double_t local[3];
2195 local[0]=GetX(isector,iplan,cl->GetLocalTimeBin());
2196 local[1]=cl->GetY();
2197 local[2]=cl->GetZ();
2198 Double_t global[3];
2199 fGeom->RotateBack(idet,local,global);
2200 p.SetXYZ(global[0],global[1],global[2]);
2201 AliAlignObj::ELayerID iLayer = AliAlignObj::kTRD1;
2202 switch (iplan) {
2203 case 0:
2204 iLayer = AliAlignObj::kTRD1;
2205 break;
2206 case 1:
2207 iLayer = AliAlignObj::kTRD2;
2208 break;
2209 case 2:
2210 iLayer = AliAlignObj::kTRD3;
2211 break;
2212 case 3:
2213 iLayer = AliAlignObj::kTRD4;
2214 break;
2215 case 4:
2216 iLayer = AliAlignObj::kTRD5;
2217 break;
2218 case 5:
2219 iLayer = AliAlignObj::kTRD6;
2220 break;
2221 };
2222 Int_t modId = isector*fGeom->Ncham()+ichamber;
2223 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer,modId);
2224 p.SetVolumeID(volid);
2225
2226 return kTRUE;
2227
2228}
2229
75bd7f81 2230//_____________________________________________________________________________
029cd327 2231void AliTRDtracker::CookLabel(AliKalmanTrack* pt, Float_t wrong) const
2232{
2233 //
2234 // This cooks a label. Mmmmh, smells good...
2235 //
46d29e70 2236
2237 Int_t label=123456789, index, i, j;
5443e65e 2238 Int_t ncl=pt->GetNumberOfClusters();
029cd327 2239 const Int_t kRange = fTrSec[0]->GetOuterTimeBin()+1;
5443e65e 2240
029cd327 2241 Bool_t labelAdded;
46d29e70 2242
029cd327 2243 // Int_t s[kRange][2];
2244 Int_t **s = new Int_t* [kRange];
2245 for (i=0; i<kRange; i++) {
d1b06c24 2246 s[i] = new Int_t[2];
2247 }
029cd327 2248 for (i=0; i<kRange; i++) {
46d29e70 2249 s[i][0]=-1;
2250 s[i][1]=0;
2251 }
2252
2253 Int_t t0,t1,t2;
2254 for (i=0; i<ncl; i++) {
5443e65e 2255 index=pt->GetClusterIndex(i);
46d29e70 2256 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
5443e65e 2257 t0=c->GetLabel(0);
2258 t1=c->GetLabel(1);
2259 t2=c->GetLabel(2);
46d29e70 2260 }
2261
2262 for (i=0; i<ncl; i++) {
5443e65e 2263 index=pt->GetClusterIndex(i);
46d29e70 2264 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
2265 for (Int_t k=0; k<3; k++) {
5443e65e 2266 label=c->GetLabel(k);
029cd327 2267 labelAdded=kFALSE; j=0;
46d29e70 2268 if (label >= 0) {
029cd327 2269 while ( (!labelAdded) && ( j < kRange ) ) {
a9814c09 2270 if (s[j][0]==label || s[j][1]==0) {
2271 s[j][0]=label;
2272 s[j][1]=s[j][1]+1;
029cd327 2273 labelAdded=kTRUE;
a9814c09 2274 }
2275 j++;
2276 }
46d29e70 2277 }
2278 }
2279 }
2280
46d29e70 2281 Int_t max=0;
2282 label = -123456789;
2283
029cd327 2284 for (i=0; i<kRange; i++) {
46d29e70 2285 if (s[i][1]>max) {
2286 max=s[i][1]; label=s[i][0];
2287 }
2288 }
5443e65e 2289
029cd327 2290 for (i=0; i<kRange; i++) {
5443e65e 2291 delete []s[i];
2292 }
2293
d1b06c24 2294 delete []s;
5443e65e 2295
2296 if ((1.- Float_t(max)/ncl) > wrong) label=-label;
2297
2298 pt->SetLabel(label);
2299
46d29e70 2300}
2301
75bd7f81 2302//_____________________________________________________________________________
029cd327 2303void AliTRDtracker::UseClusters(const AliKalmanTrack* t, Int_t from) const
2304{
2305 //
2306 // Use clusters, but don't abuse them!
2307 //
75bd7f81 2308
69b55c55 2309 const Float_t kmaxchi2 =18;
2310 const Float_t kmincl =10;
2311 AliTRDtrack * track = (AliTRDtrack*)t;
2312 //
5443e65e 2313 Int_t ncl=t->GetNumberOfClusters();
2314 for (Int_t i=from; i<ncl; i++) {
2315 Int_t index = t->GetClusterIndex(i);
2316 AliTRDcluster *c=(AliTRDcluster*)fClusters->UncheckedAt(index);
69b55c55 2317 //
2318 Int_t iplane = fGeom->GetPlane(c->GetDetector());
2319 if (track->fTracklets[iplane].GetChi2()>kmaxchi2) continue;
2320 if (track->fTracklets[iplane].GetN()<kmincl) continue;
2321 if (!(c->IsUsed())) c->Use();
5443e65e 2322 }
5443e65e 2323
75bd7f81 2324}
5443e65e 2325
75bd7f81 2326//_____________________________________________________________________________
029cd327 2327Double_t AliTRDtracker::ExpectedSigmaY2(Double_t , Double_t , Double_t ) const
5443e65e 2328{
75bd7f81 2329 //
5443e65e 2330 // Parametrised "expected" error of the cluster reconstruction in Y
75bd7f81 2331 //
5443e65e 2332
2333 Double_t s = 0.08 * 0.08;
2334 return s;
75bd7f81 2335
5443e65e 2336}
2337
75bd7f81 2338//_____________________________________________________________________________
029cd327 2339Double_t AliTRDtracker::ExpectedSigmaZ2(Double_t , Double_t ) const
0a29d0f1 2340{
75bd7f81 2341 //
5443e65e 2342 // Parametrised "expected" error of the cluster reconstruction in Z
75bd7f81 2343 //
5443e65e 2344
a9814c09 2345 Double_t s = 9 * 9 /12.;
5443e65e 2346 return s;
75bd7f81 2347
5443e65e 2348}
2349
75bd7f81 2350//_____________________________________________________________________________
029cd327 2351Double_t AliTRDtracker::GetX(Int_t sector, Int_t plane, Int_t localTB) const
5443e65e 2352{
2353 //
029cd327 2354 // Returns radial position which corresponds to time bin <localTB>
5443e65e 2355 // in tracking sector <sector> and plane <plane>
2356 //
2357
029cd327 2358 Int_t index = fTrSec[sector]->CookTimeBinIndex(plane, localTB);
5443e65e 2359 Int_t pl = fTrSec[sector]->GetLayerNumber(index);
2360 return fTrSec[sector]->GetLayer(pl)->GetX();
2361
2362}
2363
75bd7f81 2364//_____________________________________________________________________________
2365AliTRDtracker::AliTRDpropagationLayer
2366 ::AliTRDpropagationLayer(Double_t x, Double_t dx, Double_t rho
2367 , Double_t radLength, Int_t tbIndex, Int_t plane)
5443e65e 2368{
0a29d0f1 2369 //
5443e65e 2370 // AliTRDpropagationLayer constructor
0a29d0f1 2371 //
46d29e70 2372
029cd327 2373 fN = 0; fX = x; fdX = dx; fRho = rho; fX0 = radLength;
2374 fClusters = NULL; fIndex = NULL; fTimeBinIndex = tbIndex;
59393e34 2375 fPlane = plane;
46d29e70 2376
029cd327 2377 for(Int_t i=0; i < (Int_t) kZones; i++) {
5443e65e 2378 fZc[i]=0; fZmax[i] = 0;
a819a5f7 2379 }
5443e65e 2380
2381 fYmax = 0;
2382
2383 if(fTimeBinIndex >= 0) {
029cd327 2384 fClusters = new AliTRDcluster*[kMaxClusterPerTimeBin];
2385 fIndex = new UInt_t[kMaxClusterPerTimeBin];
a819a5f7 2386 }
46d29e70 2387
3c625a9b 2388 for (Int_t i=0;i<5;i++) fIsHole[i] = kFALSE;
5443e65e 2389 fHole = kFALSE;
2390 fHoleZc = 0;
2391 fHoleZmax = 0;
2392 fHoleYc = 0;
2393 fHoleYmax = 0;
2394 fHoleRho = 0;
2395 fHoleX0 = 0;
2396
2397}
2398
75bd7f81 2399//_____________________________________________________________________________
2400void AliTRDtracker::AliTRDpropagationLayer
2401 ::SetHole(Double_t Zmax, Double_t Ymax, Double_t rho
2402 , Double_t radLength, Double_t Yc, Double_t Zc)
5443e65e 2403{
2404 //
2405 // Sets hole in the layer
2406 //
75bd7f81 2407
5443e65e 2408 fHole = kTRUE;
2409 fHoleZc = Zc;
2410 fHoleZmax = Zmax;
2411 fHoleYc = Yc;
2412 fHoleYmax = Ymax;
2413 fHoleRho = rho;
029cd327 2414 fHoleX0 = radLength;
75bd7f81 2415
5443e65e 2416}
46d29e70 2417
75bd7f81 2418//_____________________________________________________________________________
2419AliTRDtracker::AliTRDtrackingSector
2420 ::AliTRDtrackingSector(AliTRDgeometry* geo, Int_t gs)
5443e65e 2421{
2422 //
2423 // AliTRDtrackingSector Constructor
2424 //
75bd7f81 2425
a5cadd36 2426 AliTRDpadPlane *padPlane = 0;
2427
5443e65e 2428 fGeom = geo;
5443e65e 2429 fGeomSector = gs;
5443e65e 2430 fN = 0;
3c625a9b 2431 //
2432 // get holes description from geometry
2433 Bool_t holes[AliTRDgeometry::kNcham];
2434 //printf("sector\t%d\t",gs);
2435 for (Int_t icham=0; icham<AliTRDgeometry::kNcham;icham++){
2436 holes[icham] = fGeom->IsHole(0,icham,gs);
2437 //printf("%d",holes[icham]);
2438 }
2439 //printf("\n");
2440
029cd327 2441 for(UInt_t i=0; i < kMaxTimeBinIndex; i++) fTimeBinIndex[i] = -1;
5443e65e 2442
2443
2444 AliTRDpropagationLayer* ppl;
2445
59393e34 2446 Double_t x, dx, rho, radLength;
2447 // Int_t steps;
5443e65e 2448
2449 // add layers for each of the planes
5443e65e 2450 Double_t dxAmp = (Double_t) fGeom->CamHght(); // Amplification region
a305677e 2451 //Double_t dxDrift = (Double_t) fGeom->CdrHght(); // Drift region
5443e65e 2452
a305677e 2453 Int_t tbIndex;
029cd327 2454 const Int_t kNchambers = AliTRDgeometry::Ncham();
3c625a9b 2455 Double_t ymax = 0;
3c625a9b 2456 Double_t ymaxsensitive=0;
029cd327 2457 Double_t *zc = new Double_t[kNchambers];
2458 Double_t *zmax = new Double_t[kNchambers];
3c625a9b 2459 Double_t *zmaxsensitive = new Double_t[kNchambers];
5443e65e 2460
3551db50 2461 AliTRDCommonParam* commonParam = AliTRDCommonParam::Instance();
2462 if (!commonParam)
2463 {
2464 printf("<AliTRDtracker::AliTRDtrackingSector::AliTRDtrackingSector> ");
2465 printf("Could not get common params\n");
2466 return;
2467 }
2468
5443e65e 2469 for(Int_t plane = 0; plane < AliTRDgeometry::Nplan(); plane++) {
5443e65e 2470
3c625a9b 2471 ymax = fGeom->GetChamberWidth(plane)/2.;
a5cadd36 2472 // Modidified for new pad plane class, 22.04.05 (C.B.)
3551db50 2473 padPlane = commonParam->GetPadPlane(plane,0);
59393e34 2474 ymaxsensitive = (padPlane->GetColSize(1)*padPlane->GetNcols()-4)/2.;
029cd327 2475 for(Int_t ch = 0; ch < kNchambers; ch++) {
2476 zmax[ch] = fGeom->GetChamberLength(plane,ch)/2;
a5cadd36 2477 //
2478 // Modidified for new pad plane class, 22.04.05 (C.B.)
a5cadd36 2479 Float_t pad = padPlane->GetRowSize(1);
3551db50 2480 Float_t row0 = commonParam->GetRow0(plane,ch,0);
2481 Int_t nPads = commonParam->GetRowMax(plane,ch,0);
857b3eb0 2482 zmaxsensitive[ch] = Float_t(nPads)*pad/2.;
7ad19338 2483 zc[ch] = -(pad * nPads)/2 + row0;
5443e65e 2484 }
2485
59393e34 2486 dx = AliTRDcalibDB::Instance()->GetVdrift(0,0,0)
3551db50 2487 / AliTRDcalibDB::Instance()->GetSamplingFrequency();
029cd327 2488 rho = 0.00295 * 0.85; radLength = 11.0;
5443e65e 2489
3551db50 2490 Double_t x0 = (Double_t) AliTRDgeometry::GetTime0(plane);
a305677e 2491 //Double_t xbottom = x0 - dxDrift;
2492 //Double_t xtop = x0 + dxAmp;
3c625a9b 2493 //
59393e34 2494 Int_t nTimeBins = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
2495 for (Int_t iTime = 0; iTime<nTimeBins; iTime++){
2496 Double_t xlayer = iTime*dx - dxAmp;
2497 //if (xlayer<0) xlayer=dxAmp/2.;
2498 x = x0 - xlayer;
2499 //
2500 tbIndex = CookTimeBinIndex(plane, iTime);
2501 ppl = new AliTRDpropagationLayer(x,dx,rho,radLength,tbIndex, plane);
3c625a9b 2502 ppl->SetYmax(ymax,ymaxsensitive);
2503 ppl->SetZ(zc, zmax, zmaxsensitive);
2504 ppl->SetHoles(holes);
59393e34 2505 InsertLayer(ppl);
5443e65e 2506 }
2507 }
2508
5443e65e 2509 MapTimeBinLayers();
029cd327 2510 delete [] zc;
2511 delete [] zmax;
4f1c04d3 2512 delete [] zmaxsensitive;
5443e65e 2513
2514}
2515
75bd7f81 2516//_____________________________________________________________________________
2517Int_t AliTRDtracker::AliTRDtrackingSector
2518 ::CookTimeBinIndex(Int_t plane, Int_t localTB) const
5443e65e 2519{
2520 //
2521 // depending on the digitization parameters calculates "global"
029cd327 2522 // time bin index for timebin <localTB> in plane <plane>
5443e65e 2523 //
59393e34 2524 //
75bd7f81 2525
59393e34 2526 Int_t tbPerPlane = AliTRDcalibDB::Instance()->GetNumberOfTimeBins();
7b580082 2527 Int_t gtb = (plane+1) * tbPerPlane - localTB -1;
59393e34 2528 if (localTB<0) return -1;
75bd7f81 2529 if (gtb<0) return -1;
2530
5443e65e 2531 return gtb;
5443e65e 2532
75bd7f81 2533}
5443e65e 2534
75bd7f81 2535//_____________________________________________________________________________
2536void AliTRDtracker::AliTRDtrackingSector
2537 ::MapTimeBinLayers()
5443e65e 2538{
2539 //
2540 // For all sensitive time bins sets corresponding layer index
2541 // in the array fTimeBins
2542 //
2543
2544 Int_t index;
2545
2546 for(Int_t i = 0; i < fN; i++) {
2547 index = fLayers[i]->GetTimeBinIndex();
2548
2549 // printf("gtb %d -> pl %d -> x %f \n", index, i, fLayers[i]->GetX());
2550
2551 if(index < 0) continue;
029cd327 2552 if(index >= (Int_t) kMaxTimeBinIndex) {
5443e65e 2553 printf("*** AliTRDtracker::MapTimeBinLayers: \n");
2554 printf(" index %d exceeds allowed maximum of %d!\n",
029cd327 2555 index, kMaxTimeBinIndex-1);
5443e65e 2556 continue;
2557 }
2558 fTimeBinIndex[index] = i;
2559 }
5443e65e 2560
75bd7f81 2561}
5443e65e 2562
75bd7f81 2563//_____________________________________________________________________________
2564Int_t AliTRDtracker::AliTRDtrackingSector
2565 ::GetLayerNumber(Double_t x) const
5443e65e 2566{
2567 //
2568 // Returns the number of time bin which in radial position is closest to <x>
2569 //
2570
2571 if(x >= fLayers[fN-1]->GetX()) return fN-1;
2572 if(x <= fLayers[0]->GetX()) return 0;
2573
2574 Int_t b=0, e=fN-1, m=(b+e)/2;
2575 for (; b<e; m=(b+e)/2) {
2576 if (x > fLayers[m]->GetX()) b=m+1;
2577 else e=m;
2578 }
2579 if(TMath::Abs(x - fLayers[m]->GetX()) >
2580 TMath::Abs(x - fLayers[m+1]->GetX())) return m+1;
75bd7f81 2581
5443e65e 2582 else return m;
2583
2584}
2585
75bd7f81 2586//_____________________________________________________________________________
2587Int_t AliTRDtracker::AliTRDtrackingSector
2588 ::GetInnerTimeBin() const
5443e65e 2589{
2590 //
2591 // Returns number of the innermost SENSITIVE propagation layer
2592 //
2593
2594 return GetLayerNumber(0);
5443e65e 2595
75bd7f81 2596}
5443e65e 2597
75bd7f81 2598//_____________________________________________________________________________
2599Int_t AliTRDtracker::AliTRDtrackingSector
2600 ::GetOuterTimeBin() const
5443e65e 2601{
2602 //
2603 // Returns number of the outermost SENSITIVE time bin
2604 //
2605
2606 return GetLayerNumber(GetNumberOfTimeBins() - 1);
46d29e70 2607
75bd7f81 2608}
5443e65e 2609
75bd7f81 2610//_____________________________________________________________________________
2611Int_t AliTRDtracker::AliTRDtrackingSector
2612 ::GetNumberOfTimeBins() const
5443e65e 2613{
2614 //
2615 // Returns number of SENSITIVE time bins
2616 //
2617
2618 Int_t tb, layer;
029cd327 2619 for(tb = kMaxTimeBinIndex-1; tb >=0; tb--) {
5443e65e 2620 layer = GetLayerNumber(tb);
2621 if(layer>=0) break;
2622 }
75bd7f81 2623
5443e65e 2624 return tb+1;
5443e65e 2625
75bd7f81 2626}
5443e65e 2627
75bd7f81 2628//_____________________________________________________________________________
2629void AliTRDtracker::AliTRDtrackingSector
2630 ::InsertLayer(AliTRDpropagationLayer* pl)
5443e65e 2631{
2632 //
2633 // Insert layer <pl> in fLayers array.
2634 // Layers are sorted according to X coordinate.
75bd7f81 2635 //
5443e65e 2636
029cd327 2637 if ( fN == ((Int_t) kMaxLayersPerSector)) {
5443e65e 2638 printf("AliTRDtrackingSector::InsertLayer(): Too many layers !\n");
2639 return;
2640 }
2641 if (fN==0) {fLayers[fN++] = pl; return;}
2642 Int_t i=Find(pl->GetX());
2643
2644 memmove(fLayers+i+1 ,fLayers+i,(fN-i)*sizeof(AliTRDpropagationLayer*));
2645 fLayers[i]=pl; fN++;
2646
2647}
2648
75bd7f81 2649//_____________________________________________________________________________
2650Int_t AliTRDtracker::AliTRDtrackingSector
2651 ::Find(Double_t x) const
5443e65e 2652{
2653 //
2654 // Returns index of the propagation layer nearest to X
2655 //
2656
2657 if (x <= fLayers[0]->GetX()) return 0;
2658 if (x > fLayers[fN-1]->GetX()) return fN;
2659 Int_t b=0, e=fN-1, m=(b+e)/2;
2660 for (; b<e; m=(b+e)/2) {
2661 if (x > fLayers[m]->GetX()) b=m+1;
2662 else e=m;
2663 }
7ad19338 2664
75bd7f81 2665 return m;
7ad19338 2666
75bd7f81 2667}
7ad19338 2668
75bd7f81 2669//_____________________________________________________________________________
2670void AliTRDtracker::AliTRDpropagationLayer
2671 ::SetZ(Double_t* center, Double_t *w, Double_t *wsensitive )
3c625a9b 2672{
2673 //
2674 // set centers and the width of sectors
75bd7f81 2675 //
2676
3c625a9b 2677 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2678 fZc[icham] = center[icham];
2679 fZmax[icham] = w[icham];
2680 fZmaxSensitive[icham] = wsensitive[icham];
75bd7f81 2681 // printf("chamber\t%d\tzc\t%f\tzmax\t%f\tzsens\t%f\n",icham,fZc[icham],fZmax[icham],fZmaxSensitive[icham]);
3c625a9b 2682 }
75bd7f81 2683
3c625a9b 2684}
5443e65e 2685
75bd7f81 2686//_____________________________________________________________________________
3c625a9b 2687void AliTRDtracker::AliTRDpropagationLayer::SetHoles(Bool_t *holes)
2688{
2689 //
2690 // set centers and the width of sectors
75bd7f81 2691 //
2692
3c625a9b 2693 fHole = kFALSE;
2694 for (Int_t icham=0;icham< AliTRDgeometry::kNcham;icham++){
2695 fIsHole[icham] = holes[icham];
2696 if (holes[icham]) fHole = kTRUE;
2697 }
5443e65e 2698
75bd7f81 2699}
5443e65e 2700
75bd7f81 2701//_____________________________________________________________________________
2702void AliTRDtracker::AliTRDpropagationLayer
2703 ::InsertCluster(AliTRDcluster* c, UInt_t index)
2704{
2705 //
2706 // Insert cluster in cluster array.
2707 // Clusters are sorted according to Y coordinate.
2708 //
5443e65e 2709
2710 if(fTimeBinIndex < 0) {
2711 printf("*** attempt to insert cluster into non-sensitive time bin!\n");
2712 return;
2713 }
2714
029cd327 2715 if (fN== (Int_t) kMaxClusterPerTimeBin) {
5443e65e 2716 printf("AliTRDpropagationLayer::InsertCluster(): Too many clusters !\n");
2717 return;
2718 }
2719 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
2720 Int_t i=Find(c->GetY());
2721 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTRDcluster*));
2722 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
2723 fIndex[i]=index; fClusters[i]=c; fN++;
5443e65e 2724
75bd7f81 2725}
5443e65e 2726
75bd7f81 2727//_____________________________________________________________________________
2728Int_t AliTRDtracker::AliTRDpropagationLayer::Find(Float_t y) const
2729{
2730 //
2731 // Returns index of the cluster nearest in Y
2732 //
5443e65e 2733
69b55c55 2734 if (fN<=0) return 0;
5443e65e 2735 if (y <= fClusters[0]->GetY()) return 0;
2736 if (y > fClusters[fN-1]->GetY()) return fN;
2737 Int_t b=0, e=fN-1, m=(b+e)/2;
2738 for (; b<e; m=(b+e)/2) {
2739 if (y > fClusters[m]->GetY()) b=m+1;
2740 else e=m;
2741 }
75bd7f81 2742
5443e65e 2743 return m;
75bd7f81 2744
5443e65e 2745}
2746
75bd7f81 2747//_____________________________________________________________________________
2748Int_t AliTRDtracker::AliTRDpropagationLayer
2749 ::FindNearestCluster(Float_t y, Float_t z, Float_t maxroad
2750 , Float_t maxroadz) const
7ad19338 2751{
2752 //
2753 // Returns index of the cluster nearest to the given y,z
2754 //
75bd7f81 2755
7ad19338 2756 Int_t index = -1;
2757 Int_t maxn = fN;
69b55c55 2758 Float_t mindist = maxroad;
7ad19338 2759 //
2760 for (Int_t i=Find(y-maxroad); i<maxn; i++) {
2761 AliTRDcluster* c=(AliTRDcluster*)(fClusters[i]);
69b55c55 2762 Float_t ycl = c->GetY();
7ad19338 2763 //
69b55c55 2764 if (ycl > y+maxroad) break;
2765 if (TMath::Abs(c->GetZ()-z) > maxroadz) continue;
2766 if (TMath::Abs(ycl-y)<mindist){
2767 mindist = TMath::Abs(ycl-y);
2768 index = fIndex[i];
7ad19338 2769 }
2770 }
75bd7f81 2771
7ad19338 2772 return index;
7ad19338 2773
75bd7f81 2774}
7ad19338 2775
75bd7f81 2776//_____________________________________________________________________________
2777Double_t AliTRDtracker::GetTiltFactor(const AliTRDcluster* c)
2778{
2779 //
2780 // Returns correction factor for tilted pads geometry
2781 //
5443e65e 2782
fd621f36 2783 Int_t det = c->GetDetector();
2784 Int_t plane = fGeom->GetPlane(det);
3551db50 2785 AliTRDpadPlane *padPlane = AliTRDCommonParam::Instance()->GetPadPlane(plane,0);
de4b10e5 2786 Double_t h01 = TMath::Tan(-TMath::Pi() / 180.0 * padPlane->GetTiltingAngle());
b8dc2353 2787
2788 if(fNoTilt) h01 = 0;
75bd7f81 2789
fd621f36 2790 return h01;
5443e65e 2791
75bd7f81 2792}
c630aafd 2793
75bd7f81 2794//_____________________________________________________________________________
eab5961e 2795void AliTRDtracker::CookdEdxTimBin(AliTRDtrack& TRDtrack)
2796{
75bd7f81 2797 //
eab5961e 2798 // *** ADDED TO GET MORE INFORMATION FOR TRD PID ---- PS
2799 // This is setting fdEdxPlane and fTimBinPlane
2800 // Sums up the charge in each plane for track TRDtrack and also get the
2801 // Time bin for Max. Cluster
2802 // Prashant Shukla (shukla@physi.uni-heidelberg.de)
75bd7f81 2803 //
eab5961e 2804
6d45eaef 2805 Double_t clscharge[AliESDtrack::kNPlane][AliESDtrack::kNSlice];
2806 Double_t maxclscharge[AliESDtrack::kNPlane];
2807 Int_t nCluster[AliESDtrack::kNPlane][AliESDtrack::kNSlice];
2808 Int_t timebin[AliESDtrack::kNPlane];
2809
2810 //Initialization of cluster charge per plane.
2811 for (Int_t iPlane = 0; iPlane < AliESDtrack::kNPlane; iPlane++) {
2812 for (Int_t iSlice = 0; iSlice < AliESDtrack::kNSlice; iSlice++) {
2813 clscharge[iPlane][iSlice] = 0.0;
2814 nCluster[iPlane][iSlice] = 0;
2815 }
2816 }
eab5961e 2817
2818 //Initialization of cluster charge per plane.
f122c485 2819 for (Int_t iPlane = 0; iPlane < AliESDtrack::kNPlane; iPlane++) {
eab5961e 2820 timebin[iPlane] = -1;
2821 maxclscharge[iPlane] = 0.0;
2822 }
2823
2824 // Loop through all clusters associated to track TRDtrack
2825 Int_t nClus = TRDtrack.GetNumberOfClusters(); // from Kalmantrack
2826 for (Int_t iClus = 0; iClus < nClus; iClus++) {
2827 Double_t charge = TRDtrack.GetClusterdQdl(iClus);
2828 Int_t index = TRDtrack.GetClusterIndex(iClus);
c6f438c0 2829 AliTRDcluster *pTRDcluster = (AliTRDcluster *) GetCluster(index);
2830 if (!pTRDcluster) continue;
2831 Int_t tb = pTRDcluster->GetLocalTimeBin();
eab5961e 2832 if (!tb) continue;
c6f438c0 2833 Int_t detector = pTRDcluster->GetDetector();
eab5961e 2834 Int_t iPlane = fGeom->GetPlane(detector);
6d45eaef 2835 Int_t iSlice = tb*AliESDtrack::kNSlice/AliTRDtrack::kNtimeBins;
2836 clscharge[iPlane][iSlice] = clscharge[iPlane][iSlice]+charge;
eab5961e 2837 if(charge > maxclscharge[iPlane]) {
2838 maxclscharge[iPlane] = charge;
2839 timebin[iPlane] = tb;
2840 }
6d45eaef 2841 nCluster[iPlane][iSlice]++;
eab5961e 2842 } // end of loop over cluster
2843
2844 // Setting the fdEdxPlane and fTimBinPlane variabales
c6f438c0 2845 Double_t totalCharge = 0;
6d45eaef 2846
f122c485 2847 for (Int_t iPlane = 0; iPlane < AliESDtrack::kNPlane; iPlane++) {
6d45eaef 2848 for (Int_t iSlice = 0; iSlice < AliESDtrack::kNSlice; iSlice++) {
2849 if (nCluster[iPlane][iSlice]) clscharge[iPlane][iSlice] /= nCluster[iPlane][iSlice];
2850 TRDtrack.SetPIDsignals(clscharge[iPlane][iSlice], iPlane, iSlice);
2851 totalCharge= totalCharge+clscharge[iPlane][iSlice];
bd50219c 2852 }
6d45eaef 2853 TRDtrack.SetPIDTimBin(timebin[iPlane], iPlane);
eab5961e 2854 }
6d45eaef 2855
eab5961e 2856 // Int_t i;
2857 // Int_t nc=TRDtrack.GetNumberOfClusters();
2858 // Float_t dedx=0;
2859 // for (i=0; i<nc; i++) dedx += TRDtrack.GetClusterdQdl(i);
2860 // dedx /= nc;
2861 // for (Int_t iPlane = 0; iPlane < kNPlane; iPlane++) {
2862 // TRDtrack.SetPIDsignals(dedx, iPlane);
2863 // TRDtrack.SetPIDTimBin(timbin[iPlane], iPlane);
2864 // }
2865
75bd7f81 2866}
c630aafd 2867
75bd7f81 2868//_____________________________________________________________________________
2869Int_t AliTRDtracker::FindClusters(Int_t sector, Int_t t0, Int_t t1
2870 , AliTRDtrack * track
2871 , Int_t *clusters,AliTRDtracklet&tracklet)
4f1c04d3 2872{
2873 //
2874 //
75bd7f81 2875 // Try to find nearest clusters to the track in timebins from t0 to t1
4f1c04d3 2876 //
2877 //
7ad19338 2878 //
2879 // correction coeficients - depends on TRD parameters - to be changed according it
2880 //
2881
2882 Double_t x[100],yt[100],zt[100];
2883 Double_t xmean=0; //reference x
2884 Double_t dz[10][100],dy[10][100];
2885 Float_t zmean[100], nmean[100];
2886 Int_t clfound=0;
2887 Int_t indexes[10][100]; // indexes of the clusters in the road
2888 AliTRDcluster *cl[10][100]; // pointers to the clusters in the road
2889 Int_t best[10][100]; // index of best matching cluster
2890 //
2891 //
69b55c55 2892
75bd7f81 2893 for (Int_t it=0;it<100; it++){
4f1c04d3 2894 x[it]=0;
2895 yt[it]=0;
2896 zt[it]=0;
75bd7f81 2897 clusters[it]=-2;
7ad19338 2898 zmean[it]=0;
2899 nmean[it]=0;
2900 //
2901 for (Int_t ih=0;ih<10;ih++){
2902 indexes[ih][it]=-2; //reset indexes1
2903 cl[ih][it]=0;
2904 dz[ih][it]=-100;
2905 dy[ih][it]=-100;
2906 best[ih][it]=0;
2907 }
4f1c04d3 2908 }
2909 //
2910 Double_t x0 = track->GetX();
69b55c55 2911 Double_t sigmaz = TMath::Sqrt(TMath::Abs(track->GetSigmaZ2()));
4f1c04d3 2912 Int_t nall=0;
2913 Int_t nfound=0;
7ad19338 2914 Double_t h01 =0;
2915 Int_t plane =-1;
c6f438c0 2916 Int_t detector =-1;
7ad19338 2917 Float_t padlength=0;
2918 AliTRDtrack track2(*track);
2919 Float_t snpy = track->GetSnp();
2920 Float_t tany = TMath::Sqrt(snpy*snpy/(1.-snpy*snpy));
2921 if (snpy<0) tany*=-1;
2922 //
2923 Double_t sy2=ExpectedSigmaY2(x0,track->GetTgl(),track->GetPt());
2924 Double_t sz2=ExpectedSigmaZ2(x0,track->GetTgl());
2925 Double_t road = 15.*sqrt(track->GetSigmaY2() + sy2);
2926 if (road>6.) road=6.;
4f1c04d3 2927
7ad19338 2928 //
2929 for (Int_t it=0;it<t1-t0;it++){
2930 Double_t maxChi2[2]={fgkMaxChi2,fgkMaxChi2};
2931 AliTRDpropagationLayer& timeBin=*(fTrSec[sector]->GetLayer(it+t0));
4f1c04d3 2932 if (timeBin==0) continue; // no indexes1
2933 Int_t maxn = timeBin;
2934 x[it] = timeBin.GetX();
7ad19338 2935 track2.PropagateTo(x[it]);
2936 yt[it] = track2.GetY();
2937 zt[it] = track2.GetZ();
2938
2939 Double_t y=yt[it],z=zt[it];
4f1c04d3 2940 Double_t chi2 =1000000;
2941 nall++;
2942 //
7ad19338 2943 // find 2 nearest cluster at given time bin
2944 //
2945 //
4f1c04d3 2946 for (Int_t i=timeBin.Find(y-road); i<maxn; i++) {
2947 AliTRDcluster* c=(AliTRDcluster*)(timeBin[i]);
7ad19338 2948 h01 = GetTiltFactor(c);
2949 if (plane<0){
c6f438c0 2950 Int_t det = c->GetDetector();
7ad19338 2951 plane = fGeom->GetPlane(det);
2952 padlength = TMath::Sqrt(c->GetSigmaZ2()*12.);
2953 }
2954 // if (c->GetLocalTimeBin()==0) continue;
4f1c04d3 2955 if (c->GetY() > y+road) break;
7ad19338 2956 if((c->GetZ()-z)*(c->GetZ()-z) > 12. * sz2) continue;
2957
2958 Double_t dist = TMath::Abs(c->GetZ()-z);
2959 if (dist> (0.5*padlength+6.*sigmaz)) continue; // 6 sigma boundary cut
2960 Double_t cost = 0;
2961 //
2962 if (dist> (0.5*padlength-sigmaz)){ // sigma boundary cost function
2963 cost = (dist-0.5*padlength)/(2.*sigmaz);
2964 if (cost>-1) cost= (cost+1.)*(cost+1.);
2965 else cost=0;
2966 }
2967 // Int_t label = TMath::Abs(track->GetLabel());
2968 // if (c->GetLabel(0)!=label && c->GetLabel(1)!=label&&c->GetLabel(2)!=label) continue;
2969 chi2=track2.GetPredictedChi2(c,h01)+cost;
2970 //
2971 clfound++;
2972 if (chi2 > maxChi2[1]) continue;
c6f438c0 2973 detector = c->GetDetector();
7ad19338 2974
2975 for (Int_t ih=2;ih<9; ih++){ //store the clusters in the road
2976 if (cl[ih][it]==0){
2977 cl[ih][it] = c;
2978 indexes[ih][it] =timeBin.GetIndex(i); // index - 9 - reserved for outliers
2979 break;
2980 }
4f1c04d3 2981 }
7ad19338 2982 //
2983 if (chi2 <maxChi2[0]){
2984 maxChi2[1] = maxChi2[0];
2985 maxChi2[0] = chi2;
2986 indexes[1][it] = indexes[0][it];
2987 cl[1][it] = cl[0][it];
2988 indexes[0][it] = timeBin.GetIndex(i);
2989 cl[0][it] = c;
2990 continue;
2991 }
2992 maxChi2[1]=chi2;
2993 cl[1][it] = c;
2994 indexes[1][it] =timeBin.GetIndex(i);
2995 }
2996 if (cl[0][it]){
2997 nfound++;
2998 xmean += x[it];
2999 }
4f1c04d3 3000 }
3001 //
7ad19338 3002 if (nfound<4) return 0;
3003 xmean /=Float_t(nfound); // middle x
3004 track2.PropagateTo(xmean); // propagate track to the center
4f1c04d3 3005 //
3006 // choose one of the variants
3007 //
7ad19338 3008 Int_t changes[10];
3009 Float_t sumz = 0;
3010 Float_t sum = 0;
3011 Double_t sumdy = 0;
3012 Double_t sumdy2 = 0;
3013 Double_t sumx = 0;
3014 Double_t sumxy = 0;
3015 Double_t sumx2 = 0;
3016 Double_t mpads = 0;
3017 //
3018 Int_t ngood[10];
3019 Int_t nbad[10];
3020 //
3021 Double_t meanz[10];
3022 Double_t moffset[10]; // mean offset
3023 Double_t mean[10]; // mean value
3024 Double_t angle[10]; // angle
3025 //
3026 Double_t smoffset[10]; // sigma of mean offset
3027 Double_t smean[10]; // sigma of mean value
3028 Double_t sangle[10]; // sigma of angle
3029 Double_t smeanangle[10]; // correlation
3030 //
3031 Double_t sigmas[10];
3032 Double_t tchi2s[10]; // chi2s for tracklet
75bd7f81 3033
3034 for (Int_t it=0;it<10;it++) {
3035
3036 ngood[it] = 0;
3037 nbad[it] = 0;
3038 //
3039 meanz[it] = 0;
3040 moffset[it] = 0; // mean offset
3041 mean[it] = 0; // mean value
3042 angle[it] = 0; // angle
3043 //
3044 smoffset[it] = 1e10; // sigma of mean offset
3045 smean[it] = 1e10; // sigma of mean value
3046 sangle[it] = 1e10; // sigma of angle
3047 smeanangle[it] = 0; // correlation
3048 //
3049 sigmas[it] = 1e10;
3050 tchi2s[it] = 1e10; // chi2s for tracklet
3051
3052 }
3053
7ad19338 3054 //
3055 // calculate zmean
3056 //
3057 for (Int_t it=0;it<t1-t0;it++){
3058 if (!cl[0][it]) continue;
3059 for (Int_t dt=-3;dt<=3;dt++){
3060 if (it+dt<0) continue;
8979685e 3061 if (it+dt>t1-t0) continue;
7ad19338 3062 if (!cl[0][it+dt]) continue;
3063 zmean[it]+=cl[0][it+dt]->GetZ();
3064 nmean[it]+=1.;
3065 }
3066 zmean[it]/=nmean[it];
3067 }
3068 //
3069 for (Int_t it=0; it<t1-t0;it++){
3070 best[0][it]=0;
3071 for (Int_t ih=0;ih<10;ih++){
3072 dz[ih][it]=-100;
3073 dy[ih][it]=-100;
4f1c04d3 3074 if (!cl[ih][it]) continue;
59393e34 3075 Double_t xcluster = cl[ih][it]->GetX();
3076 Double_t ytrack,ztrack;
3077 track2.GetProlongation(xcluster, ytrack, ztrack );
3078 dz[ih][it] = cl[ih][it]->GetZ()- ztrack; // calculate distance from track in z
3079 dy[ih][it] = cl[ih][it]->GetY()+ dz[ih][it]*h01 -ytrack; // in y
7ad19338 3080 }
3081 // minimize changes
3082 if (!cl[0][it]) continue;
3083 if (TMath::Abs(cl[0][it]->GetZ()-zmean[it])> padlength*0.8 &&cl[1][it])
3084 if (TMath::Abs(cl[1][it]->GetZ()-zmean[it])< padlength*0.5){
3085 best[0][it]=1;
4f1c04d3 3086 }
7ad19338 3087 }
3088 //
3089 // iterative choosing of "best path"
3090 //
3091 //
3092 Int_t label = TMath::Abs(track->GetLabel());
3093 Int_t bestiter=0;
3094 //
3095 for (Int_t iter=0;iter<9;iter++){
3096 //
3097 changes[iter]= 0;
3098 sumz = 0; sum=0; sumdy=0;sumdy2=0;sumx=0;sumx2=0;sumxy=0;mpads=0; ngood[iter]=0; nbad[iter]=0;
3099 // linear fit
3100 for (Int_t it=0;it<t1-t0;it++){
3101 if (!cl[best[iter][it]][it]) continue;
3102 //calculates pad-row changes
3103 Double_t zbefore= cl[best[iter][it]][it]->GetZ();
3104 Double_t zafter = cl[best[iter][it]][it]->GetZ();
3105 for (Int_t itd = it-1; itd>=0;itd--) {
3106 if (cl[best[iter][itd]][itd]) {
3107 zbefore= cl[best[iter][itd]][itd]->GetZ();
3108 break;
3109 }
3110 }
3111 for (Int_t itd = it+1; itd<t1-t0;itd++) {
3112 if (cl[best[iter][itd]][itd]) {
3113 zafter= cl[best[iter][itd]][itd]->GetZ();
3114 break;
3115 }
3116 }
3117 if (TMath::Abs(cl[best[iter][it]][it]->GetZ()-zbefore)>0.1&&TMath::Abs(cl[best[iter][it]][it]->GetZ()-zafter)>0.1) changes[iter]++;
3118 //
3119 Double_t dx = x[it]-xmean; // distance to reference x
3120 sumz += cl[best[iter][it]][it]->GetZ();
4f1c04d3 3121 sum++;
7ad19338 3122 sumdy += dy[best[iter][it]][it];
3123 sumdy2+= dy[best[iter][it]][it]*dy[best[iter][it]][it];
3124 sumx += dx;
3125 sumx2 += dx*dx;
3126 sumxy += dx*dy[best[iter][it]][it];
3127 mpads += cl[best[iter][it]][it]->GetNPads();
3128 if (cl[best[iter][it]][it]->GetLabel(0)==label || cl[best[iter][it]][it]->GetLabel(1)==label||cl[best[iter][it]][it]->GetLabel(2)==label){
3129 ngood[iter]++;
4f1c04d3 3130 }
3131 else{
7ad19338 3132 nbad[iter]++;
4f1c04d3 3133 }
3134 }
7ad19338 3135 //
3136 // calculates line parameters
3137 //
3138 Double_t det = sum*sumx2-sumx*sumx;
3139 angle[iter] = (sum*sumxy-sumx*sumdy)/det;
3140 mean[iter] = (sumx2*sumdy-sumx*sumxy)/det;
3141 meanz[iter] = sumz/sum;
3142 moffset[iter] = sumdy/sum;
3143 mpads /= sum; // mean number of pads
3144 //
3145 //
3146 Double_t sigma2 = 0; // normalized residuals - for line fit
3147 Double_t sigma1 = 0; // normalized residuals - constant fit
3148 //
3149 for (Int_t it=0;it<t1-t0;it++){
3150 if (!cl[best[iter][it]][it]) continue;
3151 Double_t dx = x[it]-xmean;
3152 Double_t ytr = mean[iter]+angle[iter]*dx;
3153 sigma2 += (dy[best[iter][it]][it]-ytr)*(dy[best[iter][it]][it]-ytr);
3154 sigma1 += (dy[best[iter][it]][it]-moffset[iter])*(dy[best[iter][it]][it]-moffset[iter]);
3155 sum++;
4f1c04d3 3156 }
7ad19338 3157 sigma2 /=(sum-2); // normalized residuals
3158 sigma1 /=(sum-1); // normalized residuals
3159 //
3160 smean[iter] = sigma2*(sumx2/det); // estimated error2 of mean
3161 sangle[iter] = sigma2*(sum/det); // estimated error2 of angle
3162 smeanangle[iter] = sigma2*(-sumx/det); // correlation
3163 //
3164 //
3165 sigmas[iter] = TMath::Sqrt(sigma1); //
3166 smoffset[iter]= (sigma1/sum)+0.01*0.01; // sigma of mean offset + unisochronity sigma
3167 //
3168 // iterative choosing of "better path"
3169 //
3170 for (Int_t it=0;it<t1-t0;it++){
3171 if (!cl[best[iter][it]][it]) continue;
3172 //
3173 Double_t sigmatr2 = smoffset[iter]+0.5*tany*tany; //add unisochronity + angular effect contribution
69b55c55 3174 Double_t sweight = 1./sigmatr2+1./track->GetSigmaY2();
7ad19338 3175 Double_t weighty = (moffset[iter]/sigmatr2)/sweight; // weighted mean
69b55c55 3176 Double_t sigmacl = TMath::Sqrt(sigma1*sigma1+track->GetSigmaY2()); //
7ad19338 3177 Double_t mindist=100000;
3178 Int_t ihbest=0;
3179 for (Int_t ih=0;ih<10;ih++){
3180 if (!cl[ih][it]) break;
3181 Double_t dist2 = (dy[ih][it]-weighty)/sigmacl;
3182 dist2*=dist2; //chi2 distance
3183 if (dist2<mindist){
3184 mindist = dist2;
3185 ihbest =ih;
3186 }
3187 }
3188 best[iter+1][it]=ihbest;
4f1c04d3 3189 }
4f1c04d3 3190 //
7ad19338 3191 // update best hypothesy if better chi2 according tracklet position and angle
3192 //
69b55c55 3193 Double_t sy2 = smean[iter] + track->GetSigmaY2();
7ad19338 3194 Double_t sa2 = sangle[iter] + track->fCee;
3195 Double_t say = track->fCey;
3196 // Double_t chi20 = mean[bestiter]*mean[bestiter]/sy2+angle[bestiter]*angle[bestiter]/sa2;
3197 // Double_t chi21 = mean[iter]*mean[iter]/sy2+angle[iter]*angle[iter]/sa2;
3198
3199 Double_t detchi = sy2*sa2-say*say;
3200 Double_t invers[3] = {sa2/detchi, sy2/detchi, -say/detchi}; //inverse value of covariance matrix
4f1c04d3 3201
7ad19338 3202 Double_t chi20 = mean[bestiter]*mean[bestiter]*invers[0]+angle[bestiter]*angle[bestiter]*invers[1]+
3203 2.*mean[bestiter]*angle[bestiter]*invers[2];
3204 Double_t chi21 = mean[iter]*mean[iter]*invers[0]+angle[iter]*angle[iter]*invers[1]+
3205 2*mean[iter]*angle[iter]*invers[2];
3206 tchi2s[iter] =chi21;
3207 //
3208 if (changes[iter]<=changes[bestiter] && chi21<chi20) {
3209 bestiter =iter;
3210 }
3211 }
3212 //
3213 //set clusters
3214 //
3215 Double_t sigma2 = sigmas[0]; // choose as sigma from 0 iteration
8979685e 3216 Short_t maxpos = -1;
3217 Float_t maxcharge = 0;
3218 Short_t maxpos4 = -1;
3219 Float_t maxcharge4 = 0;
3220 Short_t maxpos5 = -1;
3221 Float_t maxcharge5 = 0;
3222
7ad19338 3223 //if (tchi2s[bestiter]>25.) sigma2*=tchi2s[bestiter]/25.;
3224 //if (tchi2s[bestiter]>25.) sigma2=1000.; // dont'accept
3225
59393e34 3226 Double_t exB = AliTRDcalibDB::Instance()->GetOmegaTau(AliTRDcalibDB::Instance()->GetVdrift(0,0,0));
7ad19338 3227 Double_t expectederr = sigma2*sigma2+0.01*0.01;
3228 if (mpads>3.5) expectederr += (mpads-3.5)*0.04;
3229 if (changes[bestiter]>1) expectederr+= changes[bestiter]*0.01;
59393e34 3230 expectederr+=(0.03*(tany-exB)*(tany-exB))*15;
7ad19338 3231 // if (tchi2s[bestiter]>18.) expectederr*= tchi2s[bestiter]/18.;
3232 //expectederr+=10000;
3233 for (Int_t it=0;it<t1-t0;it++){
3234 if (!cl[best[bestiter][it]][it]) continue;
7ad19338 3235 cl[best[bestiter][it]][it]->SetSigmaY2(expectederr); // set cluster error
3236 if (!cl[best[bestiter][it]][it]->IsUsed()){
59393e34 3237 cl[best[bestiter][it]][it]->SetY( cl[best[bestiter][it]][it]->GetY());
69b55c55 3238 // cl[best[bestiter][it]][it]->Use();
3239 }
3240 //
3241 // time bins with maximal charge
3242 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge){
3243 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3244 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3245 }
3246
3247 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge4){
3248 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=4){
3249 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3250 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3251 }
3252 }
3253 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge5){
3254 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=5){
3255 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3256 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3257 }
7ad19338 3258 }
8979685e 3259 //
3260 // time bins with maximal charge
3261 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge){
3262 maxcharge = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3263 maxpos = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3264 }
3265
3266 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge4){
3267 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=4){
3268 maxcharge4 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3269 maxpos4 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3270 }
3271 }
3272 if (TMath::Abs(cl[best[bestiter][it]][it]->GetQ())> maxcharge5){
3273 if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>=5){
3274 maxcharge5 = TMath::Abs(cl[best[bestiter][it]][it]->GetQ());
3275 maxpos5 = cl[best[bestiter][it]][it]->GetLocalTimeBin();
3276 }
3277 }
7ad19338 3278 clusters[it+t0] = indexes[best[bestiter][it]][it];
69b55c55 3279 //if (cl[best[bestiter][it]][it]->GetLocalTimeBin()>4 && cl[best[bestiter][it]][it]->GetLocalTimeBin()<18) clusters[it+t0] = indexes[best[bestiter][it]][it]; //Test
7ad19338 3280 }
3281 //
3282 // set tracklet parameters
3283 //
3284 Double_t trackleterr2 = smoffset[bestiter]+0.01*0.01;
3285 if (mpads>3.5) trackleterr2 += (mpads-3.5)*0.04;
3286 trackleterr2+= changes[bestiter]*0.01;
3287 trackleterr2*= TMath::Max(14.-nfound,1.);
59393e34 3288 trackleterr2+= 0.2*(tany-exB)*(tany-exB);
7ad19338 3289 //
3290 tracklet.Set(xmean, track2.GetY()+moffset[bestiter], meanz[bestiter], track2.GetAlpha(), trackleterr2); //set tracklet parameters
3291 tracklet.SetTilt(h01);
3292 tracklet.SetP0(mean[bestiter]);
3293 tracklet.SetP1(angle[bestiter]);
3294 tracklet.SetN(nfound);
3295 tracklet.SetNCross(changes[bestiter]);
3296 tracklet.SetPlane(plane);
3297 tracklet.SetSigma2(expectederr);
3298 tracklet.SetChi2(tchi2s[bestiter]);
8979685e 3299 tracklet.SetMaxPos(maxpos,maxpos4,maxpos5);
7ad19338 3300 track->fTracklets[plane] = tracklet;
3301 track->fNWrong+=nbad[0];
3302 //
3303 // Debuging part
3304 //
69b55c55 3305 TClonesArray array0("AliTRDcluster");
3306 TClonesArray array1("AliTRDcluster");
3307 array0.ExpandCreateFast(t1-t0+1);
3308 array1.ExpandCreateFast(t1-t0+1);
7ad19338 3309 TTreeSRedirector& cstream = *fDebugStreamer;
3310 AliTRDcluster dummy;
3311 Double_t dy0[100];
8979685e 3312 Double_t dyb[100];
3313
7ad19338 3314 for (Int_t it=0;it<t1-t0;it++){
3315 dy0[it] = dy[0][it];
3316 dyb[it] = dy[best[bestiter][it]][it];
3317 if(cl[0][it]) {
3318 new(array0[it]) AliTRDcluster(*cl[0][it]);
3319 }
3320 else{
3321 new(array0[it]) AliTRDcluster(dummy);
3322 }
3323 if(cl[best[bestiter][it]][it]) {
3324 new(array1[it]) AliTRDcluster(*cl[best[bestiter][it]][it]);
3325 }
3326 else{
3327 new(array1[it]) AliTRDcluster(dummy);
3328 }
4f1c04d3 3329 }
7ad19338 3330 TGraph graph0(t1-t0,x,dy0);
3331 TGraph graph1(t1-t0,x,dyb);
3332 TGraph graphy(t1-t0,x,yt);
3333 TGraph graphz(t1-t0,x,zt);
3334 //
3335 //
d337ef8d 3336 if (AliTRDReconstructor::StreamLevel()>0)
7ad19338 3337 cstream<<"tracklet"<<
3338 "track.="<<track<< // track parameters
3339 "tany="<<tany<< // tangent of the local track angle
3340 "xmean="<<xmean<< // xmean - reference x of tracklet
3341 "tilt="<<h01<< // tilt angle
3342 "nall="<<nall<< // number of foundable clusters
3343 "nfound="<<nfound<< // number of found clusters
3344 "clfound="<<clfound<< // total number of found clusters in road
3345 "mpads="<<mpads<< // mean number of pads per cluster
3346 "plane="<<plane<< // plane number
c6f438c0 3347 "detector="<<detector<< // detector number
7ad19338 3348 "road="<<road<< // the width of the used road
3349 "graph0.="<<&graph0<< // x - y = dy for closest cluster
3350 "graph1.="<<&graph1<< // x - y = dy for second closest cluster
3351 "graphy.="<<&graphy<< // y position of the track
3352 "graphz.="<<&graphz<< // z position of the track
69b55c55 3353 // "fCl.="<<&array0<< // closest cluster
f6625211 3354 //"fCl2.="<<&array1<< // second closest cluster
8979685e 3355 "maxpos="<<maxpos<< // maximal charge postion
3356 "maxcharge="<<maxcharge<< // maximal charge
3357 "maxpos4="<<maxpos4<< // maximal charge postion - after bin 4
3358 "maxcharge4="<<maxcharge4<< // maximal charge - after bin 4
3359 "maxpos5="<<maxpos5<< // maximal charge postion - after bin 5
3360 "maxcharge5="<<maxcharge5<< // maximal charge - after bin 5
7ad19338 3361 //
3362 "bestiter="<<bestiter<< // best iteration number
3363 "tracklet.="<<&tracklet<< // corrspond to the best iteration
3364 "tchi20="<<tchi2s[0]<< // chi2 of cluster in the 0 iteration
3365 "tchi2b="<<tchi2s[bestiter]<< // chi2 of cluster in the best iteration
3366 "sigmas0="<<sigmas[0]<< // residuals sigma
3367 "sigmasb="<<sigmas[bestiter]<< // residulas sigma
3368 //
3369 "ngood0="<<ngood[0]<< // number of good clusters in 0 iteration
3370 "nbad0="<<nbad[0]<< // number of bad clusters in 0 iteration
3371 "ngoodb="<<ngood[bestiter]<< // in best iteration
3372 "nbadb="<<nbad[bestiter]<< // in best iteration
3373 //
3374 "changes0="<<changes[0]<< // changes of pardrows in iteration number 0
3375 "changesb="<<changes[bestiter]<< // changes of pardrows in best iteration
3376 //
3377 "moffset0="<<moffset[0]<< // offset fixing angle in iter=0
3378 "smoffset0="<<smoffset[0]<< // sigma of offset fixing angle in iter=0
3379 "moffsetb="<<moffset[bestiter]<< // offset fixing angle in iter=best
3380 "smoffsetb="<<smoffset[bestiter]<< // sigma of offset fixing angle in iter=best
3381 //
3382 "mean0="<<mean[0]<< // mean dy in iter=0;
3383 "smean0="<<smean[0]<< // sigma of mean dy in iter=0
3384 "meanb="<<mean[bestiter]<< // mean dy in iter=best
3385 "smeanb="<<smean[bestiter]<< // sigma of mean dy in iter=best
3386 //
3387 "angle0="<<angle[0]<< // angle deviation in the iteration number 0
3388 "sangle0="<<sangle[0]<< // sigma of angular deviation in iteration number 0
3389 "angleb="<<angle[bestiter]<< // angle deviation in the best iteration
3390 "sangleb="<<sangle[bestiter]<< // sigma of angle deviation in the best iteration
3391 //
3392 "expectederr="<<expectederr<< // expected error of cluster position
3393 "\n";
3394 //
3395 //
4f1c04d3 3396 return nfound;
4f1c04d3 3397
75bd7f81 3398}
4f1c04d3 3399
75bd7f81 3400//_____________________________________________________________________________
3401Int_t AliTRDtracker::Freq(Int_t n, const Int_t *inlist
3402 , Int_t *outlist, Bool_t down)
69b55c55 3403{
3404 //
3405 // Sort eleements according occurancy
3406 // The size of output array has is 2*n
3407 //
75bd7f81 3408
69b55c55 3409 Int_t * sindexS = new Int_t[n]; // temp array for sorting
3410 Int_t * sindexF = new Int_t[2*n];
3411 for (Int_t i=0;i<n;i++) sindexF[i]=0;
3412 //
3413 TMath::Sort(n,inlist, sindexS, down);
3414 Int_t last = inlist[sindexS[0]];
3415 Int_t val = last;
3416 sindexF[0] = 1;
3417 sindexF[0+n] = last;
3418 Int_t countPos = 0;
3419 //
3420 // find frequency
3421 for(Int_t i=1;i<n; i++){
3422 val = inlist[sindexS[i]];
3423 if (last == val) sindexF[countPos]++;
3424 else{
3425 countPos++;
3426 sindexF[countPos+n] = val;
3427 sindexF[countPos]++;
3428 last =val;
3429 }
3430 }
3431 if (last==val) countPos++;
3432 // sort according frequency
3433 TMath::Sort(countPos, sindexF, sindexS, kTRUE);
3434 for (Int_t i=0;i<countPos;i++){
3435 outlist[2*i ] = sindexF[sindexS[i]+n];
3436 outlist[2*i+1] = sindexF[sindexS[i]];
3437 }
3438 delete [] sindexS;
3439 delete [] sindexF;
3440
3441 return countPos;
75bd7f81 3442
69b55c55 3443}
3444
75bd7f81 3445//_____________________________________________________________________________
3446AliTRDtrack *AliTRDtracker::RegisterSeed(AliTRDseed * seeds, Double_t * params)
69b55c55 3447{
3448 //
75bd7f81 3449 // Register a seed
69b55c55 3450 //
3451 Double_t alpha=AliTRDgeometry::GetAlpha();
3452 Double_t shift=AliTRDgeometry::GetAlpha()/2.;
3453 Double_t c[15];
3454 c[0] = 0.2;
3455 c[1] = 0 ; c[2] = 2;
3456 c[3] = 0 ; c[4] = 0; c[5] = 0.02;
3457 c[6] = 0 ; c[7] = 0; c[8] = 0; c[9] = 0.1;
3458 c[10] = 0 ; c[11] = 0; c[12] = 0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
3459 //
3460 Int_t index =0;
3461 AliTRDcluster *cl =0;
3462 for (Int_t ilayer=0;ilayer<6;ilayer++){
c6f438c0 3463 if (seeds[ilayer].IsOK()){
69b55c55 3464 for (Int_t itime=22;itime>0;itime--){
3465 if (seeds[ilayer].fIndexes[itime]>0){
3466 index = seeds[ilayer].fIndexes[itime];
3467 cl = seeds[ilayer].fClusters[itime];
3468 break;
3469 }
3470 }
3471 }
3472 if (index>0) break;
3473 }
3474 if (cl==0) return 0;
3475 AliTRDtrack * track = new AliTRDtrack(cl,index,&params[1],c, params[0],params[6]*alpha+shift);
3476 track->PropagateTo(params[0]-5.);
3477 track->ResetCovariance(1);
3478 //
f6625211 3479 Int_t rc=FollowBackProlongation(*track);
69b55c55 3480 if (rc<30) {
3481 delete track;
3482 track =0;
3483 }else{
3484 track->CookdEdx();
3485 CookdEdxTimBin(*track);
3486 CookLabel(track, 0.9);
3487 }
69b55c55 3488
75bd7f81 3489 return track;
69b55c55 3490
69b55c55 3491}