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>
40 #include "AliESDEvent.h"
41 #include "AliGeomManager.h"
42 #include "AliRieman.h"
43 #include "AliTrackPointArray.h"
45 #include "AliTRDgeometry.h"
46 #include "AliTRDpadPlane.h"
47 #include "AliTRDcalibDB.h"
48 #include "AliTRDReconstructor.h"
49 #include "AliTRDCalibraFillHisto.h"
50 #include "AliTRDrecoParam.h"
52 #include "AliTRDcluster.h"
53 #include "AliTRDseedV1.h"
54 #include "AliTRDtrackV1.h"
55 #include "AliTRDtrackerV1.h"
56 #include "AliTRDtrackerDebug.h"
57 #include "AliTRDtrackingChamber.h"
58 #include "AliTRDchamberTimeBin.h"
62 ClassImp(AliTRDtrackerV1)
65 const Float_t AliTRDtrackerV1::fgkMinClustersInTrack = 0.5; //
66 const Float_t AliTRDtrackerV1::fgkLabelFraction = 0.8; //
67 const Double_t AliTRDtrackerV1::fgkMaxChi2 = 12.0; //
68 const Double_t AliTRDtrackerV1::fgkMaxSnp = 0.95; // Maximum local sine of the azimuthal angle
69 const Double_t AliTRDtrackerV1::fgkMaxStep = 2.0; // Maximal step size in propagation
70 Double_t AliTRDtrackerV1::fgTopologicQA[kNConfigs] = {
71 0.1112, 0.1112, 0.1112, 0.0786, 0.0786,
72 0.0786, 0.0786, 0.0579, 0.0579, 0.0474,
73 0.0474, 0.0408, 0.0335, 0.0335, 0.0335
75 Int_t AliTRDtrackerV1::fgNTimeBins = 0;
76 AliRieman* AliTRDtrackerV1::fgRieman = 0x0;
77 TLinearFitter* AliTRDtrackerV1::fgTiltedRieman = 0x0;
78 TLinearFitter* AliTRDtrackerV1::fgTiltedRiemanConstrained = 0x0;
80 //____________________________________________________________________
81 AliTRDtrackerV1::AliTRDtrackerV1(AliTRDReconstructor *rec)
84 ,fGeom(new AliTRDgeometry())
91 // Default constructor.
93 AliTRDcalibDB *trd = 0x0;
94 if (!(trd = AliTRDcalibDB::Instance())) {
95 AliFatal("Could not get calibration object");
98 if(!fgNTimeBins) fgNTimeBins = trd->GetNumberOfTimeBins();
100 for (Int_t isector = 0; isector < AliTRDgeometry::kNsector; isector++) new(&fTrSec[isector]) AliTRDtrackingSector(fGeom, isector);
102 for(Int_t isl =0; isl<kNSeedPlanes; isl++) fSeedTB[isl] = 0x0;
104 // Initialize debug stream
105 if(rec) SetReconstructor(rec);
108 //____________________________________________________________________
109 AliTRDtrackerV1::~AliTRDtrackerV1()
115 if(fgRieman) delete fgRieman; fgRieman = 0x0;
116 if(fgTiltedRieman) delete fgTiltedRieman; fgTiltedRieman = 0x0;
117 if(fgTiltedRiemanConstrained) delete fgTiltedRiemanConstrained; fgTiltedRiemanConstrained = 0x0;
118 for(Int_t isl =0; isl<kNSeedPlanes; isl++) if(fSeedTB[isl]) delete fSeedTB[isl];
119 if(fTracks) {fTracks->Delete(); delete fTracks;}
120 if(fTracklets) {fTracklets->Delete(); delete fTracklets;}
122 fClusters->Delete(); delete fClusters;
124 if(fGeom) delete fGeom;
127 //____________________________________________________________________
128 Int_t AliTRDtrackerV1::Clusters2Tracks(AliESDEvent *esd)
131 // Steering stand alone tracking for full TRD detector
134 // esd : The ESD event. On output it contains
135 // the ESD tracks found in TRD.
138 // Number of tracks found in the TRD detector.
140 // Detailed description
141 // 1. Launch individual SM trackers.
142 // See AliTRDtrackerV1::Clusters2TracksSM() for details.
145 if(!fReconstructor->GetRecoParam() ){
146 AliError("Reconstruction configuration not initialized. Call first AliTRDReconstructor::SetRecoParam().");
150 //AliInfo("Start Track Finder ...");
152 for(int ism=0; ism<AliTRDgeometry::kNsector; ism++){
153 // for(int ism=1; ism<2; ism++){
154 //AliInfo(Form("Processing supermodule %i ...", ism));
155 ntracks += Clusters2TracksSM(ism, esd);
157 AliInfo(Form("Number of found tracks : %d", ntracks));
162 //_____________________________________________________________________________
163 Bool_t AliTRDtrackerV1::GetTrackPoint(Int_t index, AliTrackPoint &p) const
165 //AliInfo(Form("Asking for tracklet %d", index));
167 // reset position of the point before using it
168 p.SetXYZ(0., 0., 0.);
169 AliTRDseedV1 *tracklet = GetTracklet(index);
170 if (!tracklet) return kFALSE;
172 // get detector for this tracklet
173 Int_t idet = tracklet->GetDetector();
176 local[0] = tracklet->GetX0();
177 local[1] = tracklet->GetYfit(0);
178 local[2] = tracklet->GetZfit(0);
180 fGeom->RotateBack(idet, local, global);
181 p.SetXYZ(global[0],global[1],global[2]);
185 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
186 switch (fGeom->GetLayer(idet)) {
188 iLayer = AliGeomManager::kTRD1;
191 iLayer = AliGeomManager::kTRD2;
194 iLayer = AliGeomManager::kTRD3;
197 iLayer = AliGeomManager::kTRD4;
200 iLayer = AliGeomManager::kTRD5;
203 iLayer = AliGeomManager::kTRD6;
206 Int_t modId = fGeom->GetSector(idet) * fGeom->Nstack() + fGeom->GetStack(idet);
207 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer, modId);
208 p.SetVolumeID(volid);
213 //____________________________________________________________________
214 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitter()
216 if(!fgTiltedRieman) fgTiltedRieman = new TLinearFitter(4, "hyp4");
217 return fgTiltedRieman;
220 //____________________________________________________________________
221 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitterConstraint()
223 if(!fgTiltedRiemanConstrained) fgTiltedRiemanConstrained = new TLinearFitter(2, "hyp2");
224 return fgTiltedRiemanConstrained;
227 //____________________________________________________________________
228 AliRieman* AliTRDtrackerV1::GetRiemanFitter()
230 if(!fgRieman) fgRieman = new AliRieman(AliTRDtrackingChamber::kNTimeBins * AliTRDgeometry::kNlayer);
234 //_____________________________________________________________________________
235 Int_t AliTRDtrackerV1::PropagateBack(AliESDEvent *event)
238 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
239 // backpropagated by the TPC tracker. Each seed is first propagated
240 // to the TRD, and then its prolongation is searched in the TRD.
241 // If sufficiently long continuation of the track is found in the TRD
242 // the track is updated, otherwise it's stored as originaly defined
243 // by the TPC tracker.
246 // Calibration monitor
247 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
248 if (!calibra) AliInfo("Could not get Calibra instance\n");
250 Int_t found = 0; // number of tracks found
251 Float_t foundMin = 20.0;
253 Float_t *quality = 0x0;
255 Int_t nSeed = event->GetNumberOfTracks();
257 quality = new Float_t[nSeed];
258 index = new Int_t[nSeed];
259 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
260 AliESDtrack *seed = event->GetTrack(iSeed);
261 Double_t covariance[15];
262 seed->GetExternalCovariance(covariance);
263 quality[iSeed] = covariance[0] + covariance[2];
265 // Sort tracks according to covariance of local Y and Z
266 TMath::Sort(nSeed,quality,index,kFALSE);
269 // Backpropagate all seeds
272 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
274 // Get the seeds in sorted sequence
275 AliESDtrack *seed = event->GetTrack(index[iSeed]);
277 // Check the seed status
278 ULong_t status = seed->GetStatus();
279 if ((status & AliESDtrack::kTPCout) == 0) continue;
280 if ((status & AliESDtrack::kTRDout) != 0) continue;
282 // Do the back prolongation
283 new(&track) AliTRDtrackV1(*seed);
284 track.SetReconstructor(fReconstructor);
286 //Int_t lbl = seed->GetLabel();
287 //track.SetSeedLabel(lbl);
289 // Make backup and mark entrance in the TRD
290 seed->UpdateTrackParams(&track, AliESDtrack::kTRDin);
291 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
292 Float_t p4 = track.GetC(track.GetBz());
293 expectedClr = FollowBackProlongation(track);
295 if (expectedClr<0) continue; // Back prolongation failed
299 // computes PID for track
301 // update calibration references using this track
302 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(&track);
303 // save calibration object
304 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0 /*&& quality TODO*/){
305 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(track);
306 calibTrack->SetOwner();
307 seed->AddCalibObject(calibTrack);
310 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
311 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
312 track.UpdateESDtrack(seed);
316 if ((TMath::Abs(track.GetC(track.GetBz()) - p4) / TMath::Abs(p4) < 0.2) ||(track.Pt() > 0.8)) {
318 // Make backup for back propagation
320 Int_t foundClr = track.GetNumberOfClusters();
321 if (foundClr >= foundMin) {
322 //AliInfo(Form("Making backup track ncls [%d]...", foundClr));
324 //track.CookdEdxTimBin(seed->GetID());
325 track.CookLabel(1. - fgkLabelFraction);
326 if(track.GetBackupTrack()) UseClusters(track.GetBackupTrack());
328 // Sign only gold tracks
329 if (track.GetChi2() / track.GetNumberOfClusters() < 4) {
330 if ((seed->GetKinkIndex(0) == 0) && (track.Pt() < 1.5)){
331 //UseClusters(&track);
334 Bool_t isGold = kFALSE;
337 if (track.GetChi2() / track.GetNumberOfClusters() < 5) {
338 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
344 if ((!isGold) && (track.GetNCross() == 0) && (track.GetChi2() / track.GetNumberOfClusters() < 7)) {
345 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
346 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
351 if ((!isGold) && (track.GetBackupTrack())) {
352 if ((track.GetBackupTrack()->GetNumberOfClusters() > foundMin) && ((track.GetBackupTrack()->GetChi2()/(track.GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
353 seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
358 //if ((track->StatusForTOF() > 0) && (track->GetNCross() == 0) && (Float_t(track->GetNumberOfClusters()) / Float_t(track->GetNExpected()) > 0.4)) {
359 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
364 // Propagation to the TOF (I.Belikov)
365 if (track.IsStopped() == kFALSE) {
366 Double_t xtof = 371.0;
367 Double_t xTOF0 = 370.0;
369 Double_t c2 = track.GetSnp() + track.GetC(track.GetBz()) * (xtof - track.GetX());
370 if (TMath::Abs(c2) >= 0.99) continue;
372 if (!PropagateToX(track, xTOF0, fgkMaxStep)) continue;
374 // Energy losses taken to the account - check one more time
375 c2 = track.GetSnp() + track.GetC(track.GetBz()) * (xtof - track.GetX());
376 if (TMath::Abs(c2) >= 0.99) continue;
378 //if (!PropagateToX(*track,xTOF0,fgkMaxStep)) {
379 // fHBackfit->Fill(7);
384 Double_t ymax = xtof * TMath::Tan(0.5 * AliTRDgeometry::GetAlpha());
386 track.GetYAt(xtof,GetBz(),y);
388 if (!track.Rotate( AliTRDgeometry::GetAlpha())) continue;
389 }else if (y < -ymax) {
390 if (!track.Rotate(-AliTRDgeometry::GetAlpha())) continue;
393 if (track.PropagateTo(xtof)) {
394 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
395 track.UpdateESDtrack(seed);
398 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
399 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
401 track.UpdateESDtrack(seed);
405 seed->SetTRDQuality(track.StatusForTOF());
406 seed->SetTRDBudget(track.GetBudget(0));
408 if(index) delete [] index;
409 if(quality) delete [] quality;
412 AliInfo(Form("Number of seeds: %d", nSeed));
413 AliInfo(Form("Number of back propagated TRD tracks: %d", found));
415 // run stand alone tracking
416 if (fReconstructor->IsSeeding()) Clusters2Tracks(event);
422 //____________________________________________________________________
423 Int_t AliTRDtrackerV1::RefitInward(AliESDEvent *event)
426 // Refits tracks within the TRD. The ESD event is expected to contain seeds
427 // at the outer part of the TRD.
428 // The tracks are propagated to the innermost time bin
429 // of the TRD and the ESD event is updated
430 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
433 Int_t nseed = 0; // contor for loaded seeds
434 Int_t found = 0; // contor for updated TRD tracks
438 for (Int_t itrack = 0; itrack < event->GetNumberOfTracks(); itrack++) {
439 AliESDtrack *seed = event->GetTrack(itrack);
440 new(&track) AliTRDtrackV1(*seed);
442 if (track.GetX() < 270.0) {
443 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
447 ULong_t status = seed->GetStatus();
448 // reject tracks which failed propagation in the TRD
449 if((status & AliESDtrack::kTRDout) == 0) continue;
451 // reject tracks which are produced by the TRD stand alone track finder.
452 if((status & AliESDtrack::kTRDin) == 0) continue;
455 track.ResetCovariance(50.0);
457 // do the propagation and processing
458 Bool_t kUPDATE = kFALSE;
459 Double_t xTPC = 250.0;
460 if(FollowProlongation(track)){
462 if (PropagateToX(track, xTPC, fgkMaxStep)) { // -with update
463 seed->UpdateTrackParams(&track, AliESDtrack::kTRDrefit);
469 // Prolongate to TPC without update
471 AliTRDtrackV1 tt(*seed);
472 if (PropagateToX(tt, xTPC, fgkMaxStep)) seed->UpdateTrackParams(&tt, AliESDtrack::kTRDrefit);
475 AliInfo(Form("Number of loaded seeds: %d",nseed));
476 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
481 //____________________________________________________________________
482 Int_t AliTRDtrackerV1::FollowProlongation(AliTRDtrackV1 &t)
484 // Extrapolates the TRD track in the TPC direction.
487 // t : the TRD track which has to be extrapolated
490 // number of clusters attached to the track
492 // Detailed description
494 // Starting from current radial position of track <t> this function
495 // extrapolates the track through the 6 TRD layers. The following steps
496 // are being performed for each plane:
498 // a. get plane limits in the local x direction
499 // b. check crossing sectors
500 // c. check track inclination
501 // 2. search tracklet in the tracker list (see GetTracklet() for details)
502 // 3. evaluate material budget using the geo manager
503 // 4. propagate and update track using the tracklet information.
508 Int_t nClustersExpected = 0;
509 Int_t lastplane = 5; //GetLastPlane(&t);
510 for (Int_t iplane = lastplane; iplane >= 0; iplane--) {
512 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
513 if(!tracklet) continue;
514 if(!tracklet->IsOK()) AliWarning("tracklet not OK");
516 Double_t x = tracklet->GetX0();
517 // reject tracklets which are not considered for inward refit
518 if(x > t.GetX()+fgkMaxStep) continue;
520 // append tracklet to track
521 t.SetTracklet(tracklet, index);
523 if (x < (t.GetX()-fgkMaxStep) && !PropagateToX(t, x+fgkMaxStep, fgkMaxStep)) break;
524 if (!AdjustSector(&t)) break;
526 // Start global position
530 // End global position
531 Double_t alpha = t.GetAlpha(), y, z;
532 if (!t.GetProlongation(x,y,z)) break;
534 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
535 xyz1[1] = x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
538 Double_t length = TMath::Sqrt(
539 (xyz0[0]-xyz1[0])*(xyz0[0]-xyz1[0]) +
540 (xyz0[1]-xyz1[1])*(xyz0[1]-xyz1[1]) +
541 (xyz0[2]-xyz1[2])*(xyz0[2]-xyz1[2])
544 // Get material budget
546 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) break;
547 Double_t xrho= param[0]*param[4];
548 Double_t xx0 = param[1]; // Get mean propagation parameters
550 // Propagate and update
551 t.PropagateTo(x, xx0, xrho);
552 if (!AdjustSector(&t)) break;
555 Double_t maxChi2 = t.GetPredictedChi2(tracklet);
556 if (maxChi2 < 1e+10 && t.Update(tracklet, maxChi2)){
557 nClustersExpected += tracklet->GetN();
561 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
563 for(int iplane=0; iplane<AliTRDgeometry::kNlayer; iplane++){
564 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
565 if(!tracklet) continue;
566 t.SetTracklet(tracklet, index);
569 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
570 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
571 cstreamer << "FollowProlongation"
572 << "EventNumber=" << eventNumber
573 << "ncl=" << nClustersExpected
578 return nClustersExpected;
582 //_____________________________________________________________________________
583 Int_t AliTRDtrackerV1::FollowBackProlongation(AliTRDtrackV1 &t)
585 // Extrapolates the TRD track in the TOF direction.
588 // t : the TRD track which has to be extrapolated
591 // number of clusters attached to the track
593 // Detailed description
595 // Starting from current radial position of track <t> this function
596 // extrapolates the track through the 6 TRD layers. The following steps
597 // are being performed for each plane:
599 // a. get plane limits in the local x direction
600 // b. check crossing sectors
601 // c. check track inclination
602 // 2. build tracklet (see AliTRDseed::AttachClusters() for details)
603 // 3. evaluate material budget using the geo manager
604 // 4. propagate and update track using the tracklet information.
609 Int_t nClustersExpected = 0;
610 Double_t clength = AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick();
611 AliTRDtrackingChamber *chamber = 0x0;
613 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
614 // in case of stand alone tracking we store all the pointers to the tracklets in a temporary array
615 AliTRDseedV1 *tracklets[kNPlanes];
616 memset(tracklets, 0, sizeof(AliTRDseedV1 *) * kNPlanes);
617 for(Int_t ip = 0; ip < kNPlanes; ip++){
618 tracklets[ip] = t.GetTracklet(ip);
622 // Loop through the TRD layers
623 for (Int_t ilayer = 0; ilayer < AliTRDgeometry::Nlayer(); ilayer++) {
624 // BUILD TRACKLET IF NOT ALREADY BUILT
625 Double_t x = 0., y, z, alpha;
626 ptrTracklet = tracklets[ilayer];
628 ptrTracklet = new(&tracklet) AliTRDseedV1(ilayer);
629 ptrTracklet->SetReconstructor(fReconstructor);
630 alpha = t.GetAlpha();
631 Int_t sector = Int_t(alpha/AliTRDgeometry::GetAlpha() + (alpha>0. ? 0 : AliTRDgeometry::kNsector));
633 if(!fTrSec[sector].GetNChambers()) continue;
635 if((x = fTrSec[sector].GetX(ilayer)) < 1.) continue;
637 if (!t.GetProlongation(x, y, z)) return -1/*nClustersExpected*/;
638 Int_t stack = fGeom->GetStack(z, ilayer);
639 Int_t nCandidates = stack >= 0 ? 1 : 2;
640 z -= stack >= 0 ? 0. : 4.;
642 for(int icham=0; icham<nCandidates; icham++, z+=8){
643 if((stack = fGeom->GetStack(z, ilayer)) < 0) continue;
645 if(!(chamber = fTrSec[sector].GetChamber(stack, ilayer))) continue;
647 if(chamber->GetNClusters() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
651 AliTRDpadPlane *pp = fGeom->GetPadPlane(ilayer, stack);
652 tracklet.SetTilt(TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle()));
653 tracklet.SetPadLength(pp->GetLengthIPad());
654 tracklet.SetDetector(chamber->GetDetector());
656 if(!tracklet.Init(&t)){
658 return nClustersExpected;
660 if(!tracklet.AttachClustersIter(chamber, 1000./*, kTRUE*/)) continue;
663 if(tracklet.GetN() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
667 //ptrTracklet->UseClusters();
668 } else ptrTracklet->Init(&t);
669 if(!ptrTracklet->IsOK()){
670 if(x < 1.) continue; //temporary
671 if(!PropagateToX(t, x-fgkMaxStep, fgkMaxStep)) return -1/*nClustersExpected*/;
672 if(!AdjustSector(&t)) return -1/*nClustersExpected*/;
673 if(TMath::Abs(t.GetSnp()) > fgkMaxSnp) return -1/*nClustersExpected*/;
677 // Propagate closer to the current chamber if neccessary
679 if (x > (fgkMaxStep + t.GetX()) && !PropagateToX(t, x-fgkMaxStep, fgkMaxStep)) return -1/*nClustersExpected*/;
680 if (!AdjustSector(&t)) return -1/*nClustersExpected*/;
681 if (TMath::Abs(t.GetSnp()) > fgkMaxSnp) return -1/*nClustersExpected*/;
683 // load tracklet to the tracker and the track
684 ptrTracklet = SetTracklet(ptrTracklet);
685 t.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
688 // Calculate the mean material budget along the path inside the chamber
689 //Calculate global entry and exit positions of the track in chamber (only track prolongation)
690 Double_t xyz0[3]; // entry point
692 alpha = t.GetAlpha();
693 x = ptrTracklet->GetX0();
694 if (!t.GetProlongation(x, y, z)) return -1/*nClustersExpected*/;
695 Double_t xyz1[3]; // exit point
696 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
697 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
700 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) return -1;
701 // The mean propagation parameters
702 Double_t xrho = param[0]*param[4]; // density*length
703 Double_t xx0 = param[1]; // radiation length
705 // Propagate and update track
706 if (!t.PropagateTo(x, xx0, xrho)) return -1/*nClustersExpected*/;
707 if (!AdjustSector(&t)) return -1/*nClustersExpected*/;
708 Double_t maxChi2 = t.GetPredictedChi2(ptrTracklet);
709 if (!t.Update(ptrTracklet, maxChi2)) return -1/*nClustersExpected*/;
711 nClustersExpected += ptrTracklet->GetN();
712 //t.SetTracklet(&tracklet, index);
714 // Reset material budget if 2 consecutive gold
715 if(ilayer>0 && t.GetTracklet(ilayer-1) && ptrTracklet->GetN() + t.GetTracklet(ilayer-1)->GetN() > 20) t.SetBudget(2, 0.);
717 // Make backup of the track until is gold
718 // TO DO update quality check of the track.
719 // consider comparison with fTimeBinsRange
720 Float_t ratio0 = ptrTracklet->GetN() / Float_t(fgNTimeBins);
721 //Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);
722 //printf("tracklet.GetChi2() %f [< 18.0]\n", tracklet.GetChi2());
723 //printf("ratio0 %f [> 0.8]\n", ratio0);
724 //printf("ratio1 %f [> 0.6]\n", ratio1);
725 //printf("ratio0+ratio1 %f [> 1.5]\n", ratio0+ratio1);
726 //printf("t.GetNCross() %d [== 0]\n", t.GetNCross());
727 //printf("TMath::Abs(t.GetSnp()) %f [< 0.85]\n", TMath::Abs(t.GetSnp()));
728 //printf("t.GetNumberOfClusters() %d [> 20]\n", t.GetNumberOfClusters());
730 if (//(tracklet.GetChi2() < 18.0) && TO DO check with FindClusters and move it to AliTRDseed::Update
733 //(ratio0+ratio1 > 1.5) &&
734 (t.GetNCross() == 0) &&
735 (TMath::Abs(t.GetSnp()) < 0.85) &&
736 (t.GetNumberOfClusters() > 20)) t.MakeBackupTrack();
740 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
741 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
742 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
743 //AliTRDtrackV1 *debugTrack = new AliTRDtrackV1(t);
744 //debugTrack->SetOwner();
745 cstreamer << "FollowBackProlongation"
746 << "EventNumber=" << eventNumber
747 << "ncl=" << nClustersExpected
748 //<< "track.=" << debugTrack
752 return nClustersExpected;
755 //_________________________________________________________________________
756 Float_t AliTRDtrackerV1::FitRieman(AliTRDseedV1 *tracklets, Double_t *chi2, Int_t *planes){
758 // Fits a Riemann-circle to the given points without tilting pad correction.
759 // The fit is performed using an instance of the class AliRieman (equations
760 // and transformations see documentation of this class)
761 // Afterwards all the tracklets are Updated
763 // Parameters: - Array of tracklets (AliTRDseedV1)
764 // - Storage for the chi2 values (beginning with direction z)
765 // - Seeding configuration
766 // Output: - The curvature
768 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
770 Int_t allplanes[] = {0, 1, 2, 3, 4, 5};
771 Int_t *ppl = &allplanes[0];
777 for(Int_t il = 0; il < maxLayers; il++){
778 if(!tracklets[ppl[il]].IsOK()) continue;
779 fitter->AddPoint(tracklets[ppl[il]].GetX0(), tracklets[ppl[il]].GetYfitR(0), tracklets[ppl[il]].GetZProb(),1,10);
782 // Set the reference position of the fit and calculate the chi2 values
783 memset(chi2, 0, sizeof(Double_t) * 2);
784 for(Int_t il = 0; il < maxLayers; il++){
785 // Reference positions
786 tracklets[ppl[il]].Init(fitter);
789 if((!tracklets[ppl[il]].IsOK()) && (!planes)) continue;
790 chi2[0] += tracklets[ppl[il]].GetChi2Y();
791 chi2[1] += tracklets[ppl[il]].GetChi2Z();
793 return fitter->GetC();
796 //_________________________________________________________________________
797 void AliTRDtrackerV1::FitRieman(AliTRDcluster **seedcl, Double_t chi2[2])
800 // Performs a Riemann helix fit using the seedclusters as spacepoints
801 // Afterwards the chi2 values are calculated and the seeds are updated
803 // Parameters: - The four seedclusters
804 // - The tracklet array (AliTRDseedV1)
805 // - The seeding configuration
810 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
812 for(Int_t i = 0; i < 4; i++)
813 fitter->AddPoint(seedcl[i]->GetX(), seedcl[i]->GetY(), seedcl[i]->GetZ(), 1, 10);
817 // Update the seed and calculated the chi2 value
818 chi2[0] = 0; chi2[1] = 0;
819 for(Int_t ipl = 0; ipl < kNSeedPlanes; ipl++){
821 chi2[0] += (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX()));
822 chi2[1] += (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX()));
827 //_________________________________________________________________________
828 Float_t AliTRDtrackerV1::FitTiltedRiemanConstraint(AliTRDseedV1 *tracklets, Double_t zVertex)
831 // Fits a helix to the clusters. Pad tilting is considered. As constraint it is
832 // assumed that the vertex position is set to 0.
833 // This method is very usefull for high-pt particles
834 // Basis for the fit: (x - x0)^2 + (y - y0)^2 - R^2 = 0
835 // x0, y0: Center of the circle
836 // Measured y-position: ymeas = y - tan(phiT)(zc - zt)
837 // zc: center of the pad row
838 // Equation which has to be fitted (after transformation):
839 // a + b * u + e * v + 2*(ymeas + tan(phiT)(z - zVertex))*t = 0
843 // v = 2 * x * tan(phiT) * t
844 // Parameters in the equation:
845 // a = -1/y0, b = x0/y0, e = dz/dx
847 // The Curvature is calculated by the following equation:
848 // - curv = a/Sqrt(b^2 + 1) = 1/R
849 // Parameters: - the 6 tracklets
850 // - the Vertex constraint
851 // Output: - the Chi2 value of the track
856 TLinearFitter *fitter = GetTiltedRiemanFitterConstraint();
857 fitter->StoreData(kTRUE);
858 fitter->ClearPoints();
859 AliTRDcluster *cl = 0x0;
861 Float_t x, y, z, w, t, error, tilt;
864 for(Int_t ilr = 0; ilr < AliTRDgeometry::kNlayer; ilr++){
865 if(!tracklets[ilr].IsOK()) continue;
866 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
867 if(!tracklets[ilr].IsUsable(itb)) continue;
868 cl = tracklets[ilr].GetClusters(itb);
872 tilt = tracklets[ilr].GetTilt();
874 t = 1./(x * x + y * y);
876 uvt[1] = 2. * x * t * tilt ;
877 w = 2. * (y + tilt * (z - zVertex)) * t;
878 error = 2. * 0.2 * t;
879 fitter->AddPoint(uvt, w, error);
885 // Calculate curvature
886 Double_t a = fitter->GetParameter(0);
887 Double_t b = fitter->GetParameter(1);
888 Double_t curvature = a/TMath::Sqrt(b*b + 1);
890 Float_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
891 for(Int_t ip = 0; ip < AliTRDtrackerV1::kNPlanes; ip++)
892 tracklets[ip].SetCC(curvature);
894 /* if(fReconstructor->GetStreamLevel() >= 5){
895 //Linear Model on z-direction
896 Double_t xref = CalculateReferenceX(tracklets); // Relative to the middle of the stack
897 Double_t slope = fitter->GetParameter(2);
898 Double_t zref = slope * xref;
899 Float_t chi2Z = CalculateChi2Z(tracklets, zref, slope, xref);
900 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
901 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
902 TTreeSRedirector &treeStreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
903 treeStreamer << "FitTiltedRiemanConstraint"
904 << "EventNumber=" << eventNumber
905 << "CandidateNumber=" << candidateNumber
906 << "Curvature=" << curvature
907 << "Chi2Track=" << chi2track
915 //_________________________________________________________________________
916 Float_t AliTRDtrackerV1::FitTiltedRieman(AliTRDseedV1 *tracklets, Bool_t sigError)
919 // Performs a Riemann fit taking tilting pad correction into account
920 // The equation of a Riemann circle, where the y position is substituted by the
921 // measured y-position taking pad tilting into account, has to be transformed
922 // into a 4-dimensional hyperplane equation
923 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
924 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
925 // zc: center of the pad row
926 // zt: z-position of the track
927 // The z-position of the track is assumed to be linear dependent on the x-position
928 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
929 // Transformation: u = 2 * x * t
930 // v = 2 * tan(phiT) * t
931 // w = 2 * tan(phiT) * (x - xref) * t
932 // t = 1 / (x^2 + ymeas^2)
933 // Parameters: a = -1/y0
935 // c = (R^2 -x0^2 - y0^2)/y0
938 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
939 // results from the simple riemann fit. Afterwards the fit is redone.
940 // The curvature is calculated according to the formula:
941 // curv = a/(1 + b^2 + c*a) = 1/R
943 // Paramters: - Array of tracklets (connected to the track candidate)
944 // - Flag selecting the error definition
945 // Output: - Chi2 values of the track (in Parameter list)
947 TLinearFitter *fitter = GetTiltedRiemanFitter();
948 fitter->StoreData(kTRUE);
949 fitter->ClearPoints();
950 AliTRDLeastSquare zfitter;
951 AliTRDcluster *cl = 0x0;
953 Double_t xref = CalculateReferenceX(tracklets);
954 Double_t x, y, z, t, tilt, dx, w, we;
957 // Containers for Least-square fitter
958 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
959 if(!tracklets[ipl].IsOK()) continue;
960 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
961 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
962 if (!tracklets[ipl].IsUsable(itb)) continue;
966 tilt = tracklets[ipl].GetTilt();
972 uvt[2] = 2. * tilt * t;
973 uvt[3] = 2. * tilt * dx * t;
974 w = 2. * (y + tilt*z) * t;
975 // error definition changes for the different calls
977 we *= sigError ? tracklets[ipl].GetSigmaY() : 0.2;
978 fitter->AddPoint(uvt, w, we);
979 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
986 Double_t offset = fitter->GetParameter(3);
987 Double_t slope = fitter->GetParameter(4);
989 // Linear fitter - not possible to make boundaries
990 // Do not accept non possible z and dzdx combinations
991 Bool_t acceptablez = kTRUE;
993 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
994 if(!tracklets[iLayer].IsOK()) continue;
995 zref = offset + slope * (tracklets[iLayer].GetX0() - xref);
996 if (TMath::Abs(tracklets[iLayer].GetZProb() - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
997 acceptablez = kFALSE;
1000 Double_t dzmf = zfitter.GetFunctionParameter(1);
1001 Double_t zmf = zfitter.GetFunctionValue(&xref);
1002 fgTiltedRieman->FixParameter(3, zmf);
1003 fgTiltedRieman->FixParameter(4, dzmf);
1005 fitter->ReleaseParameter(3);
1006 fitter->ReleaseParameter(4);
1007 offset = fitter->GetParameter(3);
1008 slope = fitter->GetParameter(4);
1011 // Calculate Curvarture
1012 Double_t a = fitter->GetParameter(0);
1013 Double_t b = fitter->GetParameter(1);
1014 Double_t c = fitter->GetParameter(2);
1015 Double_t curvature = 1.0 + b*b - c*a;
1016 if (curvature > 0.0)
1017 curvature = a / TMath::Sqrt(curvature);
1019 Double_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
1021 // Update the tracklets
1023 for(Int_t iLayer = 0; iLayer < AliTRDtrackerV1::kNPlanes; iLayer++) {
1025 x = tracklets[iLayer].GetX0();
1031 // y: R^2 = (x - x0)^2 + (y - y0)^2
1032 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1033 // R = Sqrt() = 1/Curvature
1034 // => y = y0 +/- Sqrt(1/Curvature^2 - (x - x0)^2)
1035 Double_t res = (x * a + b); // = (x - x0)/y0
1037 res = 1.0 - c * a + b * b - res; // = (R^2 - (x - x0)^2)/y0^2
1039 res = TMath::Sqrt(res);
1040 y = (1.0 - res) / a;
1043 // dy: R^2 = (x - x0)^2 + (y - y0)^2
1044 // => y = +/- Sqrt(R^2 - (x - x0)^2) + y0
1045 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1046 // Curvature: cr = 1/R = a/Sqrt(1 + b^2 - c*a)
1047 // => dy/dx = (x - x0)/(1/(cr^2) - (x - x0)^2)
1048 Double_t x0 = -b / a;
1049 if (-c * a + b * b + 1 > 0) {
1050 if (1.0/(curvature * curvature) - (x - x0) * (x - x0) > 0.0) {
1051 Double_t yderiv = (x - x0) / TMath::Sqrt(1.0/(curvature * curvature) - (x - x0) * (x - x0));
1052 if (a < 0) yderiv *= -1.0;
1056 z = offset + slope * (x - xref);
1058 tracklets[iLayer].SetYref(0, y);
1059 tracklets[iLayer].SetYref(1, dy);
1060 tracklets[iLayer].SetZref(0, z);
1061 tracklets[iLayer].SetZref(1, dz);
1062 tracklets[iLayer].SetC(curvature);
1063 tracklets[iLayer].SetChi2(chi2track);
1066 /* if(fReconstructor->GetStreamLevel() >=5){
1067 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
1068 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1069 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1070 Double_t chi2z = CalculateChi2Z(tracklets, offset, slope, xref);
1071 cstreamer << "FitTiltedRieman0"
1072 << "EventNumber=" << eventNumber
1073 << "CandidateNumber=" << candidateNumber
1075 << "Chi2Z=" << chi2z
1082 //____________________________________________________________________
1083 Double_t AliTRDtrackerV1::FitLine(const AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t err, Int_t np, AliTrackPoint *points)
1085 AliTRDLeastSquare yfitter, zfitter;
1086 AliTRDcluster *cl = 0x0;
1088 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1090 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1091 if(!(tracklet = track->GetTracklet(ipl))) continue;
1092 if(!tracklet->IsOK()) continue;
1093 new(&work[ipl]) AliTRDseedV1(*tracklet);
1095 tracklets = &work[0];
1098 Double_t xref = CalculateReferenceX(tracklets);
1099 Double_t x, y, z, dx, ye, yr, tilt;
1100 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1101 if(!tracklets[ipl].IsOK()) continue;
1102 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1103 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1104 if (!tracklets[ipl].IsUsable(itb)) continue;
1108 zfitter.AddPoint(&dx, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1112 Double_t z0 = zfitter.GetFunctionParameter(0);
1113 Double_t dzdx = zfitter.GetFunctionParameter(1);
1114 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1115 if(!tracklets[ipl].IsOK()) continue;
1116 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1117 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1118 if (!tracklets[ipl].IsUsable(itb)) continue;
1122 tilt = tracklets[ipl].GetTilt();
1124 yr = y + tilt*(z - z0 - dzdx*dx);
1125 // error definition changes for the different calls
1126 ye = tilt*TMath::Sqrt(cl->GetSigmaZ2());
1127 ye += err ? tracklets[ipl].GetSigmaY() : 0.2;
1128 yfitter.AddPoint(&dx, yr, ye);
1132 Double_t y0 = yfitter.GetFunctionParameter(0);
1133 Double_t dydx = yfitter.GetFunctionParameter(1);
1134 Double_t chi2 = 0.;//yfitter.GetChisquare()/Double_t(nPoints);
1136 //update track points array
1139 for(int ip=0; ip<np; ip++){
1140 points[ip].GetXYZ(xyz);
1141 xyz[1] = y0 + dydx * (xyz[0] - xref);
1142 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1143 points[ip].SetXYZ(xyz);
1150 //_________________________________________________________________________
1151 Double_t AliTRDtrackerV1::FitRiemanTilt(const AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t sigError, Int_t np, AliTrackPoint *points)
1154 // Performs a Riemann fit taking tilting pad correction into account
1155 // The equation of a Riemann circle, where the y position is substituted by the
1156 // measured y-position taking pad tilting into account, has to be transformed
1157 // into a 4-dimensional hyperplane equation
1158 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
1159 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
1160 // zc: center of the pad row
1161 // zt: z-position of the track
1162 // The z-position of the track is assumed to be linear dependent on the x-position
1163 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
1164 // Transformation: u = 2 * x * t
1165 // v = 2 * tan(phiT) * t
1166 // w = 2 * tan(phiT) * (x - xref) * t
1167 // t = 1 / (x^2 + ymeas^2)
1168 // Parameters: a = -1/y0
1170 // c = (R^2 -x0^2 - y0^2)/y0
1173 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
1174 // results from the simple riemann fit. Afterwards the fit is redone.
1175 // The curvature is calculated according to the formula:
1176 // curv = a/(1 + b^2 + c*a) = 1/R
1178 // Paramters: - Array of tracklets (connected to the track candidate)
1179 // - Flag selecting the error definition
1180 // Output: - Chi2 values of the track (in Parameter list)
1182 TLinearFitter *fitter = GetTiltedRiemanFitter();
1183 fitter->StoreData(kTRUE);
1184 fitter->ClearPoints();
1185 AliTRDLeastSquare zfitter;
1186 AliTRDcluster *cl = 0x0;
1188 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1190 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1191 if(!(tracklet = track->GetTracklet(ipl))) continue;
1192 if(!tracklet->IsOK()) continue;
1193 new(&work[ipl]) AliTRDseedV1(*tracklet);
1195 tracklets = &work[0];
1198 Double_t xref = CalculateReferenceX(tracklets);
1199 Double_t x, y, z, t, tilt, dx, w, we;
1202 // Containers for Least-square fitter
1203 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1204 if(!tracklets[ipl].IsOK()) continue;
1205 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1206 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1207 if (!tracklets[ipl].IsUsable(itb)) continue;
1211 tilt = tracklets[ipl].GetTilt();
1215 uvt[0] = 2. * x * t;
1217 uvt[2] = 2. * tilt * t;
1218 uvt[3] = 2. * tilt * dx * t;
1219 w = 2. * (y + tilt*z) * t;
1220 // error definition changes for the different calls
1222 we *= sigError ? tracklets[ipl].GetSigmaY() : 0.2;
1223 fitter->AddPoint(uvt, w, we);
1224 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1228 if(fitter->Eval()) return 1.E10;
1230 Double_t z0 = fitter->GetParameter(3);
1231 Double_t dzdx = fitter->GetParameter(4);
1234 // Linear fitter - not possible to make boundaries
1235 // Do not accept non possible z and dzdx combinations
1236 Bool_t accept = kTRUE;
1237 Double_t zref = 0.0;
1238 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
1239 if(!tracklets[iLayer].IsOK()) continue;
1240 zref = z0 + dzdx * (tracklets[iLayer].GetX0() - xref);
1241 if (TMath::Abs(tracklets[iLayer].GetZProb() - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
1246 Double_t dzmf = zfitter.GetFunctionParameter(1);
1247 Double_t zmf = zfitter.GetFunctionValue(&xref);
1248 fitter->FixParameter(3, zmf);
1249 fitter->FixParameter(4, dzmf);
1251 fitter->ReleaseParameter(3);
1252 fitter->ReleaseParameter(4);
1253 z0 = fitter->GetParameter(3); // = zmf ?
1254 dzdx = fitter->GetParameter(4); // = dzmf ?
1257 // Calculate Curvature
1258 Double_t a = fitter->GetParameter(0);
1259 Double_t b = fitter->GetParameter(1);
1260 Double_t c = fitter->GetParameter(2);
1261 Double_t y0 = 1. / a;
1262 Double_t x0 = -b * y0;
1263 Double_t tmp = y0*y0 + x0*x0 - c*y0;
1264 if(tmp<=0.) return 1.E10;
1265 Double_t R = TMath::Sqrt(tmp);
1266 Double_t C = 1.0 + b*b - c*a;
1267 if (C > 0.0) C = a / TMath::Sqrt(C);
1269 // Calculate chi2 of the fit
1270 Double_t chi2 = fitter->GetChisquare()/Double_t(nPoints);
1272 // Update the tracklets
1274 for(Int_t ip = 0; ip < kNPlanes; ip++) {
1275 x = tracklets[ip].GetX0();
1276 tmp = R*R-(x-x0)*(x-x0);
1277 if(tmp <= 0.) continue;
1278 tmp = TMath::Sqrt(tmp);
1280 // y: R^2 = (x - x0)^2 + (y - y0)^2
1281 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1282 tracklets[ip].SetYref(0, y0 - (y0>0.?1.:-1)*tmp);
1283 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1284 tracklets[ip].SetYref(1, (x - x0) / tmp);
1285 tracklets[ip].SetZref(0, z0 + dzdx * (x - xref));
1286 tracklets[ip].SetZref(1, dzdx);
1287 tracklets[ip].SetC(C);
1288 tracklets[ip].SetChi2(chi2);
1291 //update track points array
1294 for(int ip=0; ip<np; ip++){
1295 points[ip].GetXYZ(xyz);
1296 xyz[1] = TMath::Abs(xyz[0] - x0) > R ? 100. : y0 - (y0>0.?1.:-1.)*TMath::Sqrt(R*R-(xyz[0]-x0)*(xyz[0]-x0));
1297 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1298 points[ip].SetXYZ(xyz);
1306 //____________________________________________________________________
1307 Double_t AliTRDtrackerV1::FitKalman(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t up, Int_t np, AliTrackPoint *points)
1309 // Kalman filter implementation for the TRD.
1310 // It returns the positions of the fit in the array "points"
1312 // Author : A.Bercuci@gsi.de
1314 // printf("Start track @ x[%f]\n", track->GetX());
1316 //prepare marker points along the track
1317 Int_t ip = np ? 0 : 1;
1319 if((up?-1:1) * (track->GetX() - points[ip].GetX()) > 0.) break;
1320 //printf("AliTRDtrackerV1::FitKalman() : Skip track marker x[%d] = %7.3f. Before track start ( %7.3f ).\n", ip, points[ip].GetX(), track->GetX());
1323 //if(points) printf("First marker point @ x[%d] = %f\n", ip, points[ip].GetX());
1326 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
1328 //Loop through the TRD planes
1329 for (Int_t jplane = 0; jplane < kNPlanes; jplane++) {
1330 // GET TRACKLET OR BUILT IT
1331 Int_t iplane = up ? jplane : kNPlanes - 1 - jplane;
1333 if(!(ptrTracklet = &tracklets[iplane])) continue;
1335 if(!(ptrTracklet = track->GetTracklet(iplane))){
1336 /*AliTRDtrackerV1 *tracker = 0x0;
1337 if(!(tracker = dynamic_cast<AliTRDtrackerV1*>( AliTRDReconstructor::Tracker()))) continue;
1338 ptrTracklet = new(&tracklet) AliTRDseedV1(iplane);
1339 if(!tracker->MakeTracklet(ptrTracklet, track)) */
1343 if(!ptrTracklet->IsOK()) continue;
1345 Double_t x = ptrTracklet->GetX0();
1348 //don't do anything if next marker is after next update point.
1349 if((up?-1:1) * (points[ip].GetX() - x) - fgkMaxStep < 0) break;
1350 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1352 Double_t xyz[3]; // should also get the covariance
1354 track->Global2LocalPosition(xyz, track->GetAlpha());
1355 points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1358 // printf("plane[%d] tracklet[%p] x[%f]\n", iplane, ptrTracklet, x);
1360 // Propagate closer to the next update point
1361 if(((up?-1:1) * (x - track->GetX()) + fgkMaxStep < 0) && !PropagateToX(*track, x + (up?-1:1)*fgkMaxStep, fgkMaxStep)) return -1.;
1363 if(!AdjustSector(track)) return -1;
1364 if(TMath::Abs(track->GetSnp()) > fgkMaxSnp) return -1;
1366 //load tracklet to the tracker and the track
1368 if((index = FindTracklet(ptrTracklet)) < 0){
1369 ptrTracklet = SetTracklet(&tracklet);
1370 index = fTracklets->GetEntriesFast()-1;
1372 track->SetTracklet(ptrTracklet, index);*/
1375 // register tracklet to track with tracklet creation !!
1376 // PropagateBack : loaded tracklet to the tracker and update index
1377 // RefitInward : update index
1378 // MakeTrack : loaded tracklet to the tracker and update index
1379 if(!tracklets) track->SetTracklet(ptrTracklet, -1);
1382 //Calculate the mean material budget along the path inside the chamber
1383 Double_t xyz0[3]; track->GetXYZ(xyz0);
1384 Double_t alpha = track->GetAlpha();
1385 Double_t xyz1[3], y, z;
1386 if(!track->GetProlongation(x, y, z)) return -1;
1387 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
1388 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
1390 if((xyz0[0] - xyz1[9] < 1e-3) && (xyz0[0] - xyz1[9] < 1e-3)) continue; // check wheter we are at the same global x position
1392 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param) <=0.) break;
1393 Double_t xrho = param[0]*param[4]; // density*length
1394 Double_t xx0 = param[1]; // radiation length
1396 //Propagate the track
1397 track->PropagateTo(x, xx0, xrho);
1398 if (!AdjustSector(track)) break;
1401 Double_t chi2 = track->GetPredictedChi2(ptrTracklet);
1402 if(chi2<1e+10) track->Update(ptrTracklet, chi2);
1405 //Reset material budget if 2 consecutive gold
1406 if(iplane>0 && track->GetTracklet(iplane-1) && ptrTracklet->GetN() + track->GetTracklet(iplane-1)->GetN() > 20) track->SetBudget(2, 0.);
1407 } // end planes loop
1411 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1413 Double_t xyz[3]; // should also get the covariance
1415 track->Global2LocalPosition(xyz, track->GetAlpha());
1416 points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1420 return track->GetChi2();
1423 //_________________________________________________________________________
1424 Float_t AliTRDtrackerV1::CalculateChi2Z(AliTRDseedV1 *tracklets, Double_t offset, Double_t slope, Double_t xref)
1427 // Calculates the chi2-value of the track in z-Direction including tilting pad correction.
1428 // A linear dependence on the x-value serves as a model.
1429 // The parameters are related to the tilted Riemann fit.
1430 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
1431 // - the offset for the reference x
1433 // - the reference x position
1434 // Output: - The Chi2 value of the track in z-Direction
1436 Float_t chi2Z = 0, nLayers = 0;
1437 for (Int_t iLayer = 0; iLayer < AliTRDgeometry::kNlayer; iLayer++) {
1438 if(!tracklets[iLayer].IsOK()) continue;
1439 Double_t z = offset + slope * (tracklets[iLayer].GetX0() - xref);
1440 chi2Z += TMath::Abs(tracklets[iLayer].GetMeanz() - z);
1443 chi2Z /= TMath::Max((nLayers - 3.0),1.0);
1447 //_____________________________________________________________________________
1448 Int_t AliTRDtrackerV1::PropagateToX(AliTRDtrackV1 &t, Double_t xToGo, Double_t maxStep)
1451 // Starting from current X-position of track <t> this function
1452 // extrapolates the track up to radial position <xToGo>.
1453 // Returns 1 if track reaches the plane, and 0 otherwise
1456 const Double_t kEpsilon = 0.00001;
1458 // Current track X-position
1459 Double_t xpos = t.GetX();
1461 // Direction: inward or outward
1462 Double_t dir = (xpos < xToGo) ? 1.0 : -1.0;
1464 while (((xToGo - xpos) * dir) > kEpsilon) {
1473 // The next step size
1474 Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
1476 // Get the global position of the starting point
1479 // X-position after next step
1482 // Get local Y and Z at the X-position of the next step
1483 if (!t.GetProlongation(x,y,z)) {
1484 return 0; // No prolongation possible
1487 // The global position of the end point of this prolongation step
1488 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1489 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1492 // Calculate the mean material budget between start and
1493 // end point of this prolongation step
1494 if(AliTracker::MeanMaterialBudget(xyz0, xyz1, param)<=0.) return 0;
1496 // Propagate the track to the X-position after the next step
1497 if (!t.PropagateTo(x,param[1],param[0]*param[4])) {
1501 // Rotate the track if necessary
1504 // New track X-position
1514 //_____________________________________________________________________________
1515 Int_t AliTRDtrackerV1::ReadClusters(TClonesArray* &array, TTree *clusterTree) const
1518 // Reads AliTRDclusters from the file.
1519 // The names of the cluster tree and branches
1520 // should match the ones used in AliTRDclusterizer::WriteClusters()
1523 Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
1524 TObjArray *clusterArray = new TObjArray(nsize+1000);
1526 TBranch *branch = clusterTree->GetBranch("TRDcluster");
1528 AliError("Can't get the branch !");
1531 branch->SetAddress(&clusterArray);
1534 Float_t nclusters = fReconstructor->GetRecoParam()->GetNClusters();
1535 if(fReconstructor->IsHLT()) nclusters /= AliTRDgeometry::kNsector;
1536 array = new TClonesArray("AliTRDcluster", Int_t(nclusters));
1537 array->SetOwner(kTRUE);
1540 // Loop through all entries in the tree
1541 Int_t nEntries = (Int_t) clusterTree->GetEntries();
1544 AliTRDcluster *c = 0x0;
1545 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1547 nbytes += clusterTree->GetEvent(iEntry);
1549 // Get the number of points in the detector
1550 Int_t nCluster = clusterArray->GetEntriesFast();
1551 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1552 if(!(c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster))) continue;
1554 new((*fClusters)[ncl++]) AliTRDcluster(*c);
1555 delete (clusterArray->RemoveAt(iCluster));
1559 delete clusterArray;
1564 //_____________________________________________________________________________
1565 Int_t AliTRDtrackerV1::LoadClusters(TTree *cTree)
1568 // Fills clusters into TRD tracking sectors
1571 if(!fReconstructor->IsWritingClusters()){
1572 fClusters = AliTRDReconstructor::GetClusters();
1574 if (ReadClusters(fClusters, cTree)) {
1575 AliError("Problem with reading the clusters !");
1581 if(!fClusters || !fClusters->GetEntriesFast()){
1582 AliInfo("No TRD clusters");
1587 BuildTrackingContainers();
1589 //Int_t ncl = fClusters->GetEntriesFast();
1590 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1595 //_____________________________________________________________________________
1596 Int_t AliTRDtrackerV1::LoadClusters(TClonesArray *clusters)
1599 // Fills clusters into TRD tracking sectors
1600 // Function for use in the HLT
1602 if(!clusters || !clusters->GetEntriesFast()){
1603 AliInfo("No TRD clusters");
1607 fClusters = clusters;
1611 BuildTrackingContainers();
1613 //Int_t ncl = fClusters->GetEntriesFast();
1614 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1620 //____________________________________________________________________
1621 Int_t AliTRDtrackerV1::BuildTrackingContainers()
1623 // Building tracking containers for clusters
1625 Int_t nin =0, icl = fClusters->GetEntriesFast();
1627 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(icl);
1628 if(c->IsInChamber()) nin++;
1629 Int_t detector = c->GetDetector();
1630 Int_t sector = fGeom->GetSector(detector);
1631 Int_t stack = fGeom->GetStack(detector);
1632 Int_t layer = fGeom->GetLayer(detector);
1634 fTrSec[sector].GetChamber(stack, layer, kTRUE)->InsertCluster(c, icl);
1637 const AliTRDCalDet *cal = AliTRDcalibDB::Instance()->GetT0Det();
1638 for(int isector =0; isector<AliTRDgeometry::kNsector; isector++){
1639 if(!fTrSec[isector].GetNChambers()) continue;
1640 fTrSec[isector].Init(fReconstructor, cal);
1648 //____________________________________________________________________
1649 void AliTRDtrackerV1::UnloadClusters()
1652 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1655 if(fTracks) fTracks->Delete();
1656 if(fTracklets) fTracklets->Delete();
1658 if(IsClustersOwner()) fClusters->Delete();
1660 // save clusters array in the reconstructor for further use.
1661 if(!fReconstructor->IsWritingClusters()){
1662 AliTRDReconstructor::SetClusters(fClusters);
1663 SetClustersOwner(kFALSE);
1664 } else AliTRDReconstructor::SetClusters(0x0);
1667 for (int i = 0; i < AliTRDgeometry::kNsector; i++) fTrSec[i].Clear();
1669 // Increment the Event Number
1670 AliTRDtrackerDebug::SetEventNumber(AliTRDtrackerDebug::GetEventNumber() + 1);
1673 //_____________________________________________________________________________
1674 Bool_t AliTRDtrackerV1::AdjustSector(AliTRDtrackV1 *track)
1677 // Rotates the track when necessary
1680 Double_t alpha = AliTRDgeometry::GetAlpha();
1681 Double_t y = track->GetY();
1682 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
1685 if (!track->Rotate( alpha)) {
1689 else if (y < -ymax) {
1690 if (!track->Rotate(-alpha)) {
1700 //____________________________________________________________________
1701 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(AliTRDtrackV1 *track, Int_t p, Int_t &idx)
1703 // Find tracklet for TRD track <track>
1712 // Detailed description
1714 idx = track->GetTrackletIndex(p);
1715 AliTRDseedV1 *tracklet = (idx==0xffff) ? 0x0 : (AliTRDseedV1*)fTracklets->UncheckedAt(idx);
1720 //____________________________________________________________________
1721 AliTRDseedV1* AliTRDtrackerV1::SetTracklet(AliTRDseedV1 *tracklet)
1723 // Add this tracklet to the list of tracklets stored in the tracker
1726 // - tracklet : pointer to the tracklet to be added to the list
1729 // - the index of the new tracklet in the tracker tracklets list
1731 // Detailed description
1732 // Build the tracklets list if it is not yet created (late initialization)
1733 // and adds the new tracklet to the list.
1736 fTracklets = new TClonesArray("AliTRDseedV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1737 fTracklets->SetOwner(kTRUE);
1739 Int_t nentries = fTracklets->GetEntriesFast();
1740 return new ((*fTracklets)[nentries]) AliTRDseedV1(*tracklet);
1743 //____________________________________________________________________
1744 AliTRDtrackV1* AliTRDtrackerV1::SetTrack(AliTRDtrackV1 *track)
1746 // Add this track to the list of tracks stored in the tracker
1749 // - track : pointer to the track to be added to the list
1752 // - the pointer added
1754 // Detailed description
1755 // Build the tracks list if it is not yet created (late initialization)
1756 // and adds the new track to the list.
1759 fTracks = new TClonesArray("AliTRDtrackV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1760 fTracks->SetOwner(kTRUE);
1762 Int_t nentries = fTracks->GetEntriesFast();
1763 return new ((*fTracks)[nentries]) AliTRDtrackV1(*track);
1768 //____________________________________________________________________
1769 Int_t AliTRDtrackerV1::Clusters2TracksSM(Int_t sector, AliESDEvent *esd)
1772 // Steer tracking for one SM.
1775 // sector : Array of (SM) propagation layers containing clusters
1776 // esd : The current ESD event. On output it contains the also
1777 // the ESD (TRD) tracks found in this SM.
1780 // Number of tracks found in this TRD supermodule.
1782 // Detailed description
1784 // 1. Unpack AliTRDpropagationLayers objects for each stack.
1785 // 2. Launch stack tracking.
1786 // See AliTRDtrackerV1::Clusters2TracksStack() for details.
1787 // 3. Pack results in the ESD event.
1790 // allocate space for esd tracks in this SM
1791 TClonesArray esdTrackList("AliESDtrack", 2*kMaxTracksStack);
1792 esdTrackList.SetOwner();
1795 Int_t nChambers = 0;
1796 AliTRDtrackingChamber **stack = 0x0, *chamber = 0x0;
1797 for(int istack = 0; istack<AliTRDgeometry::kNstack; istack++){
1798 if(!(stack = fTrSec[sector].GetStack(istack))) continue;
1800 for(int ilayer=0; ilayer<AliTRDgeometry::kNlayer; ilayer++){
1801 if(!(chamber = stack[ilayer])) continue;
1802 if(chamber->GetNClusters() < fgNTimeBins * fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
1804 //AliInfo(Form("sector %d stack %d layer %d clusters %d", sector, istack, ilayer, chamber->GetNClusters()));
1806 if(nChambers < 4) continue;
1807 //AliInfo(Form("Doing stack %d", istack));
1808 nTracks += Clusters2TracksStack(stack, &esdTrackList);
1810 //AliInfo(Form("Found %d tracks in SM %d [%d]\n", nTracks, sector, esd->GetNumberOfTracks()));
1812 for(int itrack=0; itrack<nTracks; itrack++)
1813 esd->AddTrack((AliESDtrack*)esdTrackList[itrack]);
1815 // Reset Track and Candidate Number
1816 AliTRDtrackerDebug::SetCandidateNumber(0);
1817 AliTRDtrackerDebug::SetTrackNumber(0);
1821 //____________________________________________________________________
1822 Int_t AliTRDtrackerV1::Clusters2TracksStack(AliTRDtrackingChamber **stack, TClonesArray *esdTrackList)
1825 // Make tracks in one TRD stack.
1828 // layer : Array of stack propagation layers containing clusters
1829 // esdTrackList : Array of ESD tracks found by the stand alone tracker.
1830 // On exit the tracks found in this stack are appended.
1833 // Number of tracks found in this stack.
1835 // Detailed description
1837 // 1. Find the 3 most useful seeding chambers. See BuildSeedingConfigs() for details.
1838 // 2. Steer AliTRDtrackerV1::MakeSeeds() for 3 seeding layer configurations.
1839 // See AliTRDtrackerV1::MakeSeeds() for more details.
1840 // 3. Arrange track candidates in decreasing order of their quality
1841 // 4. Classify tracks in 5 categories according to:
1842 // a) number of layers crossed
1844 // 5. Sign clusters by tracks in decreasing order of track quality
1845 // 6. Build AliTRDtrack out of seeding tracklets
1847 // 8. Build ESD track and register it to the output list
1850 const AliTRDCalDet *cal = AliTRDcalibDB::Instance()->GetT0Det();
1851 AliTRDtrackingChamber *chamber = 0x0;
1852 AliTRDseedV1 sseed[kMaxTracksStack*6]; // to be initialized
1853 Int_t pars[4]; // MakeSeeds parameters
1855 //Double_t alpha = AliTRDgeometry::GetAlpha();
1856 //Double_t shift = .5 * alpha;
1857 Int_t configs[kNConfigs];
1859 // Build initial seeding configurations
1860 Double_t quality = BuildSeedingConfigs(stack, configs);
1861 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
1862 AliInfo(Form("Plane config %d %d %d Quality %f"
1863 , configs[0], configs[1], configs[2], quality));
1867 // Initialize contors
1868 Int_t ntracks, // number of TRD track candidates
1869 ntracks1, // number of registered TRD tracks/iter
1870 ntracks2 = 0; // number of all registered TRD tracks in stack
1874 Int_t ic = 0; AliTRDtrackingChamber **cIter = &stack[0];
1875 while(ic<kNPlanes && !(*cIter)){ic++; cIter++;}
1876 if(!(*cIter)) return ntracks2;
1877 Int_t istack = fGeom->GetStack((*cIter)->GetDetector());
1880 // Loop over seeding configurations
1881 ntracks = 0; ntracks1 = 0;
1882 for (Int_t iconf = 0; iconf<3; iconf++) {
1883 pars[0] = configs[iconf];
1886 ntracks = MakeSeeds(stack, &sseed[6*ntracks], pars);
1887 if(ntracks == kMaxTracksStack) break;
1889 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1) AliInfo(Form("Candidate TRD tracks %d in iteration %d.", ntracks, fSieveSeeding));
1893 // Sort the seeds according to their quality
1894 Int_t sort[kMaxTracksStack];
1895 TMath::Sort(ntracks, fTrackQuality, sort, kTRUE);
1897 // Initialize number of tracks so far and logic switches
1898 Int_t ntracks0 = esdTrackList->GetEntriesFast();
1899 Bool_t signedTrack[kMaxTracksStack];
1900 Bool_t fakeTrack[kMaxTracksStack];
1901 for (Int_t i=0; i<ntracks; i++){
1902 signedTrack[i] = kFALSE;
1903 fakeTrack[i] = kFALSE;
1905 //AliInfo("Selecting track candidates ...");
1907 // Sieve clusters in decreasing order of track quality
1908 Double_t trackParams[7];
1909 // AliTRDseedV1 *lseed = 0x0;
1910 Int_t jSieve = 0, candidates;
1912 //AliInfo(Form("\t\tITER = %i ", jSieve));
1914 // Check track candidates
1916 for (Int_t itrack = 0; itrack < ntracks; itrack++) {
1917 Int_t trackIndex = sort[itrack];
1918 if (signedTrack[trackIndex] || fakeTrack[trackIndex]) continue;
1921 // Calculate track parameters from tracklets seeds
1926 for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
1927 Int_t jseed = kNPlanes*trackIndex+jLayer;
1928 if(!sseed[jseed].IsOK()) continue;
1929 if (TMath::Abs(sseed[jseed].GetYref(0) / sseed[jseed].GetX0()) < 0.15) findable++;
1931 sseed[jseed].UpdateUsed();
1932 ncl += sseed[jseed].GetN2();
1933 nused += sseed[jseed].GetNUsed();
1937 // Filter duplicated tracks
1939 //printf("Skip %d nused %d\n", trackIndex, nused);
1940 fakeTrack[trackIndex] = kTRUE;
1943 if (Float_t(nused)/ncl >= .25){
1944 //printf("Skip %d nused/ncl >= .25\n", trackIndex);
1945 fakeTrack[trackIndex] = kTRUE;
1950 Bool_t skip = kFALSE;
1953 if(nlayers < 6) {skip = kTRUE; break;}
1954 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
1958 if(nlayers < findable){skip = kTRUE; break;}
1959 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -4.){skip = kTRUE; break;}
1963 if ((nlayers == findable) || (nlayers == 6)) { skip = kTRUE; break;}
1964 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -6.0){skip = kTRUE; break;}
1968 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
1972 if (nlayers == 3){skip = kTRUE; break;}
1973 //if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) - nused/(nlayers-3.0) < -15.0){skip = kTRUE; break;}
1978 //printf("REJECTED : %d [%d] nlayers %d trackQuality = %e nused %d\n", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused);
1981 signedTrack[trackIndex] = kTRUE;
1985 AliTRDcluster *cl = 0x0; Int_t clusterIndex = -1;
1986 for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
1987 Int_t jseed = kNPlanes*trackIndex+jLayer;
1988 if(!sseed[jseed].IsOK()) continue;
1989 if(TMath::Abs(sseed[jseed].GetYfit(1) - sseed[jseed].GetYfit(1)) >= .2) continue; // check this condition with Marian
1990 sseed[jseed].UseClusters();
1993 while(!(cl = sseed[jseed].GetClusters(ic))) ic++;
1994 clusterIndex = sseed[jseed].GetIndexes(ic);
2000 // Build track parameters
2001 AliTRDseedV1 *lseed =&sseed[trackIndex*6];
2003 while(idx<3 && !lseed->IsOK()) {
2007 Double_t x = lseed->GetX0();// - 3.5;
2008 trackParams[0] = x; //NEW AB
2009 trackParams[1] = lseed->GetYref(0); // lseed->GetYat(x);
2010 trackParams[2] = lseed->GetZref(0); // lseed->GetZat(x);
2011 trackParams[3] = TMath::Sin(TMath::ATan(lseed->GetYref(1)));
2012 trackParams[4] = lseed->GetZref(1) / TMath::Sqrt(1. + lseed->GetYref(1) * lseed->GetYref(1));
2013 trackParams[5] = lseed->GetC();
2014 Int_t ich = 0; while(!(chamber = stack[ich])) ich++;
2015 trackParams[6] = fGeom->GetSector(chamber->GetDetector());/* *alpha+shift; // Supermodule*/
2017 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2018 AliInfo(Form("Track %d [%d] nlayers %d trackQuality = %e nused %d, yref = %3.3f", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused, trackParams[1]));
2020 Int_t nclusters = 0;
2021 AliTRDseedV1 *dseed[6];
2023 // Build track label - what happens if measured data ???
2028 Int_t labelsall[1000];
2029 Int_t nlabelsall = 0;
2030 Int_t naccepted = 0;
2032 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
2033 Int_t jseed = kNPlanes*trackIndex+iLayer;
2034 dseed[iLayer] = new AliTRDseedV1(sseed[jseed]);
2035 dseed[iLayer]->SetOwner();
2036 nclusters += sseed[jseed].GetN2();
2037 if(!sseed[jseed].IsOK()) continue;
2038 for(int ilab=0; ilab<2; ilab++){
2039 if(sseed[jseed].GetLabels(ilab) < 0) continue;
2040 labels[nlab] = sseed[jseed].GetLabels(ilab);
2045 for (Int_t itime = 0; itime < fgNTimeBins; itime++) {
2046 if(!sseed[jseed].IsUsable(itime)) continue;
2048 Int_t tindex = 0, ilab = 0;
2049 while(ilab<3 && (tindex = sseed[jseed].GetClusters(itime)->GetLabel(ilab)) >= 0){
2050 labelsall[nlabelsall++] = tindex;
2055 Freq(nlab,labels,outlab,kFALSE);
2056 Int_t label = outlab[0];
2057 Int_t frequency = outlab[1];
2058 Freq(nlabelsall,labelsall,outlab,kFALSE);
2059 Int_t label1 = outlab[0];
2060 Int_t label2 = outlab[2];
2061 Float_t fakeratio = (naccepted - outlab[1]) / Float_t(naccepted);
2063 //Int_t eventNrInFile = esd->GetEventNumberInFile();
2064 //AliInfo(Form("Number of clusters %d.", nclusters));
2065 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2066 Int_t trackNumber = AliTRDtrackerDebug::GetTrackNumber();
2067 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2068 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2069 cstreamer << "Clusters2TracksStack"
2070 << "EventNumber=" << eventNumber
2071 << "TrackNumber=" << trackNumber
2072 << "CandidateNumber=" << candidateNumber
2073 << "Iter=" << fSieveSeeding
2074 << "Like=" << fTrackQuality[trackIndex]
2075 << "S0.=" << dseed[0]
2076 << "S1.=" << dseed[1]
2077 << "S2.=" << dseed[2]
2078 << "S3.=" << dseed[3]
2079 << "S4.=" << dseed[4]
2080 << "S5.=" << dseed[5]
2081 << "p0=" << trackParams[0]
2082 << "p1=" << trackParams[1]
2083 << "p2=" << trackParams[2]
2084 << "p3=" << trackParams[3]
2085 << "p4=" << trackParams[4]
2086 << "p5=" << trackParams[5]
2087 << "p6=" << trackParams[6]
2088 << "Label=" << label
2089 << "Label1=" << label1
2090 << "Label2=" << label2
2091 << "FakeRatio=" << fakeratio
2092 << "Freq=" << frequency
2094 << "NLayers=" << nlayers
2095 << "Findable=" << findable
2096 << "NUsed=" << nused
2100 AliTRDtrackV1 *track = MakeTrack(&sseed[trackIndex*kNPlanes], trackParams);
2102 AliWarning("Fail to build a TRD Track.");
2106 //AliInfo("End of MakeTrack()");
2107 AliESDtrack *esdTrack = new ((*esdTrackList)[ntracks0++]) AliESDtrack();
2108 esdTrack->UpdateTrackParams(track, AliESDtrack::kTRDout);
2109 esdTrack->SetLabel(track->GetLabel());
2110 track->UpdateESDtrack(esdTrack);
2111 // write ESD-friends if neccessary
2112 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0){
2113 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(*track);
2114 calibTrack->SetOwner();
2115 esdTrack->AddCalibObject(calibTrack);
2118 AliTRDtrackerDebug::SetTrackNumber(AliTRDtrackerDebug::GetTrackNumber() + 1);
2122 } while(jSieve<5 && candidates); // end track candidates sieve
2123 if(!ntracks1) break;
2125 // increment counters
2126 ntracks2 += ntracks1;
2128 if(fReconstructor->IsHLT()) break;
2131 // Rebuild plane configurations and indices taking only unused clusters into account
2132 quality = BuildSeedingConfigs(stack, configs);
2133 if(quality < 1.E-7) break; //fReconstructor->GetRecoParam() ->GetPlaneQualityThreshold()) break;
2135 for(Int_t ip = 0; ip < kNPlanes; ip++){
2136 if(!(chamber = stack[ip])) continue;
2137 chamber->Build(fGeom, cal);//Indices(fSieveSeeding);
2140 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2141 AliInfo(Form("Sieve level %d Plane config %d %d %d Quality %f", fSieveSeeding, configs[0], configs[1], configs[2], quality));
2143 } while(fSieveSeeding<10); // end stack clusters sieve
2147 //AliInfo(Form("Registered TRD tracks %d in stack %d.", ntracks2, pars[1]));
2152 //___________________________________________________________________
2153 Double_t AliTRDtrackerV1::BuildSeedingConfigs(AliTRDtrackingChamber **stack, Int_t *configs)
2156 // Assign probabilities to chambers according to their
2157 // capability of producing seeds.
2161 // layers : Array of stack propagation layers for all 6 chambers in one stack
2162 // configs : On exit array of configuration indexes (see GetSeedingConfig()
2163 // for details) in the decreasing order of their seeding probabilities.
2167 // Return top configuration quality
2169 // Detailed description:
2171 // To each chamber seeding configuration (see GetSeedingConfig() for
2172 // the list of all configurations) one defines 2 quality factors:
2173 // - an apriori topological quality (see GetSeedingConfig() for details) and
2174 // - a data quality based on the uniformity of the distribution of
2175 // clusters over the x range (time bins population). See CookChamberQA() for details.
2176 // The overall chamber quality is given by the product of this 2 contributions.
2179 Double_t chamberQ[kNPlanes];
2180 AliTRDtrackingChamber *chamber = 0x0;
2181 for(int iplane=0; iplane<kNPlanes; iplane++){
2182 if(!(chamber = stack[iplane])) continue;
2183 chamberQ[iplane] = (chamber = stack[iplane]) ? chamber->GetQuality() : 0.;
2186 Double_t tconfig[kNConfigs];
2188 for(int iconf=0; iconf<kNConfigs; iconf++){
2189 GetSeedingConfig(iconf, planes);
2190 tconfig[iconf] = fgTopologicQA[iconf];
2191 for(int iplane=0; iplane<4; iplane++) tconfig[iconf] *= chamberQ[planes[iplane]];
2194 TMath::Sort((Int_t)kNConfigs, tconfig, configs, kTRUE);
2195 // AliInfo(Form("q[%d] = %f", configs[0], tconfig[configs[0]]));
2196 // AliInfo(Form("q[%d] = %f", configs[1], tconfig[configs[1]]));
2197 // AliInfo(Form("q[%d] = %f", configs[2], tconfig[configs[2]]));
2199 return tconfig[configs[0]];
2202 //____________________________________________________________________
2203 Int_t AliTRDtrackerV1::MakeSeeds(AliTRDtrackingChamber **stack, AliTRDseedV1 *sseed, Int_t *ipar)
2206 // Make tracklet seeds in the TRD stack.
2209 // layers : Array of stack propagation layers containing clusters
2210 // sseed : Array of empty tracklet seeds. On exit they are filled.
2211 // ipar : Control parameters:
2212 // ipar[0] -> seeding chambers configuration
2213 // ipar[1] -> stack index
2214 // ipar[2] -> number of track candidates found so far
2217 // Number of tracks candidates found.
2219 // Detailed description
2221 // The following steps are performed:
2222 // 1. Select seeding layers from seeding chambers
2223 // 2. Select seeding clusters from the seeding AliTRDpropagationLayerStack.
2224 // The clusters are taken from layer 3, layer 0, layer 1 and layer 2, in
2225 // this order. The parameters controling the range of accepted clusters in
2226 // layer 0, 1, and 2 are defined in AliTRDchamberTimeBin::BuildCond().
2227 // 3. Helix fit of the cluster set. (see AliTRDtrackerFitter::FitRieman(AliTRDcluster**))
2228 // 4. Initialize seeding tracklets in the seeding chambers.
2230 // Chi2 in the Y direction less than threshold ... (1./(3. - sLayer))
2231 // Chi2 in the Z direction less than threshold ... (1./(3. - sLayer))
2232 // 6. Attach clusters to seeding tracklets and find linear approximation of
2233 // the tracklet (see AliTRDseedV1::AttachClustersIter()). The number of used
2234 // clusters used by current seeds should not exceed ... (25).
2236 // All 4 seeding tracklets should be correctly constructed (see
2237 // AliTRDseedV1::AttachClustersIter())
2238 // 8. Helix fit of the seeding tracklets
2240 // Likelihood calculation of the fit. (See AliTRDtrackerV1::CookLikelihood() for details)
2241 // 10. Extrapolation of the helix fit to the other 2 chambers:
2242 // a) Initialization of extrapolation tracklet with fit parameters
2243 // b) Helix fit of tracklets
2244 // c) Attach clusters and linear interpolation to extrapolated tracklets
2245 // d) Helix fit of tracklets
2246 // 11. Improve seeding tracklets quality by reassigning clusters.
2247 // See AliTRDtrackerV1::ImproveSeedQuality() for details.
2248 // 12. Helix fit of all 6 seeding tracklets and chi2 calculation
2249 // 13. Hyperplane fit and track quality calculation. See AliTRDtrackerFitter::FitHyperplane() for details.
2250 // 14. Cooking labels for tracklets. Should be done only for MC
2251 // 15. Register seeds.
2254 AliTRDtrackingChamber *chamber = 0x0;
2255 AliTRDcluster *c[kNSeedPlanes] = {0x0, 0x0, 0x0, 0x0}; // initilize seeding clusters
2256 AliTRDseedV1 *cseed = &sseed[0]; // initialize tracklets for first track
2257 Int_t ncl, mcl; // working variable for looping over clusters
2258 Int_t index[AliTRDchamberTimeBin::kMaxClustersLayer], jndex[AliTRDchamberTimeBin::kMaxClustersLayer];
2260 // chi2[0] = tracklet chi2 on the Z direction
2261 // chi2[1] = tracklet chi2 on the R direction
2264 // Default positions for the anode wire in all 6 Layers in case of a stack with missing clusters
2265 // Positions taken using cosmic data taken with SM3 after rebuild
2266 Double_t x_def[kNPlanes] = {300.2, 312.8, 325.4, 338.0, 350.6, 363.2};
2268 // this should be data member of AliTRDtrack
2269 Double_t seedQuality[kMaxTracksStack];
2271 // unpack control parameters
2272 Int_t config = ipar[0];
2273 Int_t ntracks = ipar[1];
2274 Int_t istack = ipar[2];
2275 Int_t planes[kNSeedPlanes]; GetSeedingConfig(config, planes);
2276 Int_t planesExt[kNPlanes-kNSeedPlanes]; GetExtrapolationConfig(config, planesExt);
2279 // Init chambers geometry
2280 Double_t hL[kNPlanes]; // Tilting angle
2281 Float_t padlength[kNPlanes]; // pad lenghts
2282 AliTRDpadPlane *pp = 0x0;
2283 for(int iplane=0; iplane<kNPlanes; iplane++){
2284 pp = fGeom->GetPadPlane(iplane, istack);
2285 hL[iplane] = TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle());
2286 padlength[iplane] = pp->GetLengthIPad();
2289 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2290 AliInfo(Form("Making seeds Stack[%d] Config[%d] Tracks[%d]...", istack, config, ntracks));
2293 // Build seeding layers
2296 for(int isl=0; isl<kNSeedPlanes; isl++){
2297 if(!(chamber = stack[planes[isl]])) continue;
2298 if(!chamber->GetSeedingLayer(fSeedTB[isl], fGeom, fReconstructor)) continue;
2301 if(nlayers < 4) return ntracks;
2304 // Start finding seeds
2305 Double_t cond0[4], cond1[4], cond2[4];
2307 while((c[3] = (*fSeedTB[3])[icl++])){
2309 fSeedTB[0]->BuildCond(c[3], cond0, 0);
2310 fSeedTB[0]->GetClusters(cond0, index, ncl);
2311 //printf("Found c[3] candidates 0 %d\n", ncl);
2314 c[0] = (*fSeedTB[0])[index[jcl++]];
2316 Double_t dx = c[3]->GetX() - c[0]->GetX();
2317 Double_t theta = (c[3]->GetZ() - c[0]->GetZ())/dx;
2318 Double_t phi = (c[3]->GetY() - c[0]->GetY())/dx;
2319 fSeedTB[1]->BuildCond(c[0], cond1, 1, theta, phi);
2320 fSeedTB[1]->GetClusters(cond1, jndex, mcl);
2321 //printf("Found c[0] candidates 1 %d\n", mcl);
2325 c[1] = (*fSeedTB[1])[jndex[kcl++]];
2327 fSeedTB[2]->BuildCond(c[1], cond2, 2, theta, phi);
2328 c[2] = fSeedTB[2]->GetNearestCluster(cond2);
2329 //printf("Found c[1] candidate 2 %p\n", c[2]);
2332 // AliInfo("Seeding clusters found. Building seeds ...");
2333 // 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());
2335 for (Int_t il = 0; il < kNPlanes; il++) cseed[il].Reset();
2339 AliTRDseedV1 *tseed = &cseed[0];
2340 AliTRDtrackingChamber **cIter = &stack[0];
2341 for(int iLayer=0; iLayer<kNPlanes; iLayer++, tseed++, cIter++){
2342 tseed->SetDetector((*cIter) ? (*cIter)->GetDetector() : -1);
2343 tseed->SetTilt(hL[iLayer]);
2344 tseed->SetPadLength(padlength[iLayer]);
2345 tseed->SetReconstructor(fReconstructor);
2346 tseed->SetX0((*cIter) ? (*cIter)->GetX() : x_def[iLayer]);
2347 tseed->Init(GetRiemanFitter());
2350 Bool_t isFake = kFALSE;
2351 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2352 if (c[0]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2353 if (c[1]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2354 if (c[2]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2357 for(Int_t l = 0; l < kNSeedPlanes; l++) xpos[l] = fSeedTB[l]->GetX();
2359 for(int il=0; il<4; il++) yref[il] = cseed[planes[il]].GetYref(0);
2360 Int_t ll = c[3]->GetLabel(0);
2361 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2362 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2363 AliRieman *rim = GetRiemanFitter();
2364 TTreeSRedirector &cs0 = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2366 <<"EventNumber=" << eventNumber
2367 <<"CandidateNumber=" << candidateNumber
2368 <<"isFake=" << isFake
2369 <<"config=" << config
2371 <<"chi2z=" << chi2[0]
2372 <<"chi2y=" << chi2[1]
2373 <<"Y2exp=" << cond2[0]
2374 <<"Z2exp=" << cond2[1]
2375 <<"X0=" << xpos[0] //layer[sLayer]->GetX()
2376 <<"X1=" << xpos[1] //layer[sLayer + 1]->GetX()
2377 <<"X2=" << xpos[2] //layer[sLayer + 2]->GetX()
2378 <<"X3=" << xpos[3] //layer[sLayer + 3]->GetX()
2379 <<"yref0=" << yref[0]
2380 <<"yref1=" << yref[1]
2381 <<"yref2=" << yref[2]
2382 <<"yref3=" << yref[3]
2387 <<"Seed0.=" << &cseed[planes[0]]
2388 <<"Seed1.=" << &cseed[planes[1]]
2389 <<"Seed2.=" << &cseed[planes[2]]
2390 <<"Seed3.=" << &cseed[planes[3]]
2391 <<"RiemanFitter.=" << rim
2394 if(chi2[0] > fReconstructor->GetRecoParam() ->GetChi2Z()/*7./(3. - sLayer)*//*iter*/){
2395 // //AliInfo(Form("Failed chi2 filter on chi2Z [%f].", chi2[0]));
2396 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2399 if(chi2[1] > fReconstructor->GetRecoParam() ->GetChi2Y()/*1./(3. - sLayer)*//*iter*/){
2400 // //AliInfo(Form("Failed chi2 filter on chi2Y [%f].", chi2[1]));
2401 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2404 //AliInfo("Passed chi2 filter.");
2406 // try attaching clusters to tracklets
2409 for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
2410 Int_t jLayer = planes[iLayer];
2411 if(!cseed[jLayer].AttachClustersIter(stack[jLayer], 5., kFALSE, c[iLayer])) continue;
2412 nUsedCl += cseed[jLayer].GetNUsed();
2413 if(nUsedCl > 25) break;
2417 if(mlayers < kNSeedPlanes){
2418 //AliInfo(Form("Failed updating all seeds %d [%d].", mlayers, kNSeedPlanes));
2419 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2423 // temporary exit door for the HLT
2424 if(fReconstructor->IsHLT()){
2425 // attach clusters to extrapolation chambers
2426 for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
2427 Int_t jLayer = planesExt[iLayer];
2428 if(!(chamber = stack[jLayer])) continue;
2429 cseed[jLayer].AttachClustersIter(chamber, 1000.);
2431 fTrackQuality[ntracks] = 1.; // dummy value
2433 if(ntracks == kMaxTracksStack) return ntracks;
2439 // fit tracklets and cook likelihood
2440 FitTiltedRieman(&cseed[0], kTRUE);// Update Seeds and calculate Likelihood
2441 Double_t like = CookLikelihood(&cseed[0], planes); // to be checked
2443 if (TMath::Log(1.E-9 + like) < fReconstructor->GetRecoParam() ->GetTrackLikelihood()){
2444 //AliInfo(Form("Failed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2445 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2448 //AliInfo(Form("Passed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2450 // book preliminary results
2451 seedQuality[ntracks] = like;
2452 fSeedLayer[ntracks] = config;/*sLayer;*/
2454 // attach clusters to the extrapolation seeds
2455 Int_t nusedf = 0; // debug value
2456 for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
2457 Int_t jLayer = planesExt[iLayer];
2458 if(!(chamber = stack[jLayer])) continue;
2460 // fit extrapolated seed
2461 if ((jLayer == 0) && !(cseed[1].IsOK())) continue;
2462 if ((jLayer == 5) && !(cseed[4].IsOK())) continue;
2463 AliTRDseedV1 pseed = cseed[jLayer];
2464 if(!pseed.AttachClustersIter(chamber, 1000.)) continue;
2465 cseed[jLayer] = pseed;
2466 nusedf += cseed[jLayer].GetNUsed(); // debug value
2467 FitTiltedRieman(cseed, kTRUE);
2470 // AliInfo("Extrapolation done.");
2471 // Debug Stream containing all the 6 tracklets
2472 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2473 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2474 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2475 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2476 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2477 cstreamer << "MakeSeeds1"
2478 << "EventNumber=" << eventNumber
2479 << "CandidateNumber=" << candidateNumber
2480 << "S0.=" << &cseed[0]
2481 << "S1.=" << &cseed[1]
2482 << "S2.=" << &cseed[2]
2483 << "S3.=" << &cseed[3]
2484 << "S4.=" << &cseed[4]
2485 << "S5.=" << &cseed[5]
2486 << "FitterT.=" << tiltedRieman
2490 if(fReconstructor->GetRecoParam()->HasImproveTracklets() && ImproveSeedQuality(stack, cseed) < 4){
2491 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2494 //AliInfo("Improve seed quality done.");
2496 // fit full track and cook likelihoods
2497 // Double_t curv = FitRieman(&cseed[0], chi2);
2498 // Double_t chi2ZF = chi2[0] / TMath::Max((mlayers - 3.), 1.);
2499 // Double_t chi2RF = chi2[1] / TMath::Max((mlayers - 3.), 1.);
2501 // do the final track fitting (Once with vertex constraint and once without vertex constraint)
2502 Double_t chi2Vals[3];
2503 chi2Vals[0] = FitTiltedRieman(&cseed[0], kFALSE);
2504 if(fReconstructor->GetRecoParam()->IsVertexConstrained())
2505 chi2Vals[1] = FitTiltedRiemanConstraint(&cseed[0], GetZ()); // Do Vertex Constrained fit if desired
2508 chi2Vals[2] = GetChi2Z(&cseed[0]) / TMath::Max((mlayers - 3.), 1.);
2509 // Chi2 definitions in testing stage
2510 //chi2Vals[2] = GetChi2ZTest(&cseed[0]);
2511 fTrackQuality[ntracks] = CalculateTrackLikelihood(&cseed[0], &chi2Vals[0]);
2512 //AliInfo("Hyperplane fit done\n");
2514 // finalize tracklets
2518 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2519 if (!cseed[iLayer].IsOK()) continue;
2521 if (cseed[iLayer].GetLabels(0) >= 0) {
2522 labels[nlab] = cseed[iLayer].GetLabels(0);
2526 if (cseed[iLayer].GetLabels(1) >= 0) {
2527 labels[nlab] = cseed[iLayer].GetLabels(1);
2531 Freq(nlab,labels,outlab,kFALSE);
2532 Int_t label = outlab[0];
2533 Int_t frequency = outlab[1];
2534 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2535 cseed[iLayer].SetFreq(frequency);
2536 cseed[iLayer].SetChi2Z(chi2[1]);
2539 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2540 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2541 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2542 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2543 TLinearFitter *fitterTC = GetTiltedRiemanFitterConstraint();
2544 TLinearFitter *fitterT = GetTiltedRiemanFitter();
2546 for(Int_t iseed = 0; iseed < kNPlanes; iseed++){
2547 ncls += cseed[iseed].IsOK() ? cseed[iseed].GetN2() : 0;
2549 cstreamer << "MakeSeeds2"
2550 << "EventNumber=" << eventNumber
2551 << "CandidateNumber=" << candidateNumber
2552 << "Chi2TR=" << chi2Vals[0]
2553 << "Chi2TC=" << chi2Vals[1]
2554 << "Nlayers=" << mlayers
2555 << "NClusters=" << ncls
2556 << "NUsedS=" << nUsedCl
2557 << "NUsed=" << nusedf
2559 << "S0.=" << &cseed[0]
2560 << "S1.=" << &cseed[1]
2561 << "S2.=" << &cseed[2]
2562 << "S3.=" << &cseed[3]
2563 << "S4.=" << &cseed[4]
2564 << "S5.=" << &cseed[5]
2565 << "Label=" << label
2566 << "Freq=" << frequency
2567 << "FitterT.=" << fitterT
2568 << "FitterTC.=" << fitterTC
2573 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2574 if(ntracks == kMaxTracksStack){
2575 AliWarning(Form("Number of seeds reached maximum allowed (%d) in stack.", kMaxTracksStack));
2586 //_____________________________________________________________________________
2587 AliTRDtrackV1* AliTRDtrackerV1::MakeTrack(AliTRDseedV1 *seeds, Double_t *params)
2590 // Build a TRD track out of tracklet candidates
2593 // seeds : array of tracklets
2594 // params : track parameters (see MakeSeeds() function body for a detailed description)
2599 // Detailed description
2601 // To be discussed with Marian !!
2605 Double_t alpha = AliTRDgeometry::GetAlpha();
2606 Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
2610 c[ 1] = 0.0; c[ 2] = 2.0;
2611 c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
2612 c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
2613 c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
2615 AliTRDtrackV1 track(seeds, ¶ms[1], c, params[0], params[6]*alpha+shift);
2616 track.PropagateTo(params[0]-5.0);
2617 if(fReconstructor->IsHLT()){
2618 AliTRDseedV1 *ptrTracklet = 0x0;
2619 for(Int_t ip=0; ip<kNPlanes; ip++){
2620 track.UnsetTracklet(ip);
2621 ptrTracklet = SetTracklet(&seeds[ip]);
2622 track.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
2624 AliTRDtrackV1 *ptrTrack = SetTrack(&track);
2625 ptrTrack->SetReconstructor(fReconstructor);
2629 track.ResetCovariance(1);
2630 Int_t nc = TMath::Abs(FollowBackProlongation(track));
2631 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 5){
2632 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2633 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2634 Double_t p[5]; // Track Params for the Debug Stream
2635 track.GetExternalParameters(params[0], p);
2636 TTreeSRedirector &cs = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2638 << "EventNumber=" << eventNumber
2639 << "CandidateNumber=" << candidateNumber
2641 << "X=" << params[0]
2647 << "Yin=" << params[1]
2648 << "Zin=" << params[2]
2649 << "snpin=" << params[3]
2650 << "tndin=" << params[4]
2651 << "crvin=" << params[5]
2652 << "track.=" << &track
2655 if (nc < 30) return 0x0;
2657 AliTRDtrackV1 *ptrTrack = SetTrack(&track);
2658 ptrTrack->SetReconstructor(fReconstructor);
2659 ptrTrack->CookLabel(.9);
2661 // computes PID for track
2662 ptrTrack->CookPID();
2663 // update calibration references using this track
2664 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
2666 AliInfo("Could not get Calibra instance\n");
2667 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(ptrTrack);
2673 //____________________________________________________________________
2674 Int_t AliTRDtrackerV1::ImproveSeedQuality(AliTRDtrackingChamber **stack, AliTRDseedV1 *cseed)
2677 // Sort tracklets according to "quality" and try to "improve" the first 4 worst
2680 // layers : Array of propagation layers for a stack/supermodule
2681 // cseed : Array of 6 seeding tracklets which has to be improved
2684 // cssed : Improved seeds
2686 // Detailed description
2688 // Iterative procedure in which new clusters are searched for each
2689 // tracklet seed such that the seed quality (see AliTRDseed::GetQuality())
2690 // can be maximized. If some optimization is found the old seeds are replaced.
2695 // make a local working copy
2696 AliTRDtrackingChamber *chamber = 0x0;
2697 AliTRDseedV1 bseed[6];
2699 for (Int_t jLayer = 0; jLayer < 6; jLayer++) bseed[jLayer] = cseed[jLayer];
2701 Float_t lastquality = 10000.0;
2702 Float_t lastchi2 = 10000.0;
2703 Float_t chi2 = 1000.0;
2705 for (Int_t iter = 0; iter < 4; iter++) {
2706 Float_t sumquality = 0.0;
2707 Float_t squality[6];
2708 Int_t sortindexes[6];
2710 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2711 squality[jLayer] = bseed[jLayer].IsOK() ? bseed[jLayer].GetQuality(kTRUE) : 1000.;
2712 sumquality += squality[jLayer];
2714 if ((sumquality >= lastquality) || (chi2 > lastchi2)) break;
2717 lastquality = sumquality;
2719 if (iter > 0) for (Int_t jLayer = 0; jLayer < 6; jLayer++) cseed[jLayer] = bseed[jLayer];
2721 TMath::Sort(6, squality, sortindexes, kFALSE);
2722 for (Int_t jLayer = 5; jLayer > 1; jLayer--) {
2723 Int_t bLayer = sortindexes[jLayer];
2724 if(!(chamber = stack[bLayer])) continue;
2725 bseed[bLayer].AttachClustersIter(chamber, squality[bLayer], kTRUE);
2726 if(bseed[bLayer].IsOK()) nLayers++;
2729 chi2 = FitTiltedRieman(bseed, kTRUE);
2730 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 7){
2731 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2732 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2733 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2734 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2735 cstreamer << "ImproveSeedQuality"
2736 << "EventNumber=" << eventNumber
2737 << "CandidateNumber=" << candidateNumber
2738 << "Iteration=" << iter
2739 << "S0.=" << &bseed[0]
2740 << "S1.=" << &bseed[1]
2741 << "S2.=" << &bseed[2]
2742 << "S3.=" << &bseed[3]
2743 << "S4.=" << &bseed[4]
2744 << "S5.=" << &bseed[5]
2745 << "FitterT.=" << tiltedRieman
2750 // we are sure that at least 2 tracklets are OK !
2754 //_________________________________________________________________________
2755 Double_t AliTRDtrackerV1::CalculateTrackLikelihood(AliTRDseedV1 *tracklets, Double_t *chi2){
2757 // Calculates the Track Likelihood value. This parameter serves as main quality criterion for
2758 // the track selection
2759 // The likelihood value containes:
2760 // - The chi2 values from the both fitters and the chi2 values in z-direction from a linear fit
2761 // - The Sum of the Parameter |slope_ref - slope_fit|/Sigma of the tracklets
2762 // For all Parameters an exponential dependency is used
2764 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
2765 // - Array of chi2 values:
2766 // * Non-Constrained Tilted Riemann fit
2767 // * Vertex-Constrained Tilted Riemann fit
2768 // * z-Direction from Linear fit
2769 // Output: - The calculated track likelihood
2774 Double_t sumdaf = 0, nLayers = 0;
2775 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
2776 if(!tracklets[iLayer].IsOK()) continue;
2777 sumdaf += TMath::Abs((tracklets[iLayer].GetYfit(1) - tracklets[iLayer].GetYref(1))/ tracklets[iLayer].GetSigmaY2());
2780 sumdaf /= Float_t (nLayers - 2.0);
2782 Double_t likeChi2Z = TMath::Exp(-chi2[2] * 0.14); // Chi2Z
2783 Double_t likeChi2TC = (fReconstructor->GetRecoParam() ->IsVertexConstrained()) ?
2784 TMath::Exp(-chi2[1] * 0.677) : 1; // Constrained Tilted Riemann
2785 Double_t likeChi2TR = TMath::Exp(-chi2[0] * 0.78); // Non-constrained Tilted Riemann
2786 Double_t likeAF = TMath::Exp(-sumdaf * 3.23);
2787 Double_t trackLikelihood = likeChi2Z * likeChi2TR * likeAF;
2789 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2790 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2791 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2792 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2793 cstreamer << "CalculateTrackLikelihood0"
2794 << "EventNumber=" << eventNumber
2795 << "CandidateNumber=" << candidateNumber
2796 << "LikeChi2Z=" << likeChi2Z
2797 << "LikeChi2TR=" << likeChi2TR
2798 << "LikeChi2TC=" << likeChi2TC
2799 << "LikeAF=" << likeAF
2800 << "TrackLikelihood=" << trackLikelihood
2804 return trackLikelihood;
2807 //____________________________________________________________________
2808 Double_t AliTRDtrackerV1::CookLikelihood(AliTRDseedV1 *cseed, Int_t planes[4])
2811 // Calculate the probability of this track candidate.
2814 // cseeds : array of candidate tracklets
2815 // planes : array of seeding planes (see seeding configuration)
2816 // chi2 : chi2 values (on the Z and Y direction) from the rieman fit of the track.
2821 // Detailed description
2823 // The track quality is estimated based on the following 4 criteria:
2824 // 1. precision of the rieman fit on the Y direction (likea)
2825 // 2. chi2 on the Y direction (likechi2y)
2826 // 3. chi2 on the Z direction (likechi2z)
2827 // 4. number of attached clusters compared to a reference value
2828 // (see AliTRDrecoParam::fkFindable) (likeN)
2830 // The distributions for each type of probabilities are given below as of
2831 // (date). They have to be checked to assure consistency of estimation.
2834 // ratio of the total number of clusters/track which are expected to be found by the tracker.
2835 const AliTRDrecoParam *fRecoPars = fReconstructor->GetRecoParam();
2837 Double_t chi2y = GetChi2Y(&cseed[0]);
2838 Double_t chi2z = GetChi2Z(&cseed[0]);
2840 Float_t nclusters = 0.;
2841 Double_t sumda = 0.;
2842 for(UChar_t ilayer = 0; ilayer < 4; ilayer++){
2843 Int_t jlayer = planes[ilayer];
2844 nclusters += cseed[jlayer].GetN2();
2845 sumda += TMath::Abs(cseed[jlayer].GetYfitR(1) - cseed[jlayer].GetYref(1));
2849 Double_t likea = TMath::Exp(-sumda * fRecoPars->GetPhiSlope());
2850 Double_t likechi2y = 0.0000000001;
2851 if (fReconstructor->IsCosmic() || chi2y < fRecoPars->GetChi2YCut()) likechi2y += TMath::Exp(-TMath::Sqrt(chi2y) * fRecoPars->GetChi2YSlope());
2852 Double_t likechi2z = TMath::Exp(-chi2z * fRecoPars->GetChi2ZSlope());
2853 Double_t likeN = TMath::Exp(-(fRecoPars->GetNMeanClusters() - nclusters) / fRecoPars->GetNSigmaClusters());
2854 Double_t like = likea * likechi2y * likechi2z * likeN;
2856 // 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));
2857 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2858 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2859 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2860 Int_t nTracklets = 0; Float_t mean_ncls = 0;
2861 for(Int_t iseed=0; iseed < kNPlanes; iseed++){
2862 if(!cseed[iseed].IsOK()) continue;
2864 mean_ncls += cseed[iseed].GetN2();
2866 if(nTracklets) mean_ncls /= nTracklets;
2867 // The Debug Stream contains the seed
2868 TTreeSRedirector &cstreamer = *fReconstructor->GetDebugStream(AliTRDReconstructor::kTracker);
2869 cstreamer << "CookLikelihood"
2870 << "EventNumber=" << eventNumber
2871 << "CandidateNumber=" << candidateNumber
2872 << "tracklet0.=" << &cseed[0]
2873 << "tracklet1.=" << &cseed[1]
2874 << "tracklet2.=" << &cseed[2]
2875 << "tracklet3.=" << &cseed[3]
2876 << "tracklet4.=" << &cseed[4]
2877 << "tracklet5.=" << &cseed[5]
2878 << "sumda=" << sumda
2879 << "chi2y=" << chi2y
2880 << "chi2z=" << chi2z
2881 << "likea=" << likea
2882 << "likechi2y=" << likechi2y
2883 << "likechi2z=" << likechi2z
2884 << "nclusters=" << nclusters
2885 << "likeN=" << likeN
2887 << "meanncls=" << mean_ncls
2894 //____________________________________________________________________
2895 void AliTRDtrackerV1::GetSeedingConfig(Int_t iconfig, Int_t planes[4])
2898 // Map seeding configurations to detector planes.
2901 // iconfig : configuration index
2902 // planes : member planes of this configuration. On input empty.
2905 // planes : contains the planes which are defining the configuration
2907 // Detailed description
2909 // Here is the list of seeding planes configurations together with
2910 // their topological classification:
2928 // The topologic quality is modeled as follows:
2929 // 1. The general model is define by the equation:
2930 // p(conf) = exp(-conf/2)
2931 // 2. According to the topologic classification, configurations from the same
2932 // class are assigned the agerage value over the model values.
2933 // 3. Quality values are normalized.
2935 // The topologic quality distribution as function of configuration is given below:
2937 // <img src="gif/topologicQA.gif">
2942 case 0: // 5432 TQ 0
2948 case 1: // 4321 TQ 0
2954 case 2: // 3210 TQ 0
2960 case 3: // 5321 TQ 1
2966 case 4: // 4210 TQ 1
2972 case 5: // 5431 TQ 1
2978 case 6: // 4320 TQ 1
2984 case 7: // 5430 TQ 2
2990 case 8: // 5210 TQ 2
2996 case 9: // 5421 TQ 3
3002 case 10: // 4310 TQ 3
3008 case 11: // 5410 TQ 4
3014 case 12: // 5420 TQ 5
3020 case 13: // 5320 TQ 5
3026 case 14: // 5310 TQ 5
3035 //____________________________________________________________________
3036 void AliTRDtrackerV1::GetExtrapolationConfig(Int_t iconfig, Int_t planes[2])
3039 // Returns the extrapolation planes for a seeding configuration.
3042 // iconfig : configuration index
3043 // planes : planes which are not in this configuration. On input empty.
3046 // planes : contains the planes which are not in the configuration
3048 // Detailed description
3052 case 0: // 5432 TQ 0
3056 case 1: // 4321 TQ 0
3060 case 2: // 3210 TQ 0
3064 case 3: // 5321 TQ 1
3068 case 4: // 4210 TQ 1
3072 case 5: // 5431 TQ 1
3076 case 6: // 4320 TQ 1
3080 case 7: // 5430 TQ 2
3084 case 8: // 5210 TQ 2
3088 case 9: // 5421 TQ 3
3092 case 10: // 4310 TQ 3
3096 case 11: // 5410 TQ 4
3100 case 12: // 5420 TQ 5
3104 case 13: // 5320 TQ 5
3108 case 14: // 5310 TQ 5
3115 //____________________________________________________________________
3116 AliCluster* AliTRDtrackerV1::GetCluster(Int_t idx) const
3118 Int_t ncls = fClusters->GetEntriesFast();
3119 return idx >= 0 && idx < ncls ? (AliCluster*)fClusters->UncheckedAt(idx) : 0x0;
3122 //____________________________________________________________________
3123 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(Int_t idx) const
3125 Int_t ntrklt = fTracklets->GetEntriesFast();
3126 return idx >= 0 && idx < ntrklt ? (AliTRDseedV1*)fTracklets->UncheckedAt(idx) : 0x0;
3129 //____________________________________________________________________
3130 AliKalmanTrack* AliTRDtrackerV1::GetTrack(Int_t idx) const
3132 Int_t ntrk = fTracks->GetEntriesFast();
3133 return idx >= 0 && idx < ntrk ? (AliKalmanTrack*)fTracks->UncheckedAt(idx) : 0x0;
3136 //____________________________________________________________________
3137 Float_t AliTRDtrackerV1::CalculateReferenceX(AliTRDseedV1 *tracklets){
3139 // Calculates the reference x-position for the tilted Rieman fit defined as middle
3140 // of the stack (middle between layers 2 and 3). For the calculation all the tracklets
3141 // are taken into account
3143 // Parameters: - Array of tracklets(AliTRDseedV1)
3145 // Output: - The reference x-position(Float_t)
3147 Int_t nDistances = 0;
3148 Float_t meanDistance = 0.;
3149 Int_t startIndex = 5;
3150 for(Int_t il =5; il > 0; il--){
3151 if(tracklets[il].IsOK() && tracklets[il -1].IsOK()){
3152 Float_t xdiff = tracklets[il].GetX0() - tracklets[il -1].GetX0();
3153 meanDistance += xdiff;
3156 if(tracklets[il].IsOK()) startIndex = il;
3158 if(tracklets[0].IsOK()) startIndex = 0;
3160 // We should normally never get here
3161 Float_t xpos[2]; memset(xpos, 0, sizeof(Float_t) * 2);
3162 Int_t iok = 0, idiff = 0;
3163 // This attempt is worse and should be avoided:
3164 // check for two chambers which are OK and repeat this without taking the mean value
3165 // Strategy avoids a division by 0;
3166 for(Int_t il = 5; il >= 0; il--){
3167 if(tracklets[il].IsOK()){
3168 xpos[iok] = tracklets[il].GetX0();
3172 if(iok) idiff++; // to get the right difference;
3176 meanDistance = (xpos[0] - xpos[1])/idiff;
3179 // we have do not even have 2 layers which are OK? The we do not need to fit at all
3184 meanDistance /= nDistances;
3186 return tracklets[startIndex].GetX0() + (2.5 - startIndex) * meanDistance - 0.5 * (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
3189 //_____________________________________________________________________________
3190 Int_t AliTRDtrackerV1::Freq(Int_t n, const Int_t *inlist
3191 , Int_t *outlist, Bool_t down)
3194 // Sort eleements according occurancy
3195 // The size of output array has is 2*n
3202 Int_t *sindexS = new Int_t[n]; // Temporary array for sorting
3203 Int_t *sindexF = new Int_t[2*n];
3204 for (Int_t i = 0; i < n; i++) {
3208 TMath::Sort(n,inlist,sindexS,down);
3210 Int_t last = inlist[sindexS[0]];
3213 sindexF[0+n] = last;
3217 for (Int_t i = 1; i < n; i++) {
3218 val = inlist[sindexS[i]];
3220 sindexF[countPos]++;
3224 sindexF[countPos+n] = val;
3225 sindexF[countPos]++;
3233 // Sort according frequency
3234 TMath::Sort(countPos,sindexF,sindexS,kTRUE);
3236 for (Int_t i = 0; i < countPos; i++) {
3237 outlist[2*i ] = sindexF[sindexS[i]+n];
3238 outlist[2*i+1] = sindexF[sindexS[i]];
3249 //____________________________________________________________________
3251 //_____________________________________________________________________________
3252 Float_t AliTRDtrackerV1::GetChi2Y(AliTRDseedV1 *tracklets) const
3254 // Chi2 definition on y-direction
3257 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
3258 if(!tracklets[ipl].IsOK()) continue;
3259 Double_t distLayer = (tracklets[ipl].GetYfit(0) - tracklets[ipl].GetYref(0));// /tracklets[ipl].GetSigmaY();
3260 chi2 += distLayer * distLayer;
3265 //____________________________________________________________________
3266 void AliTRDtrackerV1::ResetSeedTB()
3268 // reset buffer for seeding time bin layers. If the time bin
3269 // layers are not allocated this function allocates them
3271 for(Int_t isl=0; isl<kNSeedPlanes; isl++){
3272 if(!fSeedTB[isl]) fSeedTB[isl] = new AliTRDchamberTimeBin();
3273 else fSeedTB[isl]->Clear();
3277 //_____________________________________________________________________________
3278 Float_t AliTRDtrackerV1::GetChi2Z(AliTRDseedV1 *tracklets) const
3280 // Calculates normalized chi2 in z-direction
3283 // chi2 = Sum ((z - zmu)/sigma)^2
3284 // Sigma for the z direction is defined as half of the padlength
3285 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
3286 if(!tracklets[ipl].IsOK()) continue;
3287 Double_t distLayer = (tracklets[ipl].GetMeanz() - tracklets[ipl].GetZref(0)); // /(tracklets[ipl].GetPadLength()/2);
3288 chi2 += distLayer * distLayer;
3293 ///////////////////////////////////////////////////////
3295 // Resources of class AliTRDLeastSquare //
3297 ///////////////////////////////////////////////////////
3299 //_____________________________________________________________________________
3300 AliTRDtrackerV1::AliTRDLeastSquare::AliTRDLeastSquare(){
3302 // Constructor of the nested class AliTRDtrackFitterLeastSquare
3304 memset(fParams, 0, sizeof(Double_t) * 2);
3305 memset(fSums, 0, sizeof(Double_t) * 5);
3306 memset(fCovarianceMatrix, 0, sizeof(Double_t) * 3);
3310 //_____________________________________________________________________________
3311 void AliTRDtrackerV1::AliTRDLeastSquare::AddPoint(Double_t *x, Double_t y, Double_t sigmaY){
3313 // Adding Point to the fitter
3315 Double_t weight = 1/(sigmaY * sigmaY);
3317 // printf("Adding point x = %f, y = %f, sigma = %f\n", xpt, y, sigmaY);
3319 fSums[1] += weight * xpt;
3320 fSums[2] += weight * y;
3321 fSums[3] += weight * xpt * y;
3322 fSums[4] += weight * xpt * xpt;
3323 fSums[5] += weight * y * y;
3326 //_____________________________________________________________________________
3327 void AliTRDtrackerV1::AliTRDLeastSquare::RemovePoint(Double_t *x, Double_t y, Double_t sigmaY){
3329 // Remove Point from the sample
3331 Double_t weight = 1/(sigmaY * sigmaY);
3334 fSums[1] -= weight * xpt;
3335 fSums[2] -= weight * y;
3336 fSums[3] -= weight * xpt * y;
3337 fSums[4] -= weight * xpt * xpt;
3338 fSums[5] -= weight * y * y;
3341 //_____________________________________________________________________________
3342 void AliTRDtrackerV1::AliTRDLeastSquare::Eval(){
3344 // Evaluation of the fit:
3345 // Calculation of the parameters
3346 // Calculation of the covariance matrix
3349 Double_t denominator = fSums[0] * fSums[4] - fSums[1] *fSums[1];
3350 if(denominator==0) return;
3352 // for(Int_t isum = 0; isum < 5; isum++)
3353 // printf("fSums[%d] = %f\n", isum, fSums[isum]);
3354 // printf("denominator = %f\n", denominator);
3355 fParams[0] = (fSums[2] * fSums[4] - fSums[1] * fSums[3])/ denominator;
3356 fParams[1] = (fSums[0] * fSums[3] - fSums[1] * fSums[2]) / denominator;
3357 // printf("fParams[0] = %f, fParams[1] = %f\n", fParams[0], fParams[1]);
3359 // Covariance matrix
3360 fCovarianceMatrix[0] = fSums[4] - fSums[1] * fSums[1] / fSums[0];
3361 fCovarianceMatrix[1] = fSums[5] - fSums[2] * fSums[2] / fSums[0];
3362 fCovarianceMatrix[2] = fSums[3] - fSums[1] * fSums[2] / fSums[0];
3365 //_____________________________________________________________________________
3366 Double_t AliTRDtrackerV1::AliTRDLeastSquare::GetFunctionValue(Double_t *xpos) const {
3368 // Returns the Function value of the fitted function at a given x-position
3370 return fParams[0] + fParams[1] * (*xpos);
3373 //_____________________________________________________________________________
3374 void AliTRDtrackerV1::AliTRDLeastSquare::GetCovarianceMatrix(Double_t *storage) const {
3376 // Copies the values of the covariance matrix into the storage
3378 memcpy(storage, fCovarianceMatrix, sizeof(Double_t) * 3);