1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ///////////////////////////////////////////////////////////////////////////////
23 // Alex Bercuci <A.Bercuci@gsi.de> //
24 // Markus Fasel <M.Fasel@gsi.de> //
26 ///////////////////////////////////////////////////////////////////////////////
28 // #include <Riostream.h>
30 // #include <string.h>
33 #include <TDirectory.h>
34 #include <TLinearFitter.h>
36 #include <TClonesArray.h>
37 #include <TTreeStream.h>
38 #include <TGeoMatrix.h>
39 #include <TGeoManager.h>
42 #include "AliMathBase.h"
43 #include "AliESDEvent.h"
44 #include "AliGeomManager.h"
45 #include "AliRieman.h"
46 #include "AliTrackPointArray.h"
48 #include "AliTRDgeometry.h"
49 #include "AliTRDpadPlane.h"
50 #include "AliTRDcalibDB.h"
51 #include "AliTRDReconstructor.h"
52 #include "AliTRDCalibraFillHisto.h"
53 #include "AliTRDrecoParam.h"
55 #include "AliTRDcluster.h"
56 #include "AliTRDseedV1.h"
57 #include "AliTRDtrackV1.h"
58 #include "AliTRDtrackerV1.h"
59 #include "AliTRDtrackerDebug.h"
60 #include "AliTRDtrackingChamber.h"
61 #include "AliTRDchamberTimeBin.h"
65 ClassImp(AliTRDtrackerV1)
68 const Float_t AliTRDtrackerV1::fgkMinClustersInTrack = 0.5; //
69 const Float_t AliTRDtrackerV1::fgkLabelFraction = 0.8; //
70 const Double_t AliTRDtrackerV1::fgkMaxChi2 = 12.0; //
71 const Double_t AliTRDtrackerV1::fgkMaxSnp = 0.95; // Maximum local sine of the azimuthal angle
72 const Double_t AliTRDtrackerV1::fgkMaxStep = 2.0; // Maximal step size in propagation
73 Double_t AliTRDtrackerV1::fgTopologicQA[kNConfigs] = {
74 0.1112, 0.1112, 0.1112, 0.0786, 0.0786,
75 0.0786, 0.0786, 0.0579, 0.0579, 0.0474,
76 0.0474, 0.0408, 0.0335, 0.0335, 0.0335
78 const Double_t AliTRDtrackerV1::fgkX0[kNPlanes] = {
79 300.2, 312.8, 325.4, 338.0, 350.6, 363.2};
80 Int_t AliTRDtrackerV1::fgNTimeBins = 0;
81 AliRieman* AliTRDtrackerV1::fgRieman = 0x0;
82 TLinearFitter* AliTRDtrackerV1::fgTiltedRieman = 0x0;
83 TLinearFitter* AliTRDtrackerV1::fgTiltedRiemanConstrained = 0x0;
85 //____________________________________________________________________
86 AliTRDtrackerV1::AliTRDtrackerV1(AliTRDReconstructor *rec)
96 // Default constructor.
99 SetReconstructor(rec); // initialize reconstructor
101 // initialize geometry
102 if(!AliGeomManager::GetGeometry()){
103 AliFatal("Could not get geometry.");
105 fGeom = new AliTRDgeometry();
106 fGeom->CreateClusterMatrixArray();
107 TGeoHMatrix *matrix = 0x0;
108 Double_t loc[] = {0., 0., 0.};
109 Double_t glb[] = {0., 0., 0.};
110 for(Int_t ily=kNPlanes; ily--;){
112 while(!(matrix = fGeom->GetClusterMatrix(AliTRDgeometry::GetDetector(ily, 2, ism)))) ism++;
114 AliError(Form("Could not get transformation matrix for layer %d. Use default.", ily));
115 fR[ily] = fgkX0[ily];
118 matrix->LocalToMaster(loc, glb);
119 fR[ily] = glb[0]+ AliTRDgeometry::AnodePos()-.5*AliTRDgeometry::AmThick() - AliTRDgeometry::DrThick();
122 // initialize calibration values
123 AliTRDcalibDB *trd = 0x0;
124 if (!(trd = AliTRDcalibDB::Instance())) {
125 AliFatal("Could not get calibration.");
127 if(!fgNTimeBins) fgNTimeBins = trd->GetNumberOfTimeBins();
129 // initialize cluster containers
130 for (Int_t isector = 0; isector < AliTRDgeometry::kNsector; isector++) new(&fTrSec[isector]) AliTRDtrackingSector(fGeom, isector);
133 memset(fTrackQuality, 0, kMaxTracksStack*sizeof(Double_t));
134 memset(fSeedLayer, 0, kMaxTracksStack*sizeof(Int_t));
135 memset(fSeedTB, 0, kNSeedPlanes*sizeof(AliTRDchamberTimeBin*));
138 //____________________________________________________________________
139 AliTRDtrackerV1::~AliTRDtrackerV1()
145 if(fgRieman) delete fgRieman; fgRieman = 0x0;
146 if(fgTiltedRieman) delete fgTiltedRieman; fgTiltedRieman = 0x0;
147 if(fgTiltedRiemanConstrained) delete fgTiltedRiemanConstrained; fgTiltedRiemanConstrained = 0x0;
148 for(Int_t isl =0; isl<kNSeedPlanes; isl++) if(fSeedTB[isl]) delete fSeedTB[isl];
149 if(fTracks) {fTracks->Delete(); delete fTracks;}
150 if(fTracklets) {fTracklets->Delete(); delete fTracklets;}
152 fClusters->Delete(); delete fClusters;
154 if(fGeom) delete fGeom;
157 //____________________________________________________________________
158 Int_t AliTRDtrackerV1::Clusters2Tracks(AliESDEvent *esd)
161 // Steering stand alone tracking for full TRD detector
164 // esd : The ESD event. On output it contains
165 // the ESD tracks found in TRD.
168 // Number of tracks found in the TRD detector.
170 // Detailed description
171 // 1. Launch individual SM trackers.
172 // See AliTRDtrackerV1::Clusters2TracksSM() for details.
175 if(!fReconstructor->GetRecoParam() ){
176 AliError("Reconstruction configuration not initialized. Call first AliTRDReconstructor::SetRecoParam().");
180 //AliInfo("Start Track Finder ...");
182 for(int ism=0; ism<AliTRDgeometry::kNsector; ism++){
183 // for(int ism=1; ism<2; ism++){
184 //AliInfo(Form("Processing supermodule %i ...", ism));
185 ntracks += Clusters2TracksSM(ism, esd);
187 AliInfo(Form("Number of found tracks : %d", ntracks));
192 //_____________________________________________________________________________
193 Bool_t AliTRDtrackerV1::GetTrackPoint(Int_t index, AliTrackPoint &p) const
195 //AliInfo(Form("Asking for tracklet %d", index));
197 // reset position of the point before using it
198 p.SetXYZ(0., 0., 0.);
199 AliTRDseedV1 *tracklet = GetTracklet(index);
200 if (!tracklet) return kFALSE;
202 // get detector for this tracklet
203 Int_t det = tracklet->GetDetector();
204 Int_t sec = fGeom->GetSector(det);
205 Double_t alpha = (sec+.5)*AliTRDgeometry::GetAlpha(),
206 sinA = TMath::Sin(alpha),
207 cosA = TMath::Cos(alpha);
209 local[0] = tracklet->GetX();
210 local[1] = tracklet->GetY();
211 local[2] = tracklet->GetZ();
213 fGeom->RotateBack(det, local, global);
215 Double_t cov2D[3]; Float_t cov[6];
216 tracklet->GetCovAt(local[0], cov2D);
217 cov[0] = cov2D[0]*sinA*sinA;
218 cov[1] =-cov2D[0]*sinA*cosA;
219 cov[2] =-cov2D[1]*sinA;
220 cov[3] = cov2D[0]*cosA*cosA;
221 cov[4] = cov2D[1]*cosA;
223 // store the global position of the tracklet and its covariance matrix in the track point
224 p.SetXYZ(global[0],global[1],global[2], cov);
227 AliGeomManager::ELayerID iLayer = AliGeomManager::ELayerID(AliGeomManager::kTRD1+fGeom->GetLayer(det));
228 Int_t modId = fGeom->GetSector(det) * AliTRDgeometry::kNstack + fGeom->GetStack(det);
229 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer, modId);
230 p.SetVolumeID(volid);
235 //____________________________________________________________________
236 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitter()
238 if(!fgTiltedRieman) fgTiltedRieman = new TLinearFitter(4, "hyp4");
239 return fgTiltedRieman;
242 //____________________________________________________________________
243 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitterConstraint()
245 if(!fgTiltedRiemanConstrained) fgTiltedRiemanConstrained = new TLinearFitter(2, "hyp2");
246 return fgTiltedRiemanConstrained;
249 //____________________________________________________________________
250 AliRieman* AliTRDtrackerV1::GetRiemanFitter()
252 if(!fgRieman) fgRieman = new AliRieman(AliTRDseedV1::kNtb * AliTRDgeometry::kNlayer);
256 //_____________________________________________________________________________
257 Int_t AliTRDtrackerV1::PropagateBack(AliESDEvent *event)
260 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
261 // backpropagated by the TPC tracker. Each seed is first propagated
262 // to the TRD, and then its prolongation is searched in the TRD.
263 // If sufficiently long continuation of the track is found in the TRD
264 // the track is updated, otherwise it's stored as originaly defined
265 // by the TPC tracker.
268 // Calibration monitor
269 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
270 if (!calibra) AliInfo("Could not get Calibra instance\n");
273 Int_t nFound = 0, // number of tracks found
274 nSeeds = 0, // total number of ESD seeds
275 nTRDseeds= 0, // number of seeds in the TRD acceptance
276 nTPCseeds= 0; // number of TPC seeds
277 Float_t foundMin = 20.0;
279 Float_t *quality = 0x0;
281 nSeeds = event->GetNumberOfTracks();
282 // Sort tracks according to quality
283 // (covariance in the yz plane)
285 quality = new Float_t[nSeeds];
286 index = new Int_t[nSeeds];
287 for (Int_t iSeed = nSeeds; iSeed--;) {
288 AliESDtrack *seed = event->GetTrack(iSeed);
289 Double_t covariance[15];
290 seed->GetExternalCovariance(covariance);
291 quality[iSeed] = covariance[0] + covariance[2];
293 TMath::Sort(nSeeds, quality, index,kFALSE);
296 // Propagate all seeds
299 for (Int_t iSeed = 0; iSeed < nSeeds; iSeed++) {
301 // Get the seeds in sorted sequence
302 AliESDtrack *seed = event->GetTrack(index[iSeed]);
303 Float_t p4 = seed->GetC(seed->GetBz());
305 // Check the seed status
306 ULong_t status = seed->GetStatus();
307 if ((status & AliESDtrack::kTPCout) == 0) continue;
308 if ((status & AliESDtrack::kTRDout) != 0) continue;
310 // Propagate to the entrance in the TRD mother volume
311 new(&track) AliTRDtrackV1(*seed);
312 if(AliTRDgeometry::GetXtrdBeg() > (fgkMaxStep + track.GetX()) && !PropagateToX(track, AliTRDgeometry::GetXtrdBeg(), fgkMaxStep)){
313 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
316 if(!AdjustSector(&track)){
317 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
320 if(TMath::Abs(track.GetSnp()) > fgkMaxSnp) {
321 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
327 // store track status at TRD entrance
328 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
330 // prepare track and do propagation in the TRD
331 track.SetReconstructor(fReconstructor);
332 track.SetKink(Bool_t(seed->GetKinkIndex(0)));
333 expectedClr = FollowBackProlongation(track);
334 // check if track entered the TRD fiducial volume
335 if(track.GetTrackLow()){
336 seed->UpdateTrackParams(&track, AliESDtrack::kTRDin);
339 // check if track was stopped in the TRD
341 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
347 // computes PID for track
349 // update calibration references using this track
350 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(&track);
351 // save calibration object
352 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0){
353 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(track);
354 calibTrack->SetOwner();
355 seed->AddCalibObject(calibTrack);
358 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
359 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
360 track.UpdateESDtrack(seed);
364 if ((TMath::Abs(track.GetC(track.GetBz()) - p4) / TMath::Abs(p4) < 0.2) ||(track.Pt() > 0.8)) {
366 // Make backup for back propagation
367 Int_t foundClr = track.GetNumberOfClusters();
368 if (foundClr >= foundMin) {
369 track.CookLabel(1. - fgkLabelFraction);
370 //if(track.GetBackupTrack()) UseClusters(track.GetBackupTrack());
372 // Sign only gold tracks
373 if (track.GetChi2() / track.GetNumberOfClusters() < 4) {
374 //if ((seed->GetKinkIndex(0) == 0) && (track.Pt() < 1.5)) UseClusters(&track);
376 Bool_t isGold = kFALSE;
379 if (track.GetChi2() / track.GetNumberOfClusters() < 5) {
380 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
386 if ((!isGold) && (track.GetNCross() == 0) && (track.GetChi2() / track.GetNumberOfClusters() < 7)) {
387 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
388 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
393 if ((!isGold) && (track.GetBackupTrack())) {
394 if ((track.GetBackupTrack()->GetNumberOfClusters() > foundMin) && ((track.GetBackupTrack()->GetChi2()/(track.GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
395 seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
402 // Propagation to the TOF
403 if(!(seed->GetStatus()&AliESDtrack::kTRDStop)) {
404 Int_t sm = track.GetSector();
405 // default value in case we have problems with the geometry.
406 Double_t xtof = 371.;
407 //Calculate radial position of the beginning of the TOF
408 //mother volume. In order to avoid mixing of the TRD
409 //and TOF modules some hard values are needed. This are:
410 //1. The path to the TOF module.
411 //2. The width of the TOF (29.05 cm)
412 //(with the help of Annalisa de Caro Mar-17-2009)
414 gGeoManager->cd(Form("/ALIC_1/B077_1/BSEGMO%d_1/BTOF%d_1", sm, sm));
415 TGeoHMatrix *m = 0x0;
416 Double_t loc[]={0., 0., -.5*29.05}, glob[3];
418 if((m=gGeoManager->GetCurrentMatrix())){
419 m->LocalToMaster(loc, glob);
420 xtof = TMath::Sqrt(glob[0]*glob[0]+glob[1]*glob[1]);
423 if(xtof > (fgkMaxStep + track.GetX()) && !PropagateToX(track, xtof, fgkMaxStep)){
424 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
427 if(!AdjustSector(&track)){
428 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
431 if(TMath::Abs(track.GetSnp()) > fgkMaxSnp){
432 seed->UpdateTrackParams(&track, AliESDtrack::kTRDStop);
435 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
436 // TODO obsolete - delete
437 seed->SetTRDQuality(track.StatusForTOF());
439 seed->SetTRDBudget(track.GetBudget(0));
441 if(index) delete [] index;
442 if(quality) delete [] quality;
445 AliInfo(Form("Number of TPC seeds: %d (%d)", nTRDseeds, nTPCseeds));
446 AliInfo(Form("Number of propagated TRD tracks: %d", nFound));
448 // run stand alone tracking
449 if (fReconstructor->IsSeeding()) Clusters2Tracks(event);
455 //____________________________________________________________________
456 Int_t AliTRDtrackerV1::RefitInward(AliESDEvent *event)
459 // Refits tracks within the TRD. The ESD event is expected to contain seeds
460 // at the outer part of the TRD.
461 // The tracks are propagated to the innermost time bin
462 // of the TRD and the ESD event is updated
463 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
466 Int_t nseed = 0; // contor for loaded seeds
467 Int_t found = 0; // contor for updated TRD tracks
471 for (Int_t itrack = 0; itrack < event->GetNumberOfTracks(); itrack++) {
472 AliESDtrack *seed = event->GetTrack(itrack);
473 new(&track) AliTRDtrackV1(*seed);
475 if (track.GetX() < 270.0) {
476 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
480 // reject tracks which failed propagation in the TRD or
481 // are produced by the TRD stand alone tracker
482 ULong_t status = seed->GetStatus();
483 if(!(status & AliESDtrack::kTRDout)) continue;
484 if(!(status & AliESDtrack::kTRDin)) continue;
487 track.ResetCovariance(50.0);
489 // do the propagation and processing
490 Bool_t kUPDATE = kFALSE;
491 Double_t xTPC = 250.0;
492 if(FollowProlongation(track)){
494 if (PropagateToX(track, xTPC, fgkMaxStep)) { // -with update
495 seed->UpdateTrackParams(&track, AliESDtrack::kTRDrefit);
500 // Update the friend track
501 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0){
502 TObject *o = 0x0; Int_t ic = 0;
503 AliTRDtrackV1 *calibTrack = 0x0;
504 while((o = seed->GetCalibObject(ic++))){
505 if(!(calibTrack = dynamic_cast<AliTRDtrackV1*>(o))) continue;
506 calibTrack->SetTrackHigh(track.GetTrackHigh());
511 // Prolongate to TPC without update
513 AliTRDtrackV1 tt(*seed);
514 if (PropagateToX(tt, xTPC, fgkMaxStep)) seed->UpdateTrackParams(&tt, AliESDtrack::kTRDbackup);
517 AliInfo(Form("Number of loaded seeds: %d",nseed));
518 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
523 //____________________________________________________________________
524 Int_t AliTRDtrackerV1::FollowProlongation(AliTRDtrackV1 &t)
526 // Extrapolates the TRD track in the TPC direction.
529 // t : the TRD track which has to be extrapolated
532 // number of clusters attached to the track
534 // Detailed description
536 // Starting from current radial position of track <t> this function
537 // extrapolates the track through the 6 TRD layers. The following steps
538 // are being performed for each plane:
540 // a. get plane limits in the local x direction
541 // b. check crossing sectors
542 // c. check track inclination
543 // 2. search tracklet in the tracker list (see GetTracklet() for details)
544 // 3. evaluate material budget using the geo manager
545 // 4. propagate and update track using the tracklet information.
550 Bool_t kStoreIn = kTRUE;
551 Int_t nClustersExpected = 0;
552 for (Int_t iplane = kNPlanes; iplane--;) {
554 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
555 if(!tracklet) continue;
556 if(!tracklet->IsOK()) AliWarning("tracklet not OK");
558 Double_t x = tracklet->GetX();//GetX0();
559 // reject tracklets which are not considered for inward refit
560 if(x > t.GetX()+fgkMaxStep) continue;
562 // append tracklet to track
563 t.SetTracklet(tracklet, index);
565 if (x < (t.GetX()-fgkMaxStep) && !PropagateToX(t, x+fgkMaxStep, fgkMaxStep)) break;
566 if (!AdjustSector(&t)) break;
568 // Start global position
572 // End global position
573 Double_t alpha = t.GetAlpha(), y, z;
574 if (!t.GetProlongation(x,y,z)) break;
576 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
577 xyz1[1] = x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
580 Double_t length = TMath::Sqrt(
581 (xyz0[0]-xyz1[0])*(xyz0[0]-xyz1[0]) +
582 (xyz0[1]-xyz1[1])*(xyz0[1]-xyz1[1]) +
583 (xyz0[2]-xyz1[2])*(xyz0[2]-xyz1[2])
586 // Get material budget
588 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) break;
589 Double_t xrho= param[0]*param[4];
590 Double_t xx0 = param[1]; // Get mean propagation parameters
592 // Propagate and update
593 t.PropagateTo(x, xx0, xrho);
594 if (!AdjustSector(&t)) break;
601 Double_t maxChi2 = t.GetPredictedChi2(tracklet);
602 if (maxChi2 < 1e+10 && t.Update(tracklet, maxChi2)){
603 nClustersExpected += tracklet->GetN();
607 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
609 for(int iplane=0; iplane<AliTRDgeometry::kNlayer; iplane++){
610 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
611 if(!tracklet) continue;
612 t.SetTracklet(tracklet, index);
615 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
616 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
617 AliTRDtrackV1 track(t);
619 cstreamer << "FollowProlongation"
620 << "EventNumber=" << eventNumber
621 << "ncl=" << nClustersExpected
622 << "track.=" << &track
626 return nClustersExpected;
630 //_____________________________________________________________________________
631 Int_t AliTRDtrackerV1::FollowBackProlongation(AliTRDtrackV1 &t)
633 // Extrapolates the TRD track in the TOF direction.
636 // t : the TRD track which has to be extrapolated
639 // number of clusters attached to the track
641 // Detailed description
643 // Starting from current radial position of track <t> this function
644 // extrapolates the track through the 6 TRD layers. The following steps
645 // are being performed for each plane:
647 // a. get plane limits in the local x direction
648 // b. check crossing sectors
649 // c. check track inclination
650 // 2. build tracklet (see AliTRDseed::AttachClusters() for details)
651 // 3. evaluate material budget using the geo manager
652 // 4. propagate and update track using the tracklet information.
658 Double_t driftLength = .5*AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick();
659 AliTRDtrackingChamber *chamber = 0x0;
661 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
662 // in case of stand alone tracking we store all the pointers to the tracklets in a temporary array
663 AliTRDseedV1 *tracklets[kNPlanes];
664 memset(tracklets, 0, sizeof(AliTRDseedV1 *) * kNPlanes);
665 for(Int_t ip = 0; ip < kNPlanes; ip++){
666 tracklets[ip] = t.GetTracklet(ip);
669 Bool_t kStoreIn = kTRUE;
671 // Loop through the TRD layers
672 TGeoHMatrix *matrix = 0x0;
674 for (Int_t ily=0, sm=-1, stk=-1, det=-1; ily < AliTRDgeometry::kNlayer; ily++) {
675 // rough estimate of the entry point
676 if (!t.GetProlongation(fR[ily], y, z)){
678 t.SetStatus(AliTRDtrackV1::kProlongation);
682 // find sector / stack / detector
684 // TODO cross check with y value !
685 stk = fGeom->GetStack(z, ily);
686 det = stk>=0 ? AliTRDgeometry::GetDetector(ily, stk, sm) : -1;
687 matrix = det>=0 ? fGeom->GetClusterMatrix(det) : 0x0;
689 // check if supermodule/chamber is installed
690 if( !fGeom->GetSMstatus(sm) ||
692 fGeom->IsHole(ily, stk, sm) ||
694 // propagate to the default radial position
695 if(fR[ily] > (fgkMaxStep + t.GetX()) && !PropagateToX(t, fR[ily], fgkMaxStep)){
697 t.SetStatus(AliTRDtrackV1::kPropagation);
700 if(!AdjustSector(&t)){
702 t.SetStatus(AliTRDtrackV1::kAdjustSector);
705 if(TMath::Abs(t.GetSnp()) > fgkMaxSnp){
707 t.SetStatus(AliTRDtrackV1::kSnp);
710 t.SetStatus(AliTRDtrackV1::kGeometry, ily);
714 // retrieve rotation matrix for the current chamber
715 Double_t loc[] = {AliTRDgeometry::AnodePos()- driftLength, 0., 0.};
716 Double_t glb[] = {0., 0., 0.};
717 matrix->LocalToMaster(loc, glb);
719 // Propagate to the radial distance of the current layer
720 x = glb[0] - fgkMaxStep;
721 if(x > (fgkMaxStep + t.GetX()) && !PropagateToX(t, x, fgkMaxStep)){
723 t.SetStatus(AliTRDtrackV1::kPropagation);
726 if(!AdjustSector(&t)){
728 t.SetStatus(AliTRDtrackV1::kAdjustSector);
731 if(TMath::Abs(t.GetSnp()) > fgkMaxSnp) {
733 t.SetStatus(AliTRDtrackV1::kSnp);
736 Bool_t RECALCULATE = kFALSE;
737 if(sm != t.GetSector()){
741 if(stk != fGeom->GetStack(z, ily)){
742 stk = fGeom->GetStack(z, ily);
746 det = AliTRDgeometry::GetDetector(ily, stk, sm);
747 if(!(matrix = fGeom->GetClusterMatrix(det))){
748 t.SetStatus(AliTRDtrackV1::kGeometry, ily);
751 matrix->LocalToMaster(loc, glb);
752 x = glb[0] - fgkMaxStep;
755 // check if track is well inside fiducial volume
756 if (!t.GetProlongation(x+fgkMaxStep, y, z)) {
758 t.SetStatus(AliTRDtrackV1::kProlongation);
761 if(fGeom->IsOnBoundary(det, y, z, .5)){
762 t.SetStatus(AliTRDtrackV1::kBoundary, ily);
765 // mark track as entering the FIDUCIAL volume of TRD
771 ptrTracklet = tracklets[ily];
772 if(!ptrTracklet){ // BUILD TRACKLET
773 // check data in supermodule
774 if(!fTrSec[sm].GetNChambers()){
775 t.SetStatus(AliTRDtrackV1::kNoClusters, ily);
778 if(fTrSec[sm].GetX(ily) < 1.){
779 t.SetStatus(AliTRDtrackV1::kNoClusters, ily);
783 // check data in chamber
784 if(!(chamber = fTrSec[sm].GetChamber(stk, ily))){
785 t.SetStatus(AliTRDtrackV1::kNoClusters, ily);
788 if(chamber->GetNClusters() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()){
789 t.SetStatus(AliTRDtrackV1::kNoClusters, ily);
793 ptrTracklet = new(&tracklet) AliTRDseedV1(det);
794 ptrTracklet->SetReconstructor(fReconstructor);
795 ptrTracklet->SetKink(t.IsKink());
796 ptrTracklet->SetPadPlane(fGeom->GetPadPlane(ily, stk));
797 ptrTracklet->SetX0(glb[0]+driftLength);
798 if(!tracklet.Init(&t)){
800 t.SetStatus(AliTRDtrackV1::kTrackletInit);
803 if(!tracklet.AttachClusters(chamber, kTRUE)){
804 t.SetStatus(AliTRDtrackV1::kNoAttach, ily);
807 if(tracklet.GetN() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()){
808 t.SetStatus(AliTRDtrackV1::kNoClustersTracklet, ily);
811 ptrTracklet->UpdateUsed();
814 // propagate track to the radial position of the tracklet
815 ptrTracklet->UseClusters(); // TODO ? do we need this here ?
816 // fit tracklet no tilt correction
817 if(!ptrTracklet->Fit(kFALSE)){
818 t.SetStatus(AliTRDtrackV1::kNoFit, ily);
821 x = ptrTracklet->GetX(); //GetX0();
822 if(x > (fgkMaxStep + t.GetX()) && !PropagateToX(t, x, fgkMaxStep)) {
824 t.SetStatus(AliTRDtrackV1::kPropagation);
827 if(!AdjustSector(&t)) {
829 t.SetStatus(AliTRDtrackV1::kAdjustSector);
832 if(TMath::Abs(t.GetSnp()) > fgkMaxSnp) {
834 t.SetStatus(AliTRDtrackV1::kSnp);
838 // update Kalman with the TRD measurement
839 Double_t chi2 = t.GetPredictedChi2(ptrTracklet);
841 t.SetStatus(AliTRDtrackV1::kChi2, ily);
844 if(!t.Update(ptrTracklet, chi2)) {
846 t.SetStatus(AliTRDtrackV1::kUpdate);
850 // load tracklet to the tracker
851 ptrTracklet->UpDate(&t);
852 ptrTracklet = SetTracklet(ptrTracklet);
853 t.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
854 n += ptrTracklet->GetN();
856 // Reset material budget if 2 consecutive gold
857 // if(ilayer>0 && t.GetTracklet(ilayer-1) && ptrTracklet->GetN() + t.GetTracklet(ilayer-1)->GetN() > 20) t.SetBudget(2, 0.);
859 // Make backup of the track until is gold
860 // TO DO update quality check of the track.
861 // consider comparison with fTimeBinsRange
862 Float_t ratio0 = ptrTracklet->GetN() / Float_t(fgNTimeBins);
863 //Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);
868 //(ratio0+ratio1 > 1.5) &&
869 (t.GetNCross() == 0) &&
870 (TMath::Abs(t.GetSnp()) < 0.85) &&
871 (t.GetNumberOfClusters() > 20)){
875 //printf("clusters[%d] chi2[%f] x[%f] status[%d ", n, t.GetChi2(), t.GetX(), t.GetStatusTRD());
876 //for(int i=0; i<6; i++) printf("%d ", t.GetStatusTRD(i)); printf("]\n");
878 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
879 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
880 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
881 AliTRDtrackV1 track(t);
883 cstreamer << "FollowBackProlongation"
884 << "EventNumber=" << eventNumber
886 << "track.=" << &track
893 //_________________________________________________________________________
894 Float_t AliTRDtrackerV1::FitRieman(AliTRDseedV1 *tracklets, Double_t *chi2, Int_t *planes){
896 // Fits a Riemann-circle to the given points without tilting pad correction.
897 // The fit is performed using an instance of the class AliRieman (equations
898 // and transformations see documentation of this class)
899 // Afterwards all the tracklets are Updated
901 // Parameters: - Array of tracklets (AliTRDseedV1)
902 // - Storage for the chi2 values (beginning with direction z)
903 // - Seeding configuration
904 // Output: - The curvature
906 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
908 Int_t allplanes[] = {0, 1, 2, 3, 4, 5};
909 Int_t *ppl = &allplanes[0];
915 for(Int_t il = 0; il < maxLayers; il++){
916 if(!tracklets[ppl[il]].IsOK()) continue;
917 fitter->AddPoint(tracklets[ppl[il]].GetX0(), tracklets[ppl[il]].GetYfit(0), tracklets[ppl[il]].GetZfit(0),1,10);
920 // Set the reference position of the fit and calculate the chi2 values
921 memset(chi2, 0, sizeof(Double_t) * 2);
922 for(Int_t il = 0; il < maxLayers; il++){
923 // Reference positions
924 tracklets[ppl[il]].Init(fitter);
927 if((!tracklets[ppl[il]].IsOK()) && (!planes)) continue;
928 chi2[0] += tracklets[ppl[il]].GetChi2Y();
929 chi2[1] += tracklets[ppl[il]].GetChi2Z();
931 return fitter->GetC();
934 //_________________________________________________________________________
935 void AliTRDtrackerV1::FitRieman(AliTRDcluster **seedcl, Double_t chi2[2])
938 // Performs a Riemann helix fit using the seedclusters as spacepoints
939 // Afterwards the chi2 values are calculated and the seeds are updated
941 // Parameters: - The four seedclusters
942 // - The tracklet array (AliTRDseedV1)
943 // - The seeding configuration
948 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
951 for(Int_t i = 0; i < 4; i++){
952 printf(" {%f, %f, %f}\n", seedcl[i]->GetX(), seedcl[i]->GetY(), seedcl[i]->GetZ());
953 fitter->AddPoint(seedcl[i]->GetX(), seedcl[i]->GetY(), seedcl[i]->GetZ(), 1., 10.);
958 // Update the seed and calculated the chi2 value
959 chi2[0] = 0; chi2[1] = 0;
960 for(Int_t ipl = 0; ipl < kNSeedPlanes; ipl++){
962 chi2[0] += (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX()));
963 chi2[1] += (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX()));
968 //_________________________________________________________________________
969 Float_t AliTRDtrackerV1::FitTiltedRiemanConstraint(AliTRDseedV1 *tracklets, Double_t zVertex)
972 // Fits a helix to the clusters. Pad tilting is considered. As constraint it is
973 // assumed that the vertex position is set to 0.
974 // This method is very usefull for high-pt particles
975 // Basis for the fit: (x - x0)^2 + (y - y0)^2 - R^2 = 0
976 // x0, y0: Center of the circle
977 // Measured y-position: ymeas = y - tan(phiT)(zc - zt)
978 // zc: center of the pad row
979 // Equation which has to be fitted (after transformation):
980 // a + b * u + e * v + 2*(ymeas + tan(phiT)(z - zVertex))*t = 0
984 // v = 2 * x * tan(phiT) * t
985 // Parameters in the equation:
986 // a = -1/y0, b = x0/y0, e = dz/dx
988 // The Curvature is calculated by the following equation:
989 // - curv = a/Sqrt(b^2 + 1) = 1/R
990 // Parameters: - the 6 tracklets
991 // - the Vertex constraint
992 // Output: - the Chi2 value of the track
997 TLinearFitter *fitter = GetTiltedRiemanFitterConstraint();
998 fitter->StoreData(kTRUE);
999 fitter->ClearPoints();
1000 AliTRDcluster *cl = 0x0;
1002 Float_t x, y, z, w, t, error, tilt;
1005 for(Int_t ilr = 0; ilr < AliTRDgeometry::kNlayer; ilr++){
1006 if(!tracklets[ilr].IsOK()) continue;
1007 for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
1008 if(!tracklets[ilr].IsUsable(itb)) continue;
1009 cl = tracklets[ilr].GetClusters(itb);
1013 tilt = tracklets[ilr].GetTilt();
1015 t = 1./(x * x + y * y);
1016 uvt[0] = 2. * x * t;
1017 uvt[1] = 2. * x * t * tilt ;
1018 w = 2. * (y + tilt * (z - zVertex)) * t;
1019 error = 2. * TMath::Sqrt(cl->GetSigmaY2()) * t;
1020 fitter->AddPoint(uvt, w, error);
1026 // Calculate curvature
1027 Double_t a = fitter->GetParameter(0);
1028 Double_t b = fitter->GetParameter(1);
1029 Double_t curvature = a/TMath::Sqrt(b*b + 1);
1031 Float_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
1032 for(Int_t ip = 0; ip < AliTRDtrackerV1::kNPlanes; ip++)
1033 tracklets[ip].SetC(curvature);
1035 /* if(fReconstructor->GetStreamLevel() >= 5){
1036 //Linear Model on z-direction
1037 Double_t xref = CalculateReferenceX(tracklets); // Relative to the middle of the stack
1038 Double_t slope = fitter->GetParameter(2);
1039 Double_t zref = slope * xref;
1040 Float_t chi2Z = CalculateChi2Z(tracklets, zref, slope, xref);
1041 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1042 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1043 TTreeSRedirector &treeStreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
1044 treeStreamer << "FitTiltedRiemanConstraint"
1045 << "EventNumber=" << eventNumber
1046 << "CandidateNumber=" << candidateNumber
1047 << "Curvature=" << curvature
1048 << "Chi2Track=" << chi2track
1049 << "Chi2Z=" << chi2Z
1056 //_________________________________________________________________________
1057 Float_t AliTRDtrackerV1::FitTiltedRieman(AliTRDseedV1 *tracklets, Bool_t sigError)
1060 // Performs a Riemann fit taking tilting pad correction into account
1061 // The equation of a Riemann circle, where the y position is substituted by the
1062 // measured y-position taking pad tilting into account, has to be transformed
1063 // into a 4-dimensional hyperplane equation
1064 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
1065 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
1066 // zc: center of the pad row
1067 // zt: z-position of the track
1068 // The z-position of the track is assumed to be linear dependent on the x-position
1069 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
1070 // Transformation: u = 2 * x * t
1071 // v = 2 * tan(phiT) * t
1072 // w = 2 * tan(phiT) * (x - xref) * t
1073 // t = 1 / (x^2 + ymeas^2)
1074 // Parameters: a = -1/y0
1076 // c = (R^2 -x0^2 - y0^2)/y0
1079 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
1080 // results from the simple riemann fit. Afterwards the fit is redone.
1081 // The curvature is calculated according to the formula:
1082 // curv = a/(1 + b^2 + c*a) = 1/R
1084 // Paramters: - Array of tracklets (connected to the track candidate)
1085 // - Flag selecting the error definition
1086 // Output: - Chi2 values of the track (in Parameter list)
1088 TLinearFitter *fitter = GetTiltedRiemanFitter();
1089 fitter->StoreData(kTRUE);
1090 fitter->ClearPoints();
1091 AliTRDLeastSquare zfitter;
1092 AliTRDcluster *cl = 0x0;
1094 Double_t xref = CalculateReferenceX(tracklets);
1095 Double_t x, y, z, t, tilt, dx, w, we;
1098 // Containers for Least-square fitter
1099 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1100 if(!tracklets[ipl].IsOK()) continue;
1101 tilt = tracklets[ipl].GetTilt();
1102 for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
1103 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1104 if (!tracklets[ipl].IsUsable(itb)) continue;
1111 uvt[0] = 2. * x * t;
1113 uvt[2] = 2. * tilt * t;
1114 uvt[3] = 2. * tilt * dx * t;
1115 w = 2. * (y + tilt*z) * t;
1116 // error definition changes for the different calls
1118 we *= sigError ? TMath::Sqrt(cl->GetSigmaY2()) : 0.2;
1119 fitter->AddPoint(uvt, w, we);
1120 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1127 Double_t offset = fitter->GetParameter(3);
1128 Double_t slope = fitter->GetParameter(4);
1130 // Linear fitter - not possible to make boundaries
1131 // Do not accept non possible z and dzdx combinations
1132 Bool_t acceptablez = kTRUE;
1133 Double_t zref = 0.0;
1134 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
1135 if(!tracklets[iLayer].IsOK()) continue;
1136 zref = offset + slope * (tracklets[iLayer].GetX0() - xref);
1137 if (TMath::Abs(tracklets[iLayer].GetZfit(0) - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
1138 acceptablez = kFALSE;
1141 Double_t dzmf = zfitter.GetFunctionParameter(1);
1142 Double_t zmf = zfitter.GetFunctionValue(&xref);
1143 fgTiltedRieman->FixParameter(3, zmf);
1144 fgTiltedRieman->FixParameter(4, dzmf);
1146 fitter->ReleaseParameter(3);
1147 fitter->ReleaseParameter(4);
1148 offset = fitter->GetParameter(3);
1149 slope = fitter->GetParameter(4);
1152 // Calculate Curvarture
1153 Double_t a = fitter->GetParameter(0);
1154 Double_t b = fitter->GetParameter(1);
1155 Double_t c = fitter->GetParameter(2);
1156 Double_t curvature = 1.0 + b*b - c*a;
1157 if (curvature > 0.0)
1158 curvature = a / TMath::Sqrt(curvature);
1160 Double_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
1162 // Update the tracklets
1164 for(Int_t iLayer = 0; iLayer < AliTRDtrackerV1::kNPlanes; iLayer++) {
1166 x = tracklets[iLayer].GetX0();
1172 // y: R^2 = (x - x0)^2 + (y - y0)^2
1173 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1174 // R = Sqrt() = 1/Curvature
1175 // => y = y0 +/- Sqrt(1/Curvature^2 - (x - x0)^2)
1176 Double_t res = (x * a + b); // = (x - x0)/y0
1178 res = 1.0 - c * a + b * b - res; // = (R^2 - (x - x0)^2)/y0^2
1180 res = TMath::Sqrt(res);
1181 y = (1.0 - res) / a;
1184 // dy: R^2 = (x - x0)^2 + (y - y0)^2
1185 // => y = +/- Sqrt(R^2 - (x - x0)^2) + y0
1186 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1187 // Curvature: cr = 1/R = a/Sqrt(1 + b^2 - c*a)
1188 // => dy/dx = (x - x0)/(1/(cr^2) - (x - x0)^2)
1189 Double_t x0 = -b / a;
1190 if (-c * a + b * b + 1 > 0) {
1191 if (1.0/(curvature * curvature) - (x - x0) * (x - x0) > 0.0) {
1192 Double_t yderiv = (x - x0) / TMath::Sqrt(1.0/(curvature * curvature) - (x - x0) * (x - x0));
1193 if (a < 0) yderiv *= -1.0;
1197 z = offset + slope * (x - xref);
1199 tracklets[iLayer].SetYref(0, y);
1200 tracklets[iLayer].SetYref(1, dy);
1201 tracklets[iLayer].SetZref(0, z);
1202 tracklets[iLayer].SetZref(1, dz);
1203 tracklets[iLayer].SetC(curvature);
1204 tracklets[iLayer].SetChi2(chi2track);
1207 /* if(fReconstructor->GetStreamLevel() >=5){
1208 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
1209 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1210 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1211 Double_t chi2z = CalculateChi2Z(tracklets, offset, slope, xref);
1212 cstreamer << "FitTiltedRieman0"
1213 << "EventNumber=" << eventNumber
1214 << "CandidateNumber=" << candidateNumber
1216 << "Chi2Z=" << chi2z
1223 //____________________________________________________________________
1224 Double_t AliTRDtrackerV1::FitLine(const AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t err, Int_t np, AliTrackPoint *points)
1226 AliTRDLeastSquare yfitter, zfitter;
1227 AliTRDcluster *cl = 0x0;
1229 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1231 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1232 if(!(tracklet = track->GetTracklet(ipl))) continue;
1233 if(!tracklet->IsOK()) continue;
1234 new(&work[ipl]) AliTRDseedV1(*tracklet);
1236 tracklets = &work[0];
1239 Double_t xref = CalculateReferenceX(tracklets);
1240 Double_t x, y, z, dx, ye, yr, tilt;
1241 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1242 if(!tracklets[ipl].IsOK()) continue;
1243 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1244 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1245 if (!tracklets[ipl].IsUsable(itb)) continue;
1249 zfitter.AddPoint(&dx, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1253 Double_t z0 = zfitter.GetFunctionParameter(0);
1254 Double_t dzdx = zfitter.GetFunctionParameter(1);
1255 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1256 if(!tracklets[ipl].IsOK()) continue;
1257 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1258 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1259 if (!tracklets[ipl].IsUsable(itb)) continue;
1263 tilt = tracklets[ipl].GetTilt();
1265 yr = y + tilt*(z - z0 - dzdx*dx);
1266 // error definition changes for the different calls
1267 ye = tilt*TMath::Sqrt(cl->GetSigmaZ2());
1268 ye += err ? tracklets[ipl].GetSigmaY() : 0.2;
1269 yfitter.AddPoint(&dx, yr, ye);
1273 Double_t y0 = yfitter.GetFunctionParameter(0);
1274 Double_t dydx = yfitter.GetFunctionParameter(1);
1275 Double_t chi2 = 0.;//yfitter.GetChisquare()/Double_t(nPoints);
1277 //update track points array
1280 for(int ip=0; ip<np; ip++){
1281 points[ip].GetXYZ(xyz);
1282 xyz[1] = y0 + dydx * (xyz[0] - xref);
1283 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1284 points[ip].SetXYZ(xyz);
1291 //_________________________________________________________________________
1292 Double_t AliTRDtrackerV1::FitRiemanTilt(const AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t sigError, Int_t np, AliTrackPoint *points)
1295 // Performs a Riemann fit taking tilting pad correction into account
1296 // The equation of a Riemann circle, where the y position is substituted by the
1297 // measured y-position taking pad tilting into account, has to be transformed
1298 // into a 4-dimensional hyperplane equation
1299 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
1300 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
1301 // zc: center of the pad row
1302 // zt: z-position of the track
1303 // The z-position of the track is assumed to be linear dependent on the x-position
1304 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
1305 // Transformation: u = 2 * x * t
1306 // v = 2 * tan(phiT) * t
1307 // w = 2 * tan(phiT) * (x - xref) * t
1308 // t = 1 / (x^2 + ymeas^2)
1309 // Parameters: a = -1/y0
1311 // c = (R^2 -x0^2 - y0^2)/y0
1314 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
1315 // results from the simple riemann fit. Afterwards the fit is redone.
1316 // The curvature is calculated according to the formula:
1317 // curv = a/(1 + b^2 + c*a) = 1/R
1319 // Paramters: - Array of tracklets (connected to the track candidate)
1320 // - Flag selecting the error definition
1321 // Output: - Chi2 values of the track (in Parameter list)
1323 TLinearFitter *fitter = GetTiltedRiemanFitter();
1324 fitter->StoreData(kTRUE);
1325 fitter->ClearPoints();
1326 AliTRDLeastSquare zfitter;
1327 AliTRDcluster *cl = 0x0;
1329 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1331 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1332 if(!(tracklet = track->GetTracklet(ipl))) continue;
1333 if(!tracklet->IsOK()) continue;
1334 new(&work[ipl]) AliTRDseedV1(*tracklet);
1336 tracklets = &work[0];
1339 Double_t xref = CalculateReferenceX(tracklets);
1340 Double_t x, y, z, t, tilt, dx, w, we;
1343 // Containers for Least-square fitter
1344 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1345 if(!tracklets[ipl].IsOK()) continue;
1346 for(Int_t itb = 0; itb < AliTRDseedV1::kNclusters; itb++){
1347 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1348 if (!tracklets[ipl].IsUsable(itb)) continue;
1352 tilt = tracklets[ipl].GetTilt();
1356 uvt[0] = 2. * x * t;
1358 uvt[2] = 2. * tilt * t;
1359 uvt[3] = 2. * tilt * dx * t;
1360 w = 2. * (y + tilt*z) * t;
1361 // error definition changes for the different calls
1363 we *= sigError ? TMath::Sqrt(cl->GetSigmaY2()) : 0.2;
1364 fitter->AddPoint(uvt, w, we);
1365 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1369 if(fitter->Eval()) return 1.E10;
1371 Double_t z0 = fitter->GetParameter(3);
1372 Double_t dzdx = fitter->GetParameter(4);
1375 // Linear fitter - not possible to make boundaries
1376 // Do not accept non possible z and dzdx combinations
1377 Bool_t accept = kTRUE;
1378 Double_t zref = 0.0;
1379 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
1380 if(!tracklets[iLayer].IsOK()) continue;
1381 zref = z0 + dzdx * (tracklets[iLayer].GetX0() - xref);
1382 if (TMath::Abs(tracklets[iLayer].GetZfit(0) - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
1387 Double_t dzmf = zfitter.GetFunctionParameter(1);
1388 Double_t zmf = zfitter.GetFunctionValue(&xref);
1389 fitter->FixParameter(3, zmf);
1390 fitter->FixParameter(4, dzmf);
1392 fitter->ReleaseParameter(3);
1393 fitter->ReleaseParameter(4);
1394 z0 = fitter->GetParameter(3); // = zmf ?
1395 dzdx = fitter->GetParameter(4); // = dzmf ?
1398 // Calculate Curvature
1399 Double_t a = fitter->GetParameter(0);
1400 Double_t b = fitter->GetParameter(1);
1401 Double_t c = fitter->GetParameter(2);
1402 Double_t y0 = 1. / a;
1403 Double_t x0 = -b * y0;
1404 Double_t tmp = y0*y0 + x0*x0 - c*y0;
1405 if(tmp<=0.) return 1.E10;
1406 Double_t R = TMath::Sqrt(tmp);
1407 Double_t C = 1.0 + b*b - c*a;
1408 if (C > 0.0) C = a / TMath::Sqrt(C);
1410 // Calculate chi2 of the fit
1411 Double_t chi2 = fitter->GetChisquare()/Double_t(nPoints);
1413 // Update the tracklets
1415 for(Int_t ip = 0; ip < kNPlanes; ip++) {
1416 x = tracklets[ip].GetX0();
1417 tmp = R*R-(x-x0)*(x-x0);
1418 if(tmp <= 0.) continue;
1419 tmp = TMath::Sqrt(tmp);
1421 // y: R^2 = (x - x0)^2 + (y - y0)^2
1422 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1423 tracklets[ip].SetYref(0, y0 - (y0>0.?1.:-1)*tmp);
1424 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1425 tracklets[ip].SetYref(1, (x - x0) / tmp);
1426 tracklets[ip].SetZref(0, z0 + dzdx * (x - xref));
1427 tracklets[ip].SetZref(1, dzdx);
1428 tracklets[ip].SetC(C);
1429 tracklets[ip].SetChi2(chi2);
1432 //update track points array
1435 for(int ip=0; ip<np; ip++){
1436 points[ip].GetXYZ(xyz);
1437 xyz[1] = TMath::Abs(xyz[0] - x0) > R ? 100. : y0 - (y0>0.?1.:-1.)*TMath::Sqrt((R-(xyz[0]-x0))*(R+(xyz[0]-x0)));
1438 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1439 points[ip].SetXYZ(xyz);
1447 //____________________________________________________________________
1448 Double_t AliTRDtrackerV1::FitKalman(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t up, Int_t np, AliTrackPoint *points)
1450 // Kalman filter implementation for the TRD.
1451 // It returns the positions of the fit in the array "points"
1453 // Author : A.Bercuci@gsi.de
1455 // printf("Start track @ x[%f]\n", track->GetX());
1457 //prepare marker points along the track
1458 Int_t ip = np ? 0 : 1;
1460 if((up?-1:1) * (track->GetX() - points[ip].GetX()) > 0.) break;
1461 //printf("AliTRDtrackerV1::FitKalman() : Skip track marker x[%d] = %7.3f. Before track start ( %7.3f ).\n", ip, points[ip].GetX(), track->GetX());
1464 //if(points) printf("First marker point @ x[%d] = %f\n", ip, points[ip].GetX());
1467 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
1469 //Loop through the TRD planes
1470 for (Int_t jplane = 0; jplane < kNPlanes; jplane++) {
1471 // GET TRACKLET OR BUILT IT
1472 Int_t iplane = up ? jplane : kNPlanes - 1 - jplane;
1474 if(!(ptrTracklet = &tracklets[iplane])) continue;
1476 if(!(ptrTracklet = track->GetTracklet(iplane))){
1477 /*AliTRDtrackerV1 *tracker = 0x0;
1478 if(!(tracker = dynamic_cast<AliTRDtrackerV1*>( AliTRDReconstructor::Tracker()))) continue;
1479 ptrTracklet = new(&tracklet) AliTRDseedV1(iplane);
1480 if(!tracker->MakeTracklet(ptrTracklet, track)) */
1484 if(!ptrTracklet->IsOK()) continue;
1486 Double_t x = ptrTracklet->GetX0();
1489 //don't do anything if next marker is after next update point.
1490 if((up?-1:1) * (points[ip].GetX() - x) - fgkMaxStep < 0) break;
1491 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1493 Double_t xyz[3]; // should also get the covariance
1495 track->Global2LocalPosition(xyz, track->GetAlpha());
1496 points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1499 // printf("plane[%d] tracklet[%p] x[%f]\n", iplane, ptrTracklet, x);
1501 // Propagate closer to the next update point
1502 if(((up?-1:1) * (x - track->GetX()) + fgkMaxStep < 0) && !PropagateToX(*track, x + (up?-1:1)*fgkMaxStep, fgkMaxStep)) return -1.;
1504 if(!AdjustSector(track)) return -1;
1505 if(TMath::Abs(track->GetSnp()) > fgkMaxSnp) return -1;
1507 //load tracklet to the tracker and the track
1509 if((index = FindTracklet(ptrTracklet)) < 0){
1510 ptrTracklet = SetTracklet(&tracklet);
1511 index = fTracklets->GetEntriesFast()-1;
1513 track->SetTracklet(ptrTracklet, index);*/
1516 // register tracklet to track with tracklet creation !!
1517 // PropagateBack : loaded tracklet to the tracker and update index
1518 // RefitInward : update index
1519 // MakeTrack : loaded tracklet to the tracker and update index
1520 if(!tracklets) track->SetTracklet(ptrTracklet, -1);
1523 //Calculate the mean material budget along the path inside the chamber
1524 Double_t xyz0[3]; track->GetXYZ(xyz0);
1525 Double_t alpha = track->GetAlpha();
1526 Double_t xyz1[3], y, z;
1527 if(!track->GetProlongation(x, y, z)) return -1;
1528 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
1529 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
1531 if((xyz0[0] - xyz1[9] < 1e-3) && (xyz0[0] - xyz1[9] < 1e-3)) continue; // check wheter we are at the same global x position
1533 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param) <=0.) break;
1534 Double_t xrho = param[0]*param[4]; // density*length
1535 Double_t xx0 = param[1]; // radiation length
1537 //Propagate the track
1538 track->PropagateTo(x, xx0, xrho);
1539 if (!AdjustSector(track)) break;
1542 Double_t chi2 = track->GetPredictedChi2(ptrTracklet);
1543 if(chi2<1e+10) track->Update(ptrTracklet, chi2);
1546 //Reset material budget if 2 consecutive gold
1547 if(iplane>0 && track->GetTracklet(iplane-1) && ptrTracklet->GetN() + track->GetTracklet(iplane-1)->GetN() > 20) track->SetBudget(2, 0.);
1548 } // end planes loop
1552 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1554 Double_t xyz[3]; // should also get the covariance
1556 track->Global2LocalPosition(xyz, track->GetAlpha());
1557 points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1561 return track->GetChi2();
1564 //_________________________________________________________________________
1565 Float_t AliTRDtrackerV1::CalculateChi2Z(AliTRDseedV1 *tracklets, Double_t offset, Double_t slope, Double_t xref)
1568 // Calculates the chi2-value of the track in z-Direction including tilting pad correction.
1569 // A linear dependence on the x-value serves as a model.
1570 // The parameters are related to the tilted Riemann fit.
1571 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
1572 // - the offset for the reference x
1574 // - the reference x position
1575 // Output: - The Chi2 value of the track in z-Direction
1577 Float_t chi2Z = 0, nLayers = 0;
1578 for (Int_t iLayer = 0; iLayer < AliTRDgeometry::kNlayer; iLayer++) {
1579 if(!tracklets[iLayer].IsOK()) continue;
1580 Double_t z = offset + slope * (tracklets[iLayer].GetX0() - xref);
1581 chi2Z += TMath::Abs(tracklets[iLayer].GetZfit(0) - z);
1584 chi2Z /= TMath::Max((nLayers - 3.0),1.0);
1588 //_____________________________________________________________________________
1589 Int_t AliTRDtrackerV1::PropagateToX(AliTRDtrackV1 &t, Double_t xToGo, Double_t maxStep)
1592 // Starting from current X-position of track <t> this function
1593 // extrapolates the track up to radial position <xToGo>.
1594 // Returns 1 if track reaches the plane, and 0 otherwise
1597 const Double_t kEpsilon = 0.00001;
1599 // Current track X-position
1600 Double_t xpos = t.GetX();
1602 // Direction: inward or outward
1603 Double_t dir = (xpos < xToGo) ? 1.0 : -1.0;
1605 while (((xToGo - xpos) * dir) > kEpsilon) {
1614 // The next step size
1615 Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
1617 // Get the global position of the starting point
1620 // X-position after next step
1623 // Get local Y and Z at the X-position of the next step
1624 if (!t.GetProlongation(x,y,z)) {
1625 return 0; // No prolongation possible
1628 // The global position of the end point of this prolongation step
1629 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1630 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1633 // Calculate the mean material budget between start and
1634 // end point of this prolongation step
1635 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) return 0;
1637 // Propagate the track to the X-position after the next step
1638 if (!t.PropagateTo(x, param[1], param[0]*param[4])) return 0;
1640 // Rotate the track if necessary
1643 // New track X-position
1653 //_____________________________________________________________________________
1654 Int_t AliTRDtrackerV1::ReadClusters(TClonesArray* &array, TTree *clusterTree) const
1657 // Reads AliTRDclusters from the file.
1658 // The names of the cluster tree and branches
1659 // should match the ones used in AliTRDclusterizer::WriteClusters()
1662 Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
1663 TObjArray *clusterArray = new TObjArray(nsize+1000);
1665 TBranch *branch = clusterTree->GetBranch("TRDcluster");
1667 AliError("Can't get the branch !");
1670 branch->SetAddress(&clusterArray);
1673 Float_t nclusters = fReconstructor->GetRecoParam()->GetNClusters();
1674 if(fReconstructor->IsHLT()) nclusters /= AliTRDgeometry::kNsector;
1675 array = new TClonesArray("AliTRDcluster", Int_t(nclusters));
1676 array->SetOwner(kTRUE);
1679 // Loop through all entries in the tree
1680 Int_t nEntries = (Int_t) clusterTree->GetEntries();
1683 AliTRDcluster *c = 0x0;
1684 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1686 nbytes += clusterTree->GetEvent(iEntry);
1688 // Get the number of points in the detector
1689 Int_t nCluster = clusterArray->GetEntriesFast();
1690 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1691 if(!(c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster))) continue;
1693 new((*fClusters)[ncl++]) AliTRDcluster(*c);
1694 delete (clusterArray->RemoveAt(iCluster));
1698 delete clusterArray;
1703 //_____________________________________________________________________________
1704 Int_t AliTRDtrackerV1::LoadClusters(TTree *cTree)
1707 // Fills clusters into TRD tracking sectors
1710 if(!fReconstructor->IsWritingClusters()){
1711 fClusters = AliTRDReconstructor::GetClusters();
1713 if (ReadClusters(fClusters, cTree)) {
1714 AliError("Problem with reading the clusters !");
1720 if(!fClusters || !fClusters->GetEntriesFast()){
1721 AliInfo("No TRD clusters");
1726 BuildTrackingContainers();
1728 //Int_t ncl = fClusters->GetEntriesFast();
1729 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1734 //_____________________________________________________________________________
1735 Int_t AliTRDtrackerV1::LoadClusters(TClonesArray *clusters)
1738 // Fills clusters into TRD tracking sectors
1739 // Function for use in the HLT
1741 if(!clusters || !clusters->GetEntriesFast()){
1742 AliInfo("No TRD clusters");
1746 fClusters = clusters;
1750 BuildTrackingContainers();
1752 //Int_t ncl = fClusters->GetEntriesFast();
1753 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1759 //____________________________________________________________________
1760 Int_t AliTRDtrackerV1::BuildTrackingContainers()
1762 // Building tracking containers for clusters
1764 Int_t nin =0, icl = fClusters->GetEntriesFast();
1766 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(icl);
1767 if(c->IsInChamber()) nin++;
1768 Int_t detector = c->GetDetector();
1769 Int_t sector = fGeom->GetSector(detector);
1770 Int_t stack = fGeom->GetStack(detector);
1771 Int_t layer = fGeom->GetLayer(detector);
1773 fTrSec[sector].GetChamber(stack, layer, kTRUE)->InsertCluster(c, icl);
1776 const AliTRDCalDet *cal = AliTRDcalibDB::Instance()->GetT0Det();
1777 for(int isector =0; isector<AliTRDgeometry::kNsector; isector++){
1778 if(!fTrSec[isector].GetNChambers()) continue;
1779 fTrSec[isector].Init(fReconstructor, cal);
1787 //____________________________________________________________________
1788 void AliTRDtrackerV1::UnloadClusters()
1791 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1794 if(fTracks) fTracks->Delete();
1795 if(fTracklets) fTracklets->Delete();
1797 if(IsClustersOwner()) fClusters->Delete();
1799 // save clusters array in the reconstructor for further use.
1800 if(!fReconstructor->IsWritingClusters()){
1801 AliTRDReconstructor::SetClusters(fClusters);
1802 SetClustersOwner(kFALSE);
1803 } else AliTRDReconstructor::SetClusters(0x0);
1806 for (int i = 0; i < AliTRDgeometry::kNsector; i++) fTrSec[i].Clear();
1808 // Increment the Event Number
1809 AliTRDtrackerDebug::SetEventNumber(AliTRDtrackerDebug::GetEventNumber() + 1);
1812 // //____________________________________________________________________
1813 // void AliTRDtrackerV1::UseClusters(const AliKalmanTrack *t, Int_t) const
1815 // const AliTRDtrackV1 *track = dynamic_cast<const AliTRDtrackV1*>(t);
1816 // if(!track) return;
1818 // AliTRDseedV1 *tracklet = 0x0;
1819 // for(Int_t ily=AliTRDgeometry::kNlayer; ily--;){
1820 // if(!(tracklet = track->GetTracklet(ily))) continue;
1821 // AliTRDcluster *c = 0x0;
1822 // for(Int_t ic=AliTRDseed::kNclusters; ic--;){
1823 // if(!(c=tracklet->GetClusters(ic))) continue;
1830 //_____________________________________________________________________________
1831 Bool_t AliTRDtrackerV1::AdjustSector(AliTRDtrackV1 *track)
1834 // Rotates the track when necessary
1837 Double_t alpha = AliTRDgeometry::GetAlpha();
1838 Double_t y = track->GetY();
1839 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
1842 if (!track->Rotate( alpha)) {
1846 else if (y < -ymax) {
1847 if (!track->Rotate(-alpha)) {
1857 //____________________________________________________________________
1858 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(AliTRDtrackV1 *track, Int_t p, Int_t &idx)
1860 // Find tracklet for TRD track <track>
1869 // Detailed description
1871 idx = track->GetTrackletIndex(p);
1872 AliTRDseedV1 *tracklet = (idx==0xffff) ? 0x0 : (AliTRDseedV1*)fTracklets->UncheckedAt(idx);
1877 //____________________________________________________________________
1878 AliTRDseedV1* AliTRDtrackerV1::SetTracklet(AliTRDseedV1 *tracklet)
1880 // Add this tracklet to the list of tracklets stored in the tracker
1883 // - tracklet : pointer to the tracklet to be added to the list
1886 // - the index of the new tracklet in the tracker tracklets list
1888 // Detailed description
1889 // Build the tracklets list if it is not yet created (late initialization)
1890 // and adds the new tracklet to the list.
1893 fTracklets = new TClonesArray("AliTRDseedV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1894 fTracklets->SetOwner(kTRUE);
1896 Int_t nentries = fTracklets->GetEntriesFast();
1897 return new ((*fTracklets)[nentries]) AliTRDseedV1(*tracklet);
1900 //____________________________________________________________________
1901 AliTRDtrackV1* AliTRDtrackerV1::SetTrack(AliTRDtrackV1 *track)
1903 // Add this track to the list of tracks stored in the tracker
1906 // - track : pointer to the track to be added to the list
1909 // - the pointer added
1911 // Detailed description
1912 // Build the tracks list if it is not yet created (late initialization)
1913 // and adds the new track to the list.
1916 fTracks = new TClonesArray("AliTRDtrackV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1917 fTracks->SetOwner(kTRUE);
1919 Int_t nentries = fTracks->GetEntriesFast();
1920 return new ((*fTracks)[nentries]) AliTRDtrackV1(*track);
1925 //____________________________________________________________________
1926 Int_t AliTRDtrackerV1::Clusters2TracksSM(Int_t sector, AliESDEvent *esd)
1929 // Steer tracking for one SM.
1932 // sector : Array of (SM) propagation layers containing clusters
1933 // esd : The current ESD event. On output it contains the also
1934 // the ESD (TRD) tracks found in this SM.
1937 // Number of tracks found in this TRD supermodule.
1939 // Detailed description
1941 // 1. Unpack AliTRDpropagationLayers objects for each stack.
1942 // 2. Launch stack tracking.
1943 // See AliTRDtrackerV1::Clusters2TracksStack() for details.
1944 // 3. Pack results in the ESD event.
1947 // allocate space for esd tracks in this SM
1948 TClonesArray esdTrackList("AliESDtrack", 2*kMaxTracksStack);
1949 esdTrackList.SetOwner();
1952 Int_t nChambers = 0;
1953 AliTRDtrackingChamber **stack = 0x0, *chamber = 0x0;
1954 for(int istack = 0; istack<AliTRDgeometry::kNstack; istack++){
1955 if(!(stack = fTrSec[sector].GetStack(istack))) continue;
1957 for(int ilayer=0; ilayer<AliTRDgeometry::kNlayer; ilayer++){
1958 if(!(chamber = stack[ilayer])) continue;
1959 if(chamber->GetNClusters() < fgNTimeBins * fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
1961 //AliInfo(Form("sector %d stack %d layer %d clusters %d", sector, istack, ilayer, chamber->GetNClusters()));
1963 if(nChambers < 4) continue;
1964 //AliInfo(Form("Doing stack %d", istack));
1965 nTracks += Clusters2TracksStack(stack, &esdTrackList);
1967 //AliInfo(Form("Found %d tracks in SM %d [%d]\n", nTracks, sector, esd->GetNumberOfTracks()));
1969 for(int itrack=0; itrack<nTracks; itrack++)
1970 esd->AddTrack((AliESDtrack*)esdTrackList[itrack]);
1972 // Reset Track and Candidate Number
1973 AliTRDtrackerDebug::SetCandidateNumber(0);
1974 AliTRDtrackerDebug::SetTrackNumber(0);
1978 //____________________________________________________________________
1979 Int_t AliTRDtrackerV1::Clusters2TracksStack(AliTRDtrackingChamber **stack, TClonesArray *esdTrackList)
1982 // Make tracks in one TRD stack.
1985 // layer : Array of stack propagation layers containing clusters
1986 // esdTrackList : Array of ESD tracks found by the stand alone tracker.
1987 // On exit the tracks found in this stack are appended.
1990 // Number of tracks found in this stack.
1992 // Detailed description
1994 // 1. Find the 3 most useful seeding chambers. See BuildSeedingConfigs() for details.
1995 // 2. Steer AliTRDtrackerV1::MakeSeeds() for 3 seeding layer configurations.
1996 // See AliTRDtrackerV1::MakeSeeds() for more details.
1997 // 3. Arrange track candidates in decreasing order of their quality
1998 // 4. Classify tracks in 5 categories according to:
1999 // a) number of layers crossed
2001 // 5. Sign clusters by tracks in decreasing order of track quality
2002 // 6. Build AliTRDtrack out of seeding tracklets
2004 // 8. Build ESD track and register it to the output list
2007 const AliTRDCalDet *cal = AliTRDcalibDB::Instance()->GetT0Det();
2008 AliTRDtrackingChamber *chamber = 0x0;
2009 AliTRDtrackingChamber **ci = 0x0;
2010 AliTRDseedV1 sseed[kMaxTracksStack*6]; // to be initialized
2011 Int_t pars[4]; // MakeSeeds parameters
2013 //Double_t alpha = AliTRDgeometry::GetAlpha();
2014 //Double_t shift = .5 * alpha;
2015 Int_t configs[kNConfigs];
2017 // Purge used clusters from the containers
2019 for(Int_t ic = kNPlanes; ic--; ci++){
2020 if(!(*ci)) continue;
2024 // Build initial seeding configurations
2025 Double_t quality = BuildSeedingConfigs(stack, configs);
2026 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 10){
2027 AliInfo(Form("Plane config %d %d %d Quality %f"
2028 , configs[0], configs[1], configs[2], quality));
2032 // Initialize contors
2033 Int_t ntracks, // number of TRD track candidates
2034 ntracks1, // number of registered TRD tracks/iter
2035 ntracks2 = 0; // number of all registered TRD tracks in stack
2039 Int_t ic = 0; ci = &stack[0];
2040 while(ic<kNPlanes && !(*ci)){ic++; ci++;}
2041 if(!(*ci)) return ntracks2;
2042 Int_t istack = fGeom->GetStack((*ci)->GetDetector());
2045 // Loop over seeding configurations
2046 ntracks = 0; ntracks1 = 0;
2047 for (Int_t iconf = 0; iconf<3; iconf++) {
2048 pars[0] = configs[iconf];
2051 ntracks = MakeSeeds(stack, &sseed[6*ntracks], pars);
2052 if(ntracks == kMaxTracksStack) break;
2054 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 10) AliInfo(Form("Candidate TRD tracks %d in iteration %d.", ntracks, fSieveSeeding));
2058 // Sort the seeds according to their quality
2059 Int_t sort[kMaxTracksStack];
2060 TMath::Sort(ntracks, fTrackQuality, sort, kTRUE);
2062 // Initialize number of tracks so far and logic switches
2063 Int_t ntracks0 = esdTrackList->GetEntriesFast();
2064 Bool_t signedTrack[kMaxTracksStack];
2065 Bool_t fakeTrack[kMaxTracksStack];
2066 for (Int_t i=0; i<ntracks; i++){
2067 signedTrack[i] = kFALSE;
2068 fakeTrack[i] = kFALSE;
2070 //AliInfo("Selecting track candidates ...");
2072 // Sieve clusters in decreasing order of track quality
2073 Double_t trackParams[7];
2074 // AliTRDseedV1 *lseed = 0x0;
2075 Int_t jSieve = 0, candidates;
2077 //AliInfo(Form("\t\tITER = %i ", jSieve));
2079 // Check track candidates
2081 for (Int_t itrack = 0; itrack < ntracks; itrack++) {
2082 Int_t trackIndex = sort[itrack];
2083 if (signedTrack[trackIndex] || fakeTrack[trackIndex]) continue;
2086 // Calculate track parameters from tracklets seeds
2091 for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
2092 Int_t jseed = kNPlanes*trackIndex+jLayer;
2093 if(!sseed[jseed].IsOK()) continue;
2094 if (TMath::Abs(sseed[jseed].GetYref(0) / sseed[jseed].GetX0()) < 0.158) findable++;
2095 // TODO here we get a sig fault which should never happen !
2096 sseed[jseed].UpdateUsed();
2097 ncl += sseed[jseed].GetN2();
2098 nused += sseed[jseed].GetNUsed();
2102 // Filter duplicated tracks
2104 //printf("Skip %d nused %d\n", trackIndex, nused);
2105 fakeTrack[trackIndex] = kTRUE;
2108 if (Float_t(nused)/ncl >= .25){
2109 //printf("Skip %d nused/ncl >= .25\n", trackIndex);
2110 fakeTrack[trackIndex] = kTRUE;
2115 Bool_t skip = kFALSE;
2118 if(nlayers < 6) {skip = kTRUE; break;}
2119 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
2123 if(nlayers < findable){skip = kTRUE; break;}
2124 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -4.){skip = kTRUE; break;}
2128 if ((nlayers == findable) || (nlayers == 6)) { skip = kTRUE; break;}
2129 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -6.0){skip = kTRUE; break;}
2133 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
2137 if (nlayers == 3){skip = kTRUE; break;}
2138 //if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) - nused/(nlayers-3.0) < -15.0){skip = kTRUE; break;}
2143 //printf("REJECTED : %d [%d] nlayers %d trackQuality = %e nused %d\n", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused);
2146 signedTrack[trackIndex] = kTRUE;
2148 // Build track parameters
2149 AliTRDseedV1 *lseed =&sseed[trackIndex*6];
2151 while(idx<3 && !lseed->IsOK()) {
2155 Double_t x = lseed->GetX0();// - 3.5;
2156 trackParams[0] = x; //NEW AB
2157 trackParams[1] = lseed->GetYref(0); // lseed->GetYat(x);
2158 trackParams[2] = lseed->GetZref(0); // lseed->GetZat(x);
2159 trackParams[3] = TMath::Sin(TMath::ATan(lseed->GetYref(1)));
2160 trackParams[4] = lseed->GetZref(1) / TMath::Sqrt(1. + lseed->GetYref(1) * lseed->GetYref(1));
2161 trackParams[5] = lseed->GetC();
2162 Int_t ich = 0; while(!(chamber = stack[ich])) ich++;
2163 trackParams[6] = fGeom->GetSector(chamber->GetDetector());/* *alpha+shift; // Supermodule*/
2165 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2166 //AliInfo(Form("Track %d [%d] nlayers %d trackQuality = %e nused %d, yref = %3.3f", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused, trackParams[1]));
2168 AliTRDseedV1 *dseed[6];
2169 for(Int_t iseed = AliTRDgeometry::kNlayer; iseed--;) dseed[iseed] = new AliTRDseedV1(lseed[iseed]);
2171 //Int_t eventNrInFile = esd->GetEventNumberInFile();
2172 //AliInfo(Form("Number of clusters %d.", nclusters));
2173 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2174 Int_t trackNumber = AliTRDtrackerDebug::GetTrackNumber();
2175 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2176 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2177 cstreamer << "Clusters2TracksStack"
2178 << "EventNumber=" << eventNumber
2179 << "TrackNumber=" << trackNumber
2180 << "CandidateNumber=" << candidateNumber
2181 << "Iter=" << fSieveSeeding
2182 << "Like=" << fTrackQuality[trackIndex]
2183 << "S0.=" << dseed[0]
2184 << "S1.=" << dseed[1]
2185 << "S2.=" << dseed[2]
2186 << "S3.=" << dseed[3]
2187 << "S4.=" << dseed[4]
2188 << "S5.=" << dseed[5]
2189 << "p0=" << trackParams[0]
2190 << "p1=" << trackParams[1]
2191 << "p2=" << trackParams[2]
2192 << "p3=" << trackParams[3]
2193 << "p4=" << trackParams[4]
2194 << "p5=" << trackParams[5]
2195 << "p6=" << trackParams[6]
2197 << "NLayers=" << nlayers
2198 << "Findable=" << findable
2199 << "NUsed=" << nused
2203 AliTRDtrackV1 *track = MakeTrack(&sseed[trackIndex*kNPlanes], trackParams);
2205 AliWarning("Fail to build a TRD Track.");
2209 //AliInfo("End of MakeTrack()");
2210 AliESDtrack *esdTrack = new ((*esdTrackList)[ntracks0++]) AliESDtrack();
2211 esdTrack->UpdateTrackParams(track, AliESDtrack::kTRDout);
2212 esdTrack->SetLabel(track->GetLabel());
2213 track->UpdateESDtrack(esdTrack);
2214 // write ESD-friends if neccessary
2215 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0){
2216 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(*track);
2217 calibTrack->SetOwner();
2218 esdTrack->AddCalibObject(calibTrack);
2221 AliTRDtrackerDebug::SetTrackNumber(AliTRDtrackerDebug::GetTrackNumber() + 1);
2225 } while(jSieve<5 && candidates); // end track candidates sieve
2226 if(!ntracks1) break;
2228 // increment counters
2229 ntracks2 += ntracks1;
2231 if(fReconstructor->IsHLT()) break;
2234 // Rebuild plane configurations and indices taking only unused clusters into account
2235 quality = BuildSeedingConfigs(stack, configs);
2236 if(quality < 1.E-7) break; //fReconstructor->GetRecoParam() ->GetPlaneQualityThreshold()) break;
2238 for(Int_t ip = 0; ip < kNPlanes; ip++){
2239 if(!(chamber = stack[ip])) continue;
2240 chamber->Build(fGeom, cal);//Indices(fSieveSeeding);
2243 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 10){
2244 AliInfo(Form("Sieve level %d Plane config %d %d %d Quality %f", fSieveSeeding, configs[0], configs[1], configs[2], quality));
2246 } while(fSieveSeeding<10); // end stack clusters sieve
2250 //AliInfo(Form("Registered TRD tracks %d in stack %d.", ntracks2, pars[1]));
2255 //___________________________________________________________________
2256 Double_t AliTRDtrackerV1::BuildSeedingConfigs(AliTRDtrackingChamber **stack, Int_t *configs)
2259 // Assign probabilities to chambers according to their
2260 // capability of producing seeds.
2264 // layers : Array of stack propagation layers for all 6 chambers in one stack
2265 // configs : On exit array of configuration indexes (see GetSeedingConfig()
2266 // for details) in the decreasing order of their seeding probabilities.
2270 // Return top configuration quality
2272 // Detailed description:
2274 // To each chamber seeding configuration (see GetSeedingConfig() for
2275 // the list of all configurations) one defines 2 quality factors:
2276 // - an apriori topological quality (see GetSeedingConfig() for details) and
2277 // - a data quality based on the uniformity of the distribution of
2278 // clusters over the x range (time bins population). See CookChamberQA() for details.
2279 // The overall chamber quality is given by the product of this 2 contributions.
2282 Double_t chamberQ[kNPlanes];memset(chamberQ, 0, kNPlanes*sizeof(Double_t));
2283 AliTRDtrackingChamber *chamber = 0x0;
2284 for(int iplane=0; iplane<kNPlanes; iplane++){
2285 if(!(chamber = stack[iplane])) continue;
2286 chamberQ[iplane] = (chamber = stack[iplane]) ? chamber->GetQuality() : 0.;
2289 Double_t tconfig[kNConfigs];memset(tconfig, 0, kNConfigs*sizeof(Double_t));
2290 Int_t planes[] = {0, 0, 0, 0};
2291 for(int iconf=0; iconf<kNConfigs; iconf++){
2292 GetSeedingConfig(iconf, planes);
2293 tconfig[iconf] = fgTopologicQA[iconf];
2294 for(int iplane=0; iplane<4; iplane++) tconfig[iconf] *= chamberQ[planes[iplane]];
2297 TMath::Sort((Int_t)kNConfigs, tconfig, configs, kTRUE);
2298 // AliInfo(Form("q[%d] = %f", configs[0], tconfig[configs[0]]));
2299 // AliInfo(Form("q[%d] = %f", configs[1], tconfig[configs[1]]));
2300 // AliInfo(Form("q[%d] = %f", configs[2], tconfig[configs[2]]));
2302 return tconfig[configs[0]];
2305 //____________________________________________________________________
2306 Int_t AliTRDtrackerV1::MakeSeeds(AliTRDtrackingChamber **stack, AliTRDseedV1 *sseed, Int_t *ipar)
2309 // Make tracklet seeds in the TRD stack.
2312 // layers : Array of stack propagation layers containing clusters
2313 // sseed : Array of empty tracklet seeds. On exit they are filled.
2314 // ipar : Control parameters:
2315 // ipar[0] -> seeding chambers configuration
2316 // ipar[1] -> stack index
2317 // ipar[2] -> number of track candidates found so far
2320 // Number of tracks candidates found.
2322 // Detailed description
2324 // The following steps are performed:
2325 // 1. Select seeding layers from seeding chambers
2326 // 2. Select seeding clusters from the seeding AliTRDpropagationLayerStack.
2327 // The clusters are taken from layer 3, layer 0, layer 1 and layer 2, in
2328 // this order. The parameters controling the range of accepted clusters in
2329 // layer 0, 1, and 2 are defined in AliTRDchamberTimeBin::BuildCond().
2330 // 3. Helix fit of the cluster set. (see AliTRDtrackerFitter::FitRieman(AliTRDcluster**))
2331 // 4. Initialize seeding tracklets in the seeding chambers.
2333 // Chi2 in the Y direction less than threshold ... (1./(3. - sLayer))
2334 // Chi2 in the Z direction less than threshold ... (1./(3. - sLayer))
2335 // 6. Attach clusters to seeding tracklets and find linear approximation of
2336 // the tracklet (see AliTRDseedV1::AttachClustersIter()). The number of used
2337 // clusters used by current seeds should not exceed ... (25).
2339 // All 4 seeding tracklets should be correctly constructed (see
2340 // AliTRDseedV1::AttachClustersIter())
2341 // 8. Helix fit of the seeding tracklets
2343 // Likelihood calculation of the fit. (See AliTRDtrackerV1::CookLikelihood() for details)
2344 // 10. Extrapolation of the helix fit to the other 2 chambers:
2345 // a) Initialization of extrapolation tracklet with fit parameters
2346 // b) Helix fit of tracklets
2347 // c) Attach clusters and linear interpolation to extrapolated tracklets
2348 // d) Helix fit of tracklets
2349 // 11. Improve seeding tracklets quality by reassigning clusters.
2350 // See AliTRDtrackerV1::ImproveSeedQuality() for details.
2351 // 12. Helix fit of all 6 seeding tracklets and chi2 calculation
2352 // 13. Hyperplane fit and track quality calculation. See AliTRDtrackerFitter::FitHyperplane() for details.
2353 // 14. Cooking labels for tracklets. Should be done only for MC
2354 // 15. Register seeds.
2357 AliTRDtrackingChamber *chamber = 0x0;
2358 AliTRDcluster *c[kNSeedPlanes] = {0x0, 0x0, 0x0, 0x0}; // initilize seeding clusters
2359 AliTRDseedV1 *cseed = &sseed[0]; // initialize tracklets for first track
2360 Int_t ncl, mcl; // working variable for looping over clusters
2361 Int_t index[AliTRDchamberTimeBin::kMaxClustersLayer], jndex[AliTRDchamberTimeBin::kMaxClustersLayer];
2363 // chi2[0] = tracklet chi2 on the Z direction
2364 // chi2[1] = tracklet chi2 on the R direction
2367 // this should be data member of AliTRDtrack
2368 Double_t seedQuality[kMaxTracksStack];
2370 // unpack control parameters
2371 Int_t config = ipar[0];
2372 Int_t ntracks = ipar[1];
2373 Int_t istack = ipar[2];
2374 Int_t planes[kNSeedPlanes]; GetSeedingConfig(config, planes);
2375 Int_t planesExt[kNPlanes-kNSeedPlanes]; GetExtrapolationConfig(config, planesExt);
2378 // Init chambers geometry
2379 Double_t hL[kNPlanes]; // Tilting angle
2380 Float_t padlength[kNPlanes]; // pad lenghts
2381 Float_t padwidth[kNPlanes]; // pad widths
2382 AliTRDpadPlane *pp = 0x0;
2383 for(int iplane=0; iplane<kNPlanes; iplane++){
2384 pp = fGeom->GetPadPlane(iplane, istack);
2385 hL[iplane] = TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle());
2386 padlength[iplane] = pp->GetLengthIPad();
2387 padwidth[iplane] = pp->GetWidthIPad();
2390 // Init anode wire position for chambers
2391 Double_t x0[kNPlanes], // anode wire position
2392 driftLength = .5*AliTRDgeometry::AmThick() - AliTRDgeometry::DrThick(); // drift length
2393 TGeoHMatrix *matrix = 0x0;
2394 Double_t loc[] = {AliTRDgeometry::AnodePos(), 0., 0.};
2395 Double_t glb[] = {0., 0., 0.};
2396 AliTRDtrackingChamber **cIter = &stack[0];
2397 for(int iLayer=0; iLayer<kNPlanes; iLayer++,cIter++){
2398 if(!(*cIter)) continue;
2399 if(!(matrix = fGeom->GetClusterMatrix((*cIter)->GetDetector()))){
2401 x0[iLayer] = fgkX0[iLayer];
2403 matrix->LocalToMaster(loc, glb);
2404 x0[iLayer] = glb[0];
2407 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 10){
2408 AliInfo(Form("Making seeds Stack[%d] Config[%d] Tracks[%d]...", istack, config, ntracks));
2411 // Build seeding layers
2414 for(int isl=0; isl<kNSeedPlanes; isl++){
2415 if(!(chamber = stack[planes[isl]])) continue;
2416 if(!chamber->GetSeedingLayer(fSeedTB[isl], fGeom, fReconstructor)) continue;
2419 if(nlayers < kNSeedPlanes) return ntracks;
2422 // Start finding seeds
2423 Double_t cond0[4], cond1[4], cond2[4];
2425 while((c[3] = (*fSeedTB[3])[icl++])){
2427 fSeedTB[0]->BuildCond(c[3], cond0, 0);
2428 fSeedTB[0]->GetClusters(cond0, index, ncl);
2429 //printf("Found c[3] candidates 0 %d\n", ncl);
2432 c[0] = (*fSeedTB[0])[index[jcl++]];
2434 Double_t dx = c[3]->GetX() - c[0]->GetX();
2435 Double_t theta = (c[3]->GetZ() - c[0]->GetZ())/dx;
2436 Double_t phi = (c[3]->GetY() - c[0]->GetY())/dx;
2437 fSeedTB[1]->BuildCond(c[0], cond1, 1, theta, phi);
2438 fSeedTB[1]->GetClusters(cond1, jndex, mcl);
2439 //printf("Found c[0] candidates 1 %d\n", mcl);
2443 c[1] = (*fSeedTB[1])[jndex[kcl++]];
2445 fSeedTB[2]->BuildCond(c[1], cond2, 2, theta, phi);
2446 c[2] = fSeedTB[2]->GetNearestCluster(cond2);
2447 //printf("Found c[1] candidate 2 %p\n", c[2]);
2450 // AliInfo("Seeding clusters found. Building seeds ...");
2451 // for(Int_t i = 0; i < kNSeedPlanes; i++) printf("%i. coordinates: x = %6.3f, y = %6.3f, z = %6.3f\n", i, c[i]->GetX(), c[i]->GetY(), c[i]->GetZ());
2453 for (Int_t il = 0; il < kNPlanes; il++) cseed[il].Reset();
2457 AliTRDseedV1 *tseed = &cseed[0];
2459 for(int iLayer=0; iLayer<kNPlanes; iLayer++, tseed++, cIter++){
2460 Int_t det = (*cIter) ? (*cIter)->GetDetector() : -1;
2461 tseed->SetDetector(det);
2462 tseed->SetTilt(hL[iLayer]);
2463 tseed->SetPadLength(padlength[iLayer]);
2464 tseed->SetPadWidth(padwidth[iLayer]);
2465 tseed->SetReconstructor(fReconstructor);
2466 tseed->SetX0(det<0 ? fR[iLayer]+driftLength : x0[iLayer]);
2467 tseed->Init(GetRiemanFitter());
2468 tseed->SetStandAlone(kTRUE);
2471 Bool_t isFake = kFALSE;
2472 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2473 if (c[0]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2474 if (c[1]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2475 if (c[2]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2478 for(Int_t l = 0; l < kNSeedPlanes; l++) xpos[l] = fSeedTB[l]->GetX();
2480 for(int il=0; il<4; il++) yref[il] = cseed[planes[il]].GetYref(0);
2481 Int_t ll = c[3]->GetLabel(0);
2482 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2483 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2484 AliRieman *rim = GetRiemanFitter();
2485 TTreeSRedirector &cs0 = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2487 <<"EventNumber=" << eventNumber
2488 <<"CandidateNumber=" << candidateNumber
2489 <<"isFake=" << isFake
2490 <<"config=" << config
2492 <<"chi2z=" << chi2[0]
2493 <<"chi2y=" << chi2[1]
2494 <<"Y2exp=" << cond2[0]
2495 <<"Z2exp=" << cond2[1]
2496 <<"X0=" << xpos[0] //layer[sLayer]->GetX()
2497 <<"X1=" << xpos[1] //layer[sLayer + 1]->GetX()
2498 <<"X2=" << xpos[2] //layer[sLayer + 2]->GetX()
2499 <<"X3=" << xpos[3] //layer[sLayer + 3]->GetX()
2500 <<"yref0=" << yref[0]
2501 <<"yref1=" << yref[1]
2502 <<"yref2=" << yref[2]
2503 <<"yref3=" << yref[3]
2508 <<"Seed0.=" << &cseed[planes[0]]
2509 <<"Seed1.=" << &cseed[planes[1]]
2510 <<"Seed2.=" << &cseed[planes[2]]
2511 <<"Seed3.=" << &cseed[planes[3]]
2512 <<"RiemanFitter.=" << rim
2515 if(chi2[0] > fReconstructor->GetRecoParam() ->GetChi2Z()/*7./(3. - sLayer)*//*iter*/){
2516 // //AliInfo(Form("Failed chi2 filter on chi2Z [%f].", chi2[0]));
2517 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2520 if(chi2[1] > fReconstructor->GetRecoParam() ->GetChi2Y()/*1./(3. - sLayer)*//*iter*/){
2521 // //AliInfo(Form("Failed chi2 filter on chi2Y [%f].", chi2[1]));
2522 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2525 //AliInfo("Passed chi2 filter.");
2527 // try attaching clusters to tracklets
2529 for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
2530 Int_t jLayer = planes[iLayer];
2531 if(!cseed[jLayer].AttachClusters(stack[jLayer], kTRUE)) continue;
2532 cseed[jLayer].UpdateUsed();
2533 if(!cseed[jLayer].IsOK()) continue;
2537 if(mlayers < kNSeedPlanes){
2538 //AliInfo(Form("Failed updating all seeds %d [%d].", mlayers, kNSeedPlanes));
2539 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2543 // temporary exit door for the HLT
2544 if(fReconstructor->IsHLT()){
2545 // attach clusters to extrapolation chambers
2546 for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
2547 Int_t jLayer = planesExt[iLayer];
2548 if(!(chamber = stack[jLayer])) continue;
2549 cseed[jLayer].AttachClusters(chamber, kTRUE);
2551 fTrackQuality[ntracks] = 1.; // dummy value
2553 if(ntracks == kMaxTracksStack) return ntracks;
2559 // Update Seeds and calculate Likelihood
2560 // fit tracklets and cook likelihood
2561 FitTiltedRieman(&cseed[0], kTRUE);
2562 for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
2563 Int_t jLayer = planes[iLayer];
2564 cseed[jLayer].Fit(kTRUE);
2566 Double_t like = CookLikelihood(&cseed[0], planes); // to be checked
2568 if (TMath::Log(1.E-9 + like) < fReconstructor->GetRecoParam() ->GetTrackLikelihood()){
2569 //AliInfo(Form("Failed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2570 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2573 //AliInfo(Form("Passed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2575 // book preliminary results
2576 seedQuality[ntracks] = like;
2577 fSeedLayer[ntracks] = config;/*sLayer;*/
2579 // attach clusters to the extrapolation seeds
2580 for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
2581 Int_t jLayer = planesExt[iLayer];
2582 if(!(chamber = stack[jLayer])) continue;
2584 // fit extrapolated seed
2585 if ((jLayer == 0) && !(cseed[1].IsOK())) continue;
2586 if ((jLayer == 5) && !(cseed[4].IsOK())) continue;
2587 AliTRDseedV1 pseed = cseed[jLayer];
2588 if(!pseed.AttachClusters(chamber, kTRUE)) continue;
2590 cseed[jLayer] = pseed;
2591 FitTiltedRieman(cseed, kTRUE);
2592 cseed[jLayer].Fit(kTRUE);
2595 // AliInfo("Extrapolation done.");
2596 // Debug Stream containing all the 6 tracklets
2597 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2598 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2599 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2600 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2601 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2602 cstreamer << "MakeSeeds1"
2603 << "EventNumber=" << eventNumber
2604 << "CandidateNumber=" << candidateNumber
2605 << "S0.=" << &cseed[0]
2606 << "S1.=" << &cseed[1]
2607 << "S2.=" << &cseed[2]
2608 << "S3.=" << &cseed[3]
2609 << "S4.=" << &cseed[4]
2610 << "S5.=" << &cseed[5]
2611 << "FitterT.=" << tiltedRieman
2615 if(fReconstructor->GetRecoParam()->HasImproveTracklets() && ImproveSeedQuality(stack, cseed) < 4){
2616 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2619 //AliInfo("Improve seed quality done.");
2621 // fit full track and cook likelihoods
2622 // Double_t curv = FitRieman(&cseed[0], chi2);
2623 // Double_t chi2ZF = chi2[0] / TMath::Max((mlayers - 3.), 1.);
2624 // Double_t chi2RF = chi2[1] / TMath::Max((mlayers - 3.), 1.);
2626 // do the final track fitting (Once with vertex constraint and once without vertex constraint)
2627 Double_t chi2Vals[3];
2628 chi2Vals[0] = FitTiltedRieman(&cseed[0], kFALSE);
2629 if(fReconstructor->GetRecoParam()->IsVertexConstrained())
2630 chi2Vals[1] = FitTiltedRiemanConstraint(&cseed[0], GetZ()); // Do Vertex Constrained fit if desired
2633 chi2Vals[2] = GetChi2Z(&cseed[0]) / TMath::Max((mlayers - 3.), 1.);
2634 // Chi2 definitions in testing stage
2635 //chi2Vals[2] = GetChi2ZTest(&cseed[0]);
2636 fTrackQuality[ntracks] = CalculateTrackLikelihood(&cseed[0], &chi2Vals[0]);
2637 //AliInfo("Hyperplane fit done\n");
2639 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2640 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2641 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2642 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2643 TLinearFitter *fitterTC = GetTiltedRiemanFitterConstraint();
2644 TLinearFitter *fitterT = GetTiltedRiemanFitter();
2646 for(Int_t iseed = 0; iseed < kNPlanes; iseed++){
2647 ncls += cseed[iseed].IsOK() ? cseed[iseed].GetN2() : 0;
2649 cstreamer << "MakeSeeds2"
2650 << "EventNumber=" << eventNumber
2651 << "CandidateNumber=" << candidateNumber
2652 << "Chi2TR=" << chi2Vals[0]
2653 << "Chi2TC=" << chi2Vals[1]
2654 << "Nlayers=" << mlayers
2655 << "NClusters=" << ncls
2657 << "S0.=" << &cseed[0]
2658 << "S1.=" << &cseed[1]
2659 << "S2.=" << &cseed[2]
2660 << "S3.=" << &cseed[3]
2661 << "S4.=" << &cseed[4]
2662 << "S5.=" << &cseed[5]
2663 << "FitterT.=" << fitterT
2664 << "FitterTC.=" << fitterTC
2669 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2670 if(ntracks == kMaxTracksStack){
2671 AliWarning(Form("Number of seeds reached maximum allowed (%d) in stack.", kMaxTracksStack));
2682 //_____________________________________________________________________________
2683 AliTRDtrackV1* AliTRDtrackerV1::MakeTrack(AliTRDseedV1 *seeds, Double_t *params)
2686 // Build a TRD track out of tracklet candidates
2689 // seeds : array of tracklets
2690 // params : track parameters (see MakeSeeds() function body for a detailed description)
2695 // Detailed description
2697 // To be discussed with Marian !!
2701 Double_t alpha = AliTRDgeometry::GetAlpha();
2702 Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
2706 c[ 1] = 0.0; c[ 2] = 2.0;
2707 c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
2708 c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
2709 c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
2711 AliTRDtrackV1 track(seeds, ¶ms[1], c, params[0], params[6]*alpha+shift);
2712 track.PropagateTo(params[0]-5.0);
2713 AliTRDseedV1 *ptrTracklet = 0x0;
2715 for (Int_t jLayer = 0; jLayer < AliTRDgeometry::kNlayer; jLayer++) {
2716 ptrTracklet = &seeds[jLayer];
2717 if(!ptrTracklet->IsOK()) continue;
2718 if(TMath::Abs(ptrTracklet->GetYref(1) - ptrTracklet->GetYfit(1)) >= .2) continue; // check this condition with Marian
2721 if(fReconstructor->IsHLT()){
2722 for(Int_t ip=0; ip<kNPlanes; ip++){
2723 track.UnsetTracklet(ip);
2724 ptrTracklet = SetTracklet(&seeds[ip]);
2725 ptrTracklet->UseClusters();
2726 track.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
2728 AliTRDtrackV1 *ptrTrack = SetTrack(&track);
2729 ptrTrack->SetReconstructor(fReconstructor);
2733 track.ResetCovariance(1);
2734 Int_t nc = TMath::Abs(FollowBackProlongation(track));
2735 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 5){
2736 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2737 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2738 Double_t p[5]; // Track Params for the Debug Stream
2739 track.GetExternalParameters(params[0], p);
2740 TTreeSRedirector &cs = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2742 << "EventNumber=" << eventNumber
2743 << "CandidateNumber=" << candidateNumber
2745 << "X=" << params[0]
2751 << "Yin=" << params[1]
2752 << "Zin=" << params[2]
2753 << "snpin=" << params[3]
2754 << "tndin=" << params[4]
2755 << "crvin=" << params[5]
2756 << "track.=" << &track
2759 if (nc < 30) return 0x0;
2761 AliTRDtrackV1 *ptrTrack = SetTrack(&track);
2762 ptrTrack->SetReconstructor(fReconstructor);
2763 ptrTrack->CookLabel(.9);
2765 // computes PID for track
2766 ptrTrack->CookPID();
2767 // update calibration references using this track
2768 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
2770 AliInfo("Could not get Calibra instance\n");
2771 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(ptrTrack);
2777 //____________________________________________________________________
2778 Int_t AliTRDtrackerV1::ImproveSeedQuality(AliTRDtrackingChamber **stack, AliTRDseedV1 *cseed)
2781 // Sort tracklets according to "quality" and try to "improve" the first 4 worst
2784 // layers : Array of propagation layers for a stack/supermodule
2785 // cseed : Array of 6 seeding tracklets which has to be improved
2788 // cssed : Improved seeds
2790 // Detailed description
2792 // Iterative procedure in which new clusters are searched for each
2793 // tracklet seed such that the seed quality (see AliTRDseed::GetQuality())
2794 // can be maximized. If some optimization is found the old seeds are replaced.
2799 // make a local working copy
2800 AliTRDtrackingChamber *chamber = 0x0;
2801 AliTRDseedV1 bseed[6];
2803 for (Int_t jLayer = 0; jLayer < 6; jLayer++) bseed[jLayer] = cseed[jLayer];
2805 Float_t lastquality = 10000.0;
2806 Float_t lastchi2 = 10000.0;
2807 Float_t chi2 = 1000.0;
2809 for (Int_t iter = 0; iter < 4; iter++) {
2810 Float_t sumquality = 0.0;
2811 Float_t squality[6];
2812 Int_t sortindexes[6];
2814 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2815 squality[jLayer] = bseed[jLayer].IsOK() ? bseed[jLayer].GetQuality(kTRUE) : 1000.;
2816 sumquality += squality[jLayer];
2818 if ((sumquality >= lastquality) || (chi2 > lastchi2)) break;
2821 lastquality = sumquality;
2823 if (iter > 0) for (Int_t jLayer = 0; jLayer < 6; jLayer++) cseed[jLayer] = bseed[jLayer];
2825 TMath::Sort(6, squality, sortindexes, kFALSE);
2826 for (Int_t jLayer = 5; jLayer > 1; jLayer--) {
2827 Int_t bLayer = sortindexes[jLayer];
2828 if(!(chamber = stack[bLayer])) continue;
2829 bseed[bLayer].AttachClusters(chamber, kTRUE);
2830 bseed[bLayer].Fit(kTRUE);
2831 if(bseed[bLayer].IsOK()) nLayers++;
2834 chi2 = FitTiltedRieman(bseed, kTRUE);
2835 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 7){
2836 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2837 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2838 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2839 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2840 cstreamer << "ImproveSeedQuality"
2841 << "EventNumber=" << eventNumber
2842 << "CandidateNumber=" << candidateNumber
2843 << "Iteration=" << iter
2844 << "S0.=" << &bseed[0]
2845 << "S1.=" << &bseed[1]
2846 << "S2.=" << &bseed[2]
2847 << "S3.=" << &bseed[3]
2848 << "S4.=" << &bseed[4]
2849 << "S5.=" << &bseed[5]
2850 << "FitterT.=" << tiltedRieman
2854 // we are sure that at least 2 tracklets are OK !
2858 //_________________________________________________________________________
2859 Double_t AliTRDtrackerV1::CalculateTrackLikelihood(AliTRDseedV1 *tracklets, Double_t *chi2){
2861 // Calculates the Track Likelihood value. This parameter serves as main quality criterion for
2862 // the track selection
2863 // The likelihood value containes:
2864 // - The chi2 values from the both fitters and the chi2 values in z-direction from a linear fit
2865 // - The Sum of the Parameter |slope_ref - slope_fit|/Sigma of the tracklets
2866 // For all Parameters an exponential dependency is used
2868 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
2869 // - Array of chi2 values:
2870 // * Non-Constrained Tilted Riemann fit
2871 // * Vertex-Constrained Tilted Riemann fit
2872 // * z-Direction from Linear fit
2873 // Output: - The calculated track likelihood
2878 Double_t chi2phi = 0, nLayers = 0;
2879 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
2880 if(!tracklets[iLayer].IsOK()) continue;
2881 chi2phi += tracklets[iLayer].GetChi2Phi();
2884 chi2phi /= Float_t (nLayers - 2.0);
2886 Double_t likeChi2Z = TMath::Exp(-chi2[2] * 0.14); // Chi2Z
2887 Double_t likeChi2TC = (fReconstructor->GetRecoParam() ->IsVertexConstrained()) ?
2888 TMath::Exp(-chi2[1] * 0.677) : 1; // Constrained Tilted Riemann
2889 Double_t likeChi2TR = TMath::Exp(-chi2[0] * 0.78); // Non-constrained Tilted Riemann
2890 Double_t likeChi2Phi= TMath::Exp(-chi2phi * 3.23);
2891 Double_t trackLikelihood = likeChi2Z * likeChi2TR * likeChi2Phi;
2893 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2894 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2895 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2896 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2897 cstreamer << "CalculateTrackLikelihood0"
2898 << "EventNumber=" << eventNumber
2899 << "CandidateNumber=" << candidateNumber
2900 << "LikeChi2Z=" << likeChi2Z
2901 << "LikeChi2TR=" << likeChi2TR
2902 << "LikeChi2TC=" << likeChi2TC
2903 << "LikeChi2Phi=" << likeChi2Phi
2904 << "TrackLikelihood=" << trackLikelihood
2908 return trackLikelihood;
2911 //____________________________________________________________________
2912 Double_t AliTRDtrackerV1::CookLikelihood(AliTRDseedV1 *cseed, Int_t planes[4])
2915 // Calculate the probability of this track candidate.
2918 // cseeds : array of candidate tracklets
2919 // planes : array of seeding planes (see seeding configuration)
2920 // chi2 : chi2 values (on the Z and Y direction) from the rieman fit of the track.
2925 // Detailed description
2927 // The track quality is estimated based on the following 4 criteria:
2928 // 1. precision of the rieman fit on the Y direction (likea)
2929 // 2. chi2 on the Y direction (likechi2y)
2930 // 3. chi2 on the Z direction (likechi2z)
2931 // 4. number of attached clusters compared to a reference value
2932 // (see AliTRDrecoParam::fkFindable) (likeN)
2934 // The distributions for each type of probabilities are given below as of
2935 // (date). They have to be checked to assure consistency of estimation.
2938 // ratio of the total number of clusters/track which are expected to be found by the tracker.
2939 const AliTRDrecoParam *fRecoPars = fReconstructor->GetRecoParam();
2941 Double_t chi2y = GetChi2Y(&cseed[0]);
2942 Double_t chi2z = GetChi2Z(&cseed[0]);
2944 Float_t nclusters = 0.;
2945 Double_t sumda = 0.;
2946 for(UChar_t ilayer = 0; ilayer < 4; ilayer++){
2947 Int_t jlayer = planes[ilayer];
2948 nclusters += cseed[jlayer].GetN2();
2949 sumda += TMath::Abs(cseed[jlayer].GetYfit(1) - cseed[jlayer].GetYref(1));
2953 Double_t likea = TMath::Exp(-sumda * fRecoPars->GetPhiSlope());
2954 Double_t likechi2y = 0.0000000001;
2955 if (fReconstructor->IsCosmic() || chi2y < fRecoPars->GetChi2YCut()) likechi2y += TMath::Exp(-TMath::Sqrt(chi2y) * fRecoPars->GetChi2YSlope());
2956 Double_t likechi2z = TMath::Exp(-chi2z * fRecoPars->GetChi2ZSlope());
2957 Double_t likeN = TMath::Exp(-(fRecoPars->GetNMeanClusters() - nclusters) / fRecoPars->GetNSigmaClusters());
2958 Double_t like = likea * likechi2y * likechi2z * likeN;
2960 // AliInfo(Form("sumda(%f) chi2[0](%f) chi2[1](%f) likea(%f) likechi2y(%f) likechi2z(%f) nclusters(%d) likeN(%f)", sumda, chi2[0], chi2[1], likea, likechi2y, likechi2z, nclusters, likeN));
2961 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2962 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2963 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2964 Int_t nTracklets = 0; Float_t mean_ncls = 0;
2965 for(Int_t iseed=0; iseed < kNPlanes; iseed++){
2966 if(!cseed[iseed].IsOK()) continue;
2968 mean_ncls += cseed[iseed].GetN2();
2970 if(nTracklets) mean_ncls /= nTracklets;
2971 // The Debug Stream contains the seed
2972 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2973 cstreamer << "CookLikelihood"
2974 << "EventNumber=" << eventNumber
2975 << "CandidateNumber=" << candidateNumber
2976 << "tracklet0.=" << &cseed[0]
2977 << "tracklet1.=" << &cseed[1]
2978 << "tracklet2.=" << &cseed[2]
2979 << "tracklet3.=" << &cseed[3]
2980 << "tracklet4.=" << &cseed[4]
2981 << "tracklet5.=" << &cseed[5]
2982 << "sumda=" << sumda
2983 << "chi2y=" << chi2y
2984 << "chi2z=" << chi2z
2985 << "likea=" << likea
2986 << "likechi2y=" << likechi2y
2987 << "likechi2z=" << likechi2z
2988 << "nclusters=" << nclusters
2989 << "likeN=" << likeN
2991 << "meanncls=" << mean_ncls
2998 //____________________________________________________________________
2999 void AliTRDtrackerV1::GetSeedingConfig(Int_t iconfig, Int_t planes[4])
3002 // Map seeding configurations to detector planes.
3005 // iconfig : configuration index
3006 // planes : member planes of this configuration. On input empty.
3009 // planes : contains the planes which are defining the configuration
3011 // Detailed description
3013 // Here is the list of seeding planes configurations together with
3014 // their topological classification:
3032 // The topologic quality is modeled as follows:
3033 // 1. The general model is define by the equation:
3034 // p(conf) = exp(-conf/2)
3035 // 2. According to the topologic classification, configurations from the same
3036 // class are assigned the agerage value over the model values.
3037 // 3. Quality values are normalized.
3039 // The topologic quality distribution as function of configuration is given below:
3041 // <img src="gif/topologicQA.gif">
3046 case 0: // 5432 TQ 0
3052 case 1: // 4321 TQ 0
3058 case 2: // 3210 TQ 0
3064 case 3: // 5321 TQ 1
3070 case 4: // 4210 TQ 1
3076 case 5: // 5431 TQ 1
3082 case 6: // 4320 TQ 1
3088 case 7: // 5430 TQ 2
3094 case 8: // 5210 TQ 2
3100 case 9: // 5421 TQ 3
3106 case 10: // 4310 TQ 3
3112 case 11: // 5410 TQ 4
3118 case 12: // 5420 TQ 5
3124 case 13: // 5320 TQ 5
3130 case 14: // 5310 TQ 5
3139 //____________________________________________________________________
3140 void AliTRDtrackerV1::GetExtrapolationConfig(Int_t iconfig, Int_t planes[2])
3143 // Returns the extrapolation planes for a seeding configuration.
3146 // iconfig : configuration index
3147 // planes : planes which are not in this configuration. On input empty.
3150 // planes : contains the planes which are not in the configuration
3152 // Detailed description
3156 case 0: // 5432 TQ 0
3160 case 1: // 4321 TQ 0
3164 case 2: // 3210 TQ 0
3168 case 3: // 5321 TQ 1
3172 case 4: // 4210 TQ 1
3176 case 5: // 5431 TQ 1
3180 case 6: // 4320 TQ 1
3184 case 7: // 5430 TQ 2
3188 case 8: // 5210 TQ 2
3192 case 9: // 5421 TQ 3
3196 case 10: // 4310 TQ 3
3200 case 11: // 5410 TQ 4
3204 case 12: // 5420 TQ 5
3208 case 13: // 5320 TQ 5
3212 case 14: // 5310 TQ 5
3219 //____________________________________________________________________
3220 AliCluster* AliTRDtrackerV1::GetCluster(Int_t idx) const
3222 Int_t ncls = fClusters->GetEntriesFast();
3223 return idx >= 0 && idx < ncls ? (AliCluster*)fClusters->UncheckedAt(idx) : 0x0;
3226 //____________________________________________________________________
3227 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(Int_t idx) const
3229 Int_t ntrklt = fTracklets->GetEntriesFast();
3230 return idx >= 0 && idx < ntrklt ? (AliTRDseedV1*)fTracklets->UncheckedAt(idx) : 0x0;
3233 //____________________________________________________________________
3234 AliKalmanTrack* AliTRDtrackerV1::GetTrack(Int_t idx) const
3236 Int_t ntrk = fTracks->GetEntriesFast();
3237 return idx >= 0 && idx < ntrk ? (AliKalmanTrack*)fTracks->UncheckedAt(idx) : 0x0;
3240 //____________________________________________________________________
3241 Float_t AliTRDtrackerV1::CalculateReferenceX(AliTRDseedV1 *tracklets){
3243 // Calculates the reference x-position for the tilted Rieman fit defined as middle
3244 // of the stack (middle between layers 2 and 3). For the calculation all the tracklets
3245 // are taken into account
3247 // Parameters: - Array of tracklets(AliTRDseedV1)
3249 // Output: - The reference x-position(Float_t)
3251 Int_t nDistances = 0;
3252 Float_t meanDistance = 0.;
3253 Int_t startIndex = 5;
3254 for(Int_t il =5; il > 0; il--){
3255 if(tracklets[il].IsOK() && tracklets[il -1].IsOK()){
3256 Float_t xdiff = tracklets[il].GetX0() - tracklets[il -1].GetX0();
3257 meanDistance += xdiff;
3260 if(tracklets[il].IsOK()) startIndex = il;
3262 if(tracklets[0].IsOK()) startIndex = 0;
3264 // We should normally never get here
3265 Float_t xpos[2]; memset(xpos, 0, sizeof(Float_t) * 2);
3266 Int_t iok = 0, idiff = 0;
3267 // This attempt is worse and should be avoided:
3268 // check for two chambers which are OK and repeat this without taking the mean value
3269 // Strategy avoids a division by 0;
3270 for(Int_t il = 5; il >= 0; il--){
3271 if(tracklets[il].IsOK()){
3272 xpos[iok] = tracklets[il].GetX0();
3276 if(iok) idiff++; // to get the right difference;
3280 meanDistance = (xpos[0] - xpos[1])/idiff;
3283 // we have do not even have 2 layers which are OK? The we do not need to fit at all
3288 meanDistance /= nDistances;
3290 return tracklets[startIndex].GetX0() + (2.5 - startIndex) * meanDistance - 0.5 * (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
3293 // //_____________________________________________________________________________
3294 // Int_t AliTRDtrackerV1::Freq(Int_t n, const Int_t *inlist
3295 // , Int_t *outlist, Bool_t down)
3298 // // Sort eleements according occurancy
3299 // // The size of output array has is 2*n
3306 // Int_t *sindexS = new Int_t[n]; // Temporary array for sorting
3307 // Int_t *sindexF = new Int_t[2*n];
3308 // for (Int_t i = 0; i < n; i++) {
3312 // TMath::Sort(n,inlist,sindexS,down);
3314 // Int_t last = inlist[sindexS[0]];
3315 // Int_t val = last;
3317 // sindexF[0+n] = last;
3318 // Int_t countPos = 0;
3320 // // Find frequency
3321 // for (Int_t i = 1; i < n; i++) {
3322 // val = inlist[sindexS[i]];
3323 // if (last == val) {
3324 // sindexF[countPos]++;
3328 // sindexF[countPos+n] = val;
3329 // sindexF[countPos]++;
3333 // if (last == val) {
3337 // // Sort according frequency
3338 // TMath::Sort(countPos,sindexF,sindexS,kTRUE);
3340 // for (Int_t i = 0; i < countPos; i++) {
3341 // outlist[2*i ] = sindexF[sindexS[i]+n];
3342 // outlist[2*i+1] = sindexF[sindexS[i]];
3345 // delete [] sindexS;
3346 // delete [] sindexF;
3353 //____________________________________________________________________
3354 void AliTRDtrackerV1::ResetSeedTB()
3356 // reset buffer for seeding time bin layers. If the time bin
3357 // layers are not allocated this function allocates them
3359 for(Int_t isl=0; isl<kNSeedPlanes; isl++){
3360 if(!fSeedTB[isl]) fSeedTB[isl] = new AliTRDchamberTimeBin();
3361 else fSeedTB[isl]->Clear();
3366 //_____________________________________________________________________________
3367 Float_t AliTRDtrackerV1::GetChi2Y(AliTRDseedV1 *tracklets) const
3369 // Calculates normalized chi2 in y-direction
3370 // chi2 = Sum chi2 / n_tracklets
3372 Double_t chi2 = 0.; Int_t n = 0;
3373 for(Int_t ipl = kNPlanes; ipl--;){
3374 if(!tracklets[ipl].IsOK()) continue;
3375 chi2 += tracklets[ipl].GetChi2Y();
3378 return n ? chi2/n : 0.;
3381 //_____________________________________________________________________________
3382 Float_t AliTRDtrackerV1::GetChi2Z(AliTRDseedV1 *tracklets) const
3384 // Calculates normalized chi2 in z-direction
3385 // chi2 = Sum chi2 / n_tracklets
3387 Double_t chi2 = 0; Int_t n = 0;
3388 for(Int_t ipl = kNPlanes; ipl--;){
3389 if(!tracklets[ipl].IsOK()) continue;
3390 chi2 += tracklets[ipl].GetChi2Z();
3393 return n ? chi2/n : 0.;
3396 ///////////////////////////////////////////////////////
3398 // Resources of class AliTRDLeastSquare //
3400 ///////////////////////////////////////////////////////
3402 //_____________________________________________________________________________
3403 AliTRDtrackerV1::AliTRDLeastSquare::AliTRDLeastSquare(){
3405 // Constructor of the nested class AliTRDtrackFitterLeastSquare
3407 memset(fParams, 0, sizeof(Double_t) * 2);
3408 memset(fSums, 0, sizeof(Double_t) * 5);
3409 memset(fCovarianceMatrix, 0, sizeof(Double_t) * 3);
3413 //_____________________________________________________________________________
3414 void AliTRDtrackerV1::AliTRDLeastSquare::AddPoint(Double_t *x, Double_t y, Double_t sigmaY){
3416 // Adding Point to the fitter
3418 Double_t weight = 1/(sigmaY * sigmaY);
3420 // printf("Adding point x = %f, y = %f, sigma = %f\n", xpt, y, sigmaY);
3422 fSums[1] += weight * xpt;
3423 fSums[2] += weight * y;
3424 fSums[3] += weight * xpt * y;
3425 fSums[4] += weight * xpt * xpt;
3426 fSums[5] += weight * y * y;
3429 //_____________________________________________________________________________
3430 void AliTRDtrackerV1::AliTRDLeastSquare::RemovePoint(Double_t *x, Double_t y, Double_t sigmaY){
3432 // Remove Point from the sample
3434 Double_t weight = 1/(sigmaY * sigmaY);
3437 fSums[1] -= weight * xpt;
3438 fSums[2] -= weight * y;
3439 fSums[3] -= weight * xpt * y;
3440 fSums[4] -= weight * xpt * xpt;
3441 fSums[5] -= weight * y * y;
3444 //_____________________________________________________________________________
3445 void AliTRDtrackerV1::AliTRDLeastSquare::Eval(){
3447 // Evaluation of the fit:
3448 // Calculation of the parameters
3449 // Calculation of the covariance matrix
3452 Double_t denominator = fSums[0] * fSums[4] - fSums[1] *fSums[1];
3453 if(denominator==0) return;
3455 // for(Int_t isum = 0; isum < 5; isum++)
3456 // printf("fSums[%d] = %f\n", isum, fSums[isum]);
3457 // printf("denominator = %f\n", denominator);
3458 fParams[0] = (fSums[2] * fSums[4] - fSums[1] * fSums[3])/ denominator;
3459 fParams[1] = (fSums[0] * fSums[3] - fSums[1] * fSums[2]) / denominator;
3460 // printf("fParams[0] = %f, fParams[1] = %f\n", fParams[0], fParams[1]);
3462 // Covariance matrix
3463 fCovarianceMatrix[0] = fSums[4] - fSums[1] * fSums[1] / fSums[0];
3464 fCovarianceMatrix[1] = fSums[5] - fSums[2] * fSums[2] / fSums[0];
3465 fCovarianceMatrix[2] = fSums[3] - fSums[1] * fSums[2] / fSums[0];
3468 //_____________________________________________________________________________
3469 Double_t AliTRDtrackerV1::AliTRDLeastSquare::GetFunctionValue(Double_t *xpos) const {
3471 // Returns the Function value of the fitted function at a given x-position
3473 return fParams[0] + fParams[1] * (*xpos);
3476 //_____________________________________________________________________________
3477 void AliTRDtrackerV1::AliTRDLeastSquare::GetCovarianceMatrix(Double_t *storage) const {
3479 // Copies the values of the covariance matrix into the storage
3481 memcpy(storage, fCovarianceMatrix, sizeof(Double_t) * 3);