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