2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
19 ///////////////////////////////////////////////////////////////////////////////
24 // Alex Bercuci <A.Bercuci@gsi.de> //
25 // Markus Fasel <M.Fasel@gsi.de> //
27 ///////////////////////////////////////////////////////////////////////////////
29 // #include <Riostream.h>
31 // #include <string.h>
34 #include <TDirectory.h>
35 #include <TLinearFitter.h>
37 #include <TClonesArray.h>
38 #include <TTreeStream.h>
41 #include "AliESDEvent.h"
42 #include "AliGeomManager.h"
43 #include "AliRieman.h"
44 #include "AliTrackPointArray.h"
46 #include "AliTRDgeometry.h"
47 #include "AliTRDpadPlane.h"
48 #include "AliTRDcalibDB.h"
49 #include "AliTRDReconstructor.h"
50 #include "AliTRDCalibraFillHisto.h"
51 #include "AliTRDrecoParam.h"
53 #include "AliTRDcluster.h"
54 #include "AliTRDseedV1.h"
55 #include "AliTRDtrackV1.h"
56 #include "AliTRDtrackerV1.h"
57 #include "AliTRDtrackerDebug.h"
58 #include "AliTRDtrackingChamber.h"
59 #include "AliTRDchamberTimeBin.h"
63 ClassImp(AliTRDtrackerV1)
66 const Float_t AliTRDtrackerV1::fgkMinClustersInTrack = 0.5; //
67 const Float_t AliTRDtrackerV1::fgkLabelFraction = 0.8; //
68 const Double_t AliTRDtrackerV1::fgkMaxChi2 = 12.0; //
69 const Double_t AliTRDtrackerV1::fgkMaxSnp = 0.95; // Maximum local sine of the azimuthal angle
70 const Double_t AliTRDtrackerV1::fgkMaxStep = 2.0; // Maximal step size in propagation
71 Double_t AliTRDtrackerV1::fgTopologicQA[kNConfigs] = {
72 0.1112, 0.1112, 0.1112, 0.0786, 0.0786,
73 0.0786, 0.0786, 0.0579, 0.0579, 0.0474,
74 0.0474, 0.0408, 0.0335, 0.0335, 0.0335
76 Int_t AliTRDtrackerV1::fgNTimeBins = 0;
77 TTreeSRedirector *AliTRDtrackerV1::fgDebugStreamer = 0x0;
78 AliRieman* AliTRDtrackerV1::fgRieman = 0x0;
79 TLinearFitter* AliTRDtrackerV1::fgTiltedRieman = 0x0;
80 TLinearFitter* AliTRDtrackerV1::fgTiltedRiemanConstrained = 0x0;
82 //____________________________________________________________________
83 AliTRDtrackerV1::AliTRDtrackerV1(AliTRDReconstructor *rec)
86 ,fGeom(new AliTRDgeometry())
93 // Default constructor.
95 AliTRDcalibDB *trd = 0x0;
96 if (!(trd = AliTRDcalibDB::Instance())) {
97 AliFatal("Could not get calibration object");
100 if(!fgNTimeBins) fgNTimeBins = trd->GetNumberOfTimeBins();
102 for (Int_t isector = 0; isector < AliTRDgeometry::kNsector; isector++) new(&fTrSec[isector]) AliTRDtrackingSector(fGeom, isector);
104 for(Int_t isl =0; isl<kNSeedPlanes; isl++) fSeedTB[isl] = 0x0;
106 // Initialize debug stream
107 if(rec) SetReconstructor(rec);
110 //____________________________________________________________________
111 AliTRDtrackerV1::~AliTRDtrackerV1()
117 if(fgDebugStreamer) delete fgDebugStreamer;
118 if(fgRieman) delete fgRieman;
119 if(fgTiltedRieman) delete fgTiltedRieman;
120 if(fgTiltedRiemanConstrained) delete fgTiltedRiemanConstrained;
121 for(Int_t isl =0; isl<kNSeedPlanes; isl++) if(fSeedTB[isl]) delete fSeedTB[isl];
122 if(fTracks) {fTracks->Delete(); delete fTracks;}
123 if(fTracklets) {fTracklets->Delete(); delete fTracklets;}
124 if(fClusters) {fClusters->Delete(); delete fClusters;}
125 if(fGeom) delete fGeom;
128 //____________________________________________________________________
129 Int_t AliTRDtrackerV1::Clusters2Tracks(AliESDEvent *esd)
132 // Steering stand alone tracking for full TRD detector
135 // esd : The ESD event. On output it contains
136 // the ESD tracks found in TRD.
139 // Number of tracks found in the TRD detector.
141 // Detailed description
142 // 1. Launch individual SM trackers.
143 // See AliTRDtrackerV1::Clusters2TracksSM() for details.
146 if(!fReconstructor->GetRecoParam() ){
147 AliError("Reconstruction configuration not initialized. Call first AliTRDReconstructor::SetRecoParam().");
151 //AliInfo("Start Track Finder ...");
153 for(int ism=0; ism<AliTRDgeometry::kNsector; ism++){
154 // for(int ism=1; ism<2; ism++){
155 //AliInfo(Form("Processing supermodule %i ...", ism));
156 ntracks += Clusters2TracksSM(ism, esd);
158 AliInfo(Form("Number of found tracks : %d", ntracks));
163 //_____________________________________________________________________________
164 Bool_t AliTRDtrackerV1::GetTrackPoint(Int_t index, AliTrackPoint &p) const
166 //AliInfo(Form("Asking for tracklet %d", index));
168 AliTRDseedV1 *tracklet = GetTracklet(index);
169 if (!tracklet) return kFALSE;
171 // get detector for this tracklet
172 AliTRDcluster *cl = 0x0;
173 Int_t ic = 0; do; while(!(cl = tracklet->GetClusters(ic++)));
174 Int_t idet = cl->GetDetector();
177 local[0] = tracklet->GetX0();
178 local[1] = tracklet->GetYfit(0);
179 local[2] = tracklet->GetZfit(0);
181 fGeom->RotateBack(idet, local, global);
182 p.SetXYZ(global[0],global[1],global[2]);
186 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
187 switch (fGeom->GetLayer(idet)) {
189 iLayer = AliGeomManager::kTRD1;
192 iLayer = AliGeomManager::kTRD2;
195 iLayer = AliGeomManager::kTRD3;
198 iLayer = AliGeomManager::kTRD4;
201 iLayer = AliGeomManager::kTRD5;
204 iLayer = AliGeomManager::kTRD6;
207 Int_t modId = fGeom->GetSector(idet) * fGeom->Nstack() + fGeom->GetStack(idet);
208 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer, modId);
209 p.SetVolumeID(volid);
214 //____________________________________________________________________
215 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitter()
217 if(!fgTiltedRieman) fgTiltedRieman = new TLinearFitter(4, "hyp4");
218 return fgTiltedRieman;
221 //____________________________________________________________________
222 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitterConstraint()
224 if(!fgTiltedRiemanConstrained) fgTiltedRiemanConstrained = new TLinearFitter(2, "hyp2");
225 return fgTiltedRiemanConstrained;
228 //____________________________________________________________________
229 AliRieman* AliTRDtrackerV1::GetRiemanFitter()
231 if(!fgRieman) fgRieman = new AliRieman(AliTRDtrackingChamber::kNTimeBins * AliTRDgeometry::kNlayer);
235 //_____________________________________________________________________________
236 Int_t AliTRDtrackerV1::PropagateBack(AliESDEvent *event)
239 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
240 // backpropagated by the TPC tracker. Each seed is first propagated
241 // to the TRD, and then its prolongation is searched in the TRD.
242 // If sufficiently long continuation of the track is found in the TRD
243 // the track is updated, otherwise it's stored as originaly defined
244 // by the TPC tracker.
247 // Calibration monitor
248 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
249 if (!calibra) AliInfo("Could not get Calibra instance\n");
251 Int_t found = 0; // number of tracks found
252 Float_t foundMin = 20.0;
254 Float_t *quality = 0x0;
256 Int_t nSeed = event->GetNumberOfTracks();
258 quality = new Float_t[nSeed];
259 index = new Int_t[nSeed];
260 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
261 AliESDtrack *seed = event->GetTrack(iSeed);
262 Double_t covariance[15];
263 seed->GetExternalCovariance(covariance);
264 quality[iSeed] = covariance[0] + covariance[2];
266 // Sort tracks according to covariance of local Y and Z
267 TMath::Sort(nSeed,quality,index,kFALSE);
270 // Backpropagate all seeds
273 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
275 // Get the seeds in sorted sequence
276 AliESDtrack *seed = event->GetTrack(index[iSeed]);
278 // Check the seed status
279 ULong_t status = seed->GetStatus();
280 if ((status & AliESDtrack::kTPCout) == 0) continue;
281 if ((status & AliESDtrack::kTRDout) != 0) continue;
283 // Do the back prolongation
284 new(&track) AliTRDtrackV1(*seed);
286 //Int_t lbl = seed->GetLabel();
287 //track.SetSeedLabel(lbl);
288 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup); // Make backup
289 Float_t p4 = track.GetC();
290 expectedClr = FollowBackProlongation(track);
292 if (expectedClr<0) continue; // Back prolongation failed
296 // computes PID for track
298 // update calibration references using this track
299 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(&track);
300 // save calibration object
301 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
302 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
304 track.UpdateESDtrack(seed);
306 // Add TRD track to ESDfriendTrack
307 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0 /*&& quality TODO*/){
308 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(track);
309 calibTrack->SetOwner();
310 seed->AddCalibObject(calibTrack);
315 if ((TMath::Abs(track.GetC() - p4) / TMath::Abs(p4) < 0.2) ||(track.Pt() > 0.8)) {
317 // Make backup for back propagation
319 Int_t foundClr = track.GetNumberOfClusters();
320 if (foundClr >= foundMin) {
321 //AliInfo(Form("Making backup track ncls [%d]...", foundClr));
323 //track.CookdEdxTimBin(seed->GetID());
324 track.CookLabel(1. - fgkLabelFraction);
325 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)) UseClusters(&track);
332 Bool_t isGold = kFALSE;
335 if (track.GetChi2() / track.GetNumberOfClusters() < 5) {
336 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
342 if ((!isGold) && (track.GetNCross() == 0) && (track.GetChi2() / track.GetNumberOfClusters() < 7)) {
343 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
344 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
349 if ((!isGold) && (track.GetBackupTrack())) {
350 if ((track.GetBackupTrack()->GetNumberOfClusters() > foundMin) && ((track.GetBackupTrack()->GetChi2()/(track.GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
351 seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
356 //if ((track->StatusForTOF() > 0) && (track->GetNCross() == 0) && (Float_t(track->GetNumberOfClusters()) / Float_t(track->GetNExpected()) > 0.4)) {
357 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
362 // Propagation to the TOF (I.Belikov)
363 if (track.IsStopped() == kFALSE) {
364 Double_t xtof = 371.0;
365 Double_t xTOF0 = 370.0;
367 Double_t c2 = track.GetSnp() + track.GetC() * (xtof - track.GetX());
368 if (TMath::Abs(c2) >= 0.99) continue;
370 if (!PropagateToX(track, xTOF0, fgkMaxStep)) continue;
372 // Energy losses taken to the account - check one more time
373 c2 = track.GetSnp() + track.GetC() * (xtof - track.GetX());
374 if (TMath::Abs(c2) >= 0.99) continue;
376 //if (!PropagateToX(*track,xTOF0,fgkMaxStep)) {
377 // fHBackfit->Fill(7);
382 Double_t ymax = xtof * TMath::Tan(0.5 * AliTRDgeometry::GetAlpha());
384 track.GetYAt(xtof,GetBz(),y);
386 if (!track.Rotate( AliTRDgeometry::GetAlpha())) continue;
387 }else if (y < -ymax) {
388 if (!track.Rotate(-AliTRDgeometry::GetAlpha())) continue;
391 if (track.PropagateTo(xtof)) {
392 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
393 track.UpdateESDtrack(seed);
396 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
397 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
399 track.UpdateESDtrack(seed);
403 seed->SetTRDQuality(track.StatusForTOF());
404 seed->SetTRDBudget(track.GetBudget(0));
406 if(index) delete [] index;
407 if(quality) delete [] quality;
410 AliInfo(Form("Number of seeds: %d", nSeed));
411 AliInfo(Form("Number of back propagated TRD tracks: %d", found));
413 // run stand alone tracking
414 if (fReconstructor->IsSeeding()) Clusters2Tracks(event);
420 //____________________________________________________________________
421 Int_t AliTRDtrackerV1::RefitInward(AliESDEvent *event)
424 // Refits tracks within the TRD. The ESD event is expected to contain seeds
425 // at the outer part of the TRD.
426 // The tracks are propagated to the innermost time bin
427 // of the TRD and the ESD event is updated
428 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
431 Int_t nseed = 0; // contor for loaded seeds
432 Int_t found = 0; // contor for updated TRD tracks
436 for (Int_t itrack = 0; itrack < event->GetNumberOfTracks(); itrack++) {
437 AliESDtrack *seed = event->GetTrack(itrack);
438 new(&track) AliTRDtrackV1(*seed);
440 if (track.GetX() < 270.0) {
441 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
445 ULong_t status = seed->GetStatus();
446 if((status & AliESDtrack::kTRDout) == 0) continue;
447 if((status & AliESDtrack::kTRDin) != 0) continue;
450 track.ResetCovariance(50.0);
452 // do the propagation and processing
453 Bool_t kUPDATE = kFALSE;
454 Double_t xTPC = 250.0;
455 if(FollowProlongation(track)){
457 if (PropagateToX(track, xTPC, fgkMaxStep)) { // -with update
458 seed->UpdateTrackParams(&track, AliESDtrack::kTRDrefit);
464 // Prolongate to TPC without update
466 AliTRDtrackV1 tt(*seed);
467 if (PropagateToX(tt, xTPC, fgkMaxStep)) seed->UpdateTrackParams(&tt, AliESDtrack::kTRDrefit);
470 AliInfo(Form("Number of loaded seeds: %d",nseed));
471 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
476 //____________________________________________________________________
477 Int_t AliTRDtrackerV1::FollowProlongation(AliTRDtrackV1 &t)
479 // Extrapolates the TRD track in the TPC direction.
482 // t : the TRD track which has to be extrapolated
485 // number of clusters attached to the track
487 // Detailed description
489 // Starting from current radial position of track <t> this function
490 // extrapolates the track through the 6 TRD layers. The following steps
491 // are being performed for each plane:
493 // a. get plane limits in the local x direction
494 // b. check crossing sectors
495 // c. check track inclination
496 // 2. search tracklet in the tracker list (see GetTracklet() for details)
497 // 3. evaluate material budget using the geo manager
498 // 4. propagate and update track using the tracklet information.
503 Int_t nClustersExpected = 0;
504 Int_t lastplane = 5; //GetLastPlane(&t);
505 for (Int_t iplane = lastplane; iplane >= 0; iplane--) {
507 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
508 if(!tracklet) continue;
509 if(!tracklet->IsOK()) AliWarning("tracklet not OK");
511 Double_t x = tracklet->GetX0();
512 // reject tracklets which are not considered for inward refit
513 if(x > t.GetX()+fgkMaxStep) continue;
515 // append tracklet to track
516 t.SetTracklet(tracklet, index);
518 if (x < (t.GetX()-fgkMaxStep) && !PropagateToX(t, x+fgkMaxStep, fgkMaxStep)) break;
519 if (!AdjustSector(&t)) break;
521 // Start global position
525 // End global position
526 Double_t alpha = t.GetAlpha(), y, z;
527 if (!t.GetProlongation(x,y,z)) break;
529 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
530 xyz1[1] = x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
533 // Get material budget
535 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
536 Double_t xrho= param[0]*param[4];
537 Double_t xx0 = param[1]; // Get mean propagation parameters
539 // Propagate and update
540 t.PropagateTo(x, xx0, xrho);
541 if (!AdjustSector(&t)) break;
543 Double_t maxChi2 = t.GetPredictedChi2(tracklet);
544 if (maxChi2 < 1e+10 && t.Update(tracklet, maxChi2)){
545 nClustersExpected += tracklet->GetN();
549 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
551 for(int iplane=0; iplane<6; iplane++){
552 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
553 if(!tracklet) continue;
554 t.SetTracklet(tracklet, index);
557 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
558 TTreeSRedirector &cstreamer = *fgDebugStreamer;
559 cstreamer << "FollowProlongation"
560 << "EventNumber=" << eventNumber
561 << "ncl=" << nClustersExpected
566 return nClustersExpected;
570 //_____________________________________________________________________________
571 Int_t AliTRDtrackerV1::FollowBackProlongation(AliTRDtrackV1 &t)
573 // Extrapolates the TRD track in the TOF direction.
576 // t : the TRD track which has to be extrapolated
579 // number of clusters attached to the track
581 // Detailed description
583 // Starting from current radial position of track <t> this function
584 // extrapolates the track through the 6 TRD layers. The following steps
585 // are being performed for each plane:
587 // a. get plane limits in the local x direction
588 // b. check crossing sectors
589 // c. check track inclination
590 // 2. build tracklet (see AliTRDseed::AttachClusters() for details)
591 // 3. evaluate material budget using the geo manager
592 // 4. propagate and update track using the tracklet information.
597 Int_t nClustersExpected = 0;
598 Double_t clength = AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick();
599 AliTRDtrackingChamber *chamber = 0x0;
601 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
602 // in case of stand alone tracking we store all the pointers to the tracklets in a temporary array
603 AliTRDseedV1 *tracklets[kNPlanes];
604 memset(tracklets, 0, sizeof(AliTRDseedV1 *) * kNPlanes);
605 for(Int_t ip = 0; ip < kNPlanes; ip++){
606 tracklets[ip] = t.GetTracklet(ip);
610 // Loop through the TRD layers
611 for (Int_t ilayer = 0; ilayer < AliTRDgeometry::Nlayer(); ilayer++) {
612 // BUILD TRACKLET IF NOT ALREADY BUILT
613 Double_t x = 0., y, z, alpha;
614 ptrTracklet = tracklets[ilayer];
616 ptrTracklet = new(&tracklet) AliTRDseedV1(ilayer);
617 ptrTracklet->SetReconstructor(fReconstructor);
618 alpha = t.GetAlpha();
619 Int_t sector = Int_t(alpha/AliTRDgeometry::GetAlpha() + (alpha>0. ? 0 : AliTRDgeometry::kNsector));
621 if(!fTrSec[sector].GetNChambers()) continue;
623 if((x = fTrSec[sector].GetX(ilayer)) < 1.) continue;
625 if (!t.GetProlongation(x, y, z)) return -nClustersExpected;
626 Int_t stack = fGeom->GetStack(z, ilayer);
627 Int_t nCandidates = stack >= 0 ? 1 : 2;
628 z -= stack >= 0 ? 0. : 4.;
630 for(int icham=0; icham<nCandidates; icham++, z+=8){
631 if((stack = fGeom->GetStack(z, ilayer)) < 0) continue;
633 if(!(chamber = fTrSec[sector].GetChamber(stack, ilayer))) continue;
635 if(chamber->GetNClusters() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
639 AliTRDpadPlane *pp = fGeom->GetPadPlane(ilayer, stack);
640 tracklet.SetTilt(TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle()));
641 tracklet.SetPadLength(pp->GetLengthIPad());
642 tracklet.SetPlane(ilayer);
644 if(!tracklet.Init(&t)){
646 return nClustersExpected;
648 if(!tracklet.AttachClustersIter(chamber, 1000.)) continue;
651 if(tracklet.GetN() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
656 if(!ptrTracklet->IsOK()){
657 if(x < 1.) continue; //temporary
658 if(!PropagateToX(t, x-fgkMaxStep, fgkMaxStep)) return -nClustersExpected;
659 if(!AdjustSector(&t)) return -nClustersExpected;
660 if(TMath::Abs(t.GetSnp()) > fgkMaxSnp) return -nClustersExpected;
664 // Propagate closer to the current chamber if neccessary
666 if (x > (fgkMaxStep + t.GetX()) && !PropagateToX(t, x-fgkMaxStep, fgkMaxStep)) return -nClustersExpected;
667 if (!AdjustSector(&t)) return -nClustersExpected;
668 if (TMath::Abs(t.GetSnp()) > fgkMaxSnp) return -nClustersExpected;
670 // load tracklet to the tracker and the track
671 ptrTracklet = SetTracklet(ptrTracklet);
672 t.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
675 // Calculate the mean material budget along the path inside the chamber
676 //Calculate global entry and exit positions of the track in chamber (only track prolongation)
677 Double_t xyz0[3]; // entry point
679 alpha = t.GetAlpha();
680 x = ptrTracklet->GetX0();
681 if (!t.GetProlongation(x, y, z)) return -nClustersExpected;
682 Double_t xyz1[3]; // exit point
683 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
684 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
687 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
688 // The mean propagation parameters
689 Double_t xrho = param[0]*param[4]; // density*length
690 Double_t xx0 = param[1]; // radiation length
692 // Propagate and update track
693 if (!t.PropagateTo(x, xx0, xrho)) return -nClustersExpected;
694 if (!AdjustSector(&t)) return -nClustersExpected;
695 Double_t maxChi2 = t.GetPredictedChi2(ptrTracklet);
696 if (!t.Update(ptrTracklet, maxChi2)) return -nClustersExpected;
698 nClustersExpected += ptrTracklet->GetN();
699 //t.SetTracklet(&tracklet, index);
701 // Reset material budget if 2 consecutive gold
702 if(ilayer>0 && t.GetTracklet(ilayer-1) && ptrTracklet->GetN() + t.GetTracklet(ilayer-1)->GetN() > 20) t.SetBudget(2, 0.);
704 // Make backup of the track until is gold
705 // TO DO update quality check of the track.
706 // consider comparison with fTimeBinsRange
707 Float_t ratio0 = ptrTracklet->GetN() / Float_t(fgNTimeBins);
708 //Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);
709 //printf("tracklet.GetChi2() %f [< 18.0]\n", tracklet.GetChi2());
710 //printf("ratio0 %f [> 0.8]\n", ratio0);
711 //printf("ratio1 %f [> 0.6]\n", ratio1);
712 //printf("ratio0+ratio1 %f [> 1.5]\n", ratio0+ratio1);
713 //printf("t.GetNCross() %d [== 0]\n", t.GetNCross());
714 //printf("TMath::Abs(t.GetSnp()) %f [< 0.85]\n", TMath::Abs(t.GetSnp()));
715 //printf("t.GetNumberOfClusters() %d [> 20]\n", t.GetNumberOfClusters());
717 if (//(tracklet.GetChi2() < 18.0) && TO DO check with FindClusters and move it to AliTRDseed::Update
720 //(ratio0+ratio1 > 1.5) &&
721 (t.GetNCross() == 0) &&
722 (TMath::Abs(t.GetSnp()) < 0.85) &&
723 (t.GetNumberOfClusters() > 20)) t.MakeBackupTrack();
727 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
728 TTreeSRedirector &cstreamer = *fgDebugStreamer;
729 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
730 //AliTRDtrackV1 *debugTrack = new AliTRDtrackV1(t);
731 //debugTrack->SetOwner();
732 cstreamer << "FollowBackProlongation"
733 << "EventNumber=" << eventNumber
734 << "ncl=" << nClustersExpected
735 //<< "track.=" << debugTrack
739 return nClustersExpected;
742 //_________________________________________________________________________
743 Float_t AliTRDtrackerV1::FitRieman(AliTRDseedV1 *tracklets, Double_t *chi2, Int_t *planes){
745 // Fits a Riemann-circle to the given points without tilting pad correction.
746 // The fit is performed using an instance of the class AliRieman (equations
747 // and transformations see documentation of this class)
748 // Afterwards all the tracklets are Updated
750 // Parameters: - Array of tracklets (AliTRDseedV1)
751 // - Storage for the chi2 values (beginning with direction z)
752 // - Seeding configuration
753 // Output: - The curvature
755 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
757 Int_t allplanes[] = {0, 1, 2, 3, 4, 5};
758 Int_t *ppl = &allplanes[0];
764 for(Int_t il = 0; il < maxLayers; il++){
765 if(!tracklets[ppl[il]].IsOK()) continue;
766 fitter->AddPoint(tracklets[ppl[il]].GetX0(), tracklets[ppl[il]].GetYfitR(0), tracklets[ppl[il]].GetZProb(),1,10);
769 // Set the reference position of the fit and calculate the chi2 values
770 memset(chi2, 0, sizeof(Double_t) * 2);
771 for(Int_t il = 0; il < maxLayers; il++){
772 // Reference positions
773 tracklets[ppl[il]].Init(fitter);
776 if((!tracklets[ppl[il]].IsOK()) && (!planes)) continue;
777 chi2[0] += tracklets[ppl[il]].GetChi2Y();
778 chi2[1] += tracklets[ppl[il]].GetChi2Z();
780 return fitter->GetC();
783 //_________________________________________________________________________
784 void AliTRDtrackerV1::FitRieman(AliTRDcluster **seedcl, Double_t chi2[2])
787 // Performs a Riemann helix fit using the seedclusters as spacepoints
788 // Afterwards the chi2 values are calculated and the seeds are updated
790 // Parameters: - The four seedclusters
791 // - The tracklet array (AliTRDseedV1)
792 // - The seeding configuration
797 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
799 for(Int_t i = 0; i < 4; i++)
800 fitter->AddPoint(seedcl[i]->GetX(), seedcl[i]->GetY(), seedcl[i]->GetZ(), 1, 10);
804 // Update the seed and calculated the chi2 value
805 chi2[0] = 0; chi2[1] = 0;
806 for(Int_t ipl = 0; ipl < kNSeedPlanes; ipl++){
808 chi2[0] += (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX()));
809 chi2[1] += (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX()));
814 //_________________________________________________________________________
815 Float_t AliTRDtrackerV1::FitTiltedRiemanConstraint(AliTRDseedV1 *tracklets, Double_t zVertex)
818 // Fits a helix to the clusters. Pad tilting is considered. As constraint it is
819 // assumed that the vertex position is set to 0.
820 // This method is very usefull for high-pt particles
821 // Basis for the fit: (x - x0)^2 + (y - y0)^2 - R^2 = 0
822 // x0, y0: Center of the circle
823 // Measured y-position: ymeas = y - tan(phiT)(zc - zt)
824 // zc: center of the pad row
825 // Equation which has to be fitted (after transformation):
826 // a + b * u + e * v + 2*(ymeas + tan(phiT)(z - zVertex))*t = 0
830 // v = 2 * x * tan(phiT) * t
831 // Parameters in the equation:
832 // a = -1/y0, b = x0/y0, e = dz/dx
834 // The Curvature is calculated by the following equation:
835 // - curv = a/Sqrt(b^2 + 1) = 1/R
836 // Parameters: - the 6 tracklets
837 // - the Vertex constraint
838 // Output: - the Chi2 value of the track
843 TLinearFitter *fitter = GetTiltedRiemanFitterConstraint();
844 fitter->StoreData(kTRUE);
845 fitter->ClearPoints();
846 AliTRDcluster *cl = 0x0;
848 Float_t x, y, z, w, t, error, tilt;
851 for(Int_t ilr = 0; ilr < AliTRDgeometry::kNlayer; ilr++){
852 if(!tracklets[ilr].IsOK()) continue;
853 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
854 if(!tracklets[ilr].IsUsable(itb)) continue;
855 cl = tracklets[ilr].GetClusters(itb);
859 tilt = tracklets[ilr].GetTilt();
861 t = 1./(x * x + y * y);
863 uvt[1] = 2. * x * t * tilt ;
864 w = 2. * (y + tilt * (z - zVertex)) * t;
865 error = 2. * 0.2 * t;
866 fitter->AddPoint(uvt, w, error);
872 // Calculate curvature
873 Double_t a = fitter->GetParameter(0);
874 Double_t b = fitter->GetParameter(1);
875 Double_t curvature = a/TMath::Sqrt(b*b + 1);
877 Float_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
878 for(Int_t ip = 0; ip < AliTRDtrackerV1::kNPlanes; ip++)
879 tracklets[ip].SetCC(curvature);
881 /* if(fReconstructor->GetStreamLevel() >= 5){
882 //Linear Model on z-direction
883 Double_t xref = CalculateReferenceX(tracklets); // Relative to the middle of the stack
884 Double_t slope = fitter->GetParameter(2);
885 Double_t zref = slope * xref;
886 Float_t chi2Z = CalculateChi2Z(tracklets, zref, slope, xref);
887 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
888 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
889 TTreeSRedirector &treeStreamer = *fgDebugStreamer;
890 treeStreamer << "FitTiltedRiemanConstraint"
891 << "EventNumber=" << eventNumber
892 << "CandidateNumber=" << candidateNumber
893 << "Curvature=" << curvature
894 << "Chi2Track=" << chi2track
902 //_________________________________________________________________________
903 Float_t AliTRDtrackerV1::FitTiltedRieman(AliTRDseedV1 *tracklets, Bool_t sigError)
906 // Performs a Riemann fit taking tilting pad correction into account
907 // The equation of a Riemann circle, where the y position is substituted by the
908 // measured y-position taking pad tilting into account, has to be transformed
909 // into a 4-dimensional hyperplane equation
910 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
911 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
912 // zc: center of the pad row
913 // zt: z-position of the track
914 // The z-position of the track is assumed to be linear dependent on the x-position
915 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
916 // Transformation: u = 2 * x * t
917 // v = 2 * tan(phiT) * t
918 // w = 2 * tan(phiT) * (x - xref) * t
919 // t = 1 / (x^2 + ymeas^2)
920 // Parameters: a = -1/y0
922 // c = (R^2 -x0^2 - y0^2)/y0
925 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
926 // results from the simple riemann fit. Afterwards the fit is redone.
927 // The curvature is calculated according to the formula:
928 // curv = a/(1 + b^2 + c*a) = 1/R
930 // Paramters: - Array of tracklets (connected to the track candidate)
931 // - Flag selecting the error definition
932 // Output: - Chi2 values of the track (in Parameter list)
934 TLinearFitter *fitter = GetTiltedRiemanFitter();
935 fitter->StoreData(kTRUE);
936 fitter->ClearPoints();
937 AliTRDLeastSquare zfitter;
938 AliTRDcluster *cl = 0x0;
940 Double_t xref = CalculateReferenceX(tracklets);
941 Double_t x, y, z, t, tilt, dx, w, we;
944 // Containers for Least-square fitter
945 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
946 if(!tracklets[ipl].IsOK()) continue;
947 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
948 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
949 if (!tracklets[ipl].IsUsable(itb)) continue;
953 tilt = tracklets[ipl].GetTilt();
959 uvt[2] = 2. * tilt * t;
960 uvt[3] = 2. * tilt * dx * t;
961 w = 2. * (y + tilt*z) * t;
962 // error definition changes for the different calls
964 we *= sigError ? tracklets[ipl].GetSigmaY() : 0.2;
965 fitter->AddPoint(uvt, w, we);
966 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
973 Double_t offset = fitter->GetParameter(3);
974 Double_t slope = fitter->GetParameter(4);
976 // Linear fitter - not possible to make boundaries
977 // Do not accept non possible z and dzdx combinations
978 Bool_t acceptablez = kTRUE;
980 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
981 if(!tracklets[iLayer].IsOK()) continue;
982 zref = offset + slope * (tracklets[iLayer].GetX0() - xref);
983 if (TMath::Abs(tracklets[iLayer].GetZProb() - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
984 acceptablez = kFALSE;
987 Double_t dzmf = zfitter.GetFunctionParameter(1);
988 Double_t zmf = zfitter.GetFunctionValue(&xref);
989 fgTiltedRieman->FixParameter(3, zmf);
990 fgTiltedRieman->FixParameter(4, dzmf);
992 fitter->ReleaseParameter(3);
993 fitter->ReleaseParameter(4);
994 offset = fitter->GetParameter(3);
995 slope = fitter->GetParameter(4);
998 // Calculate Curvarture
999 Double_t a = fitter->GetParameter(0);
1000 Double_t b = fitter->GetParameter(1);
1001 Double_t c = fitter->GetParameter(2);
1002 Double_t curvature = 1.0 + b*b - c*a;
1003 if (curvature > 0.0)
1004 curvature = a / TMath::Sqrt(curvature);
1006 Double_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
1008 // Update the tracklets
1010 for(Int_t iLayer = 0; iLayer < AliTRDtrackerV1::kNPlanes; iLayer++) {
1012 x = tracklets[iLayer].GetX0();
1018 // y: R^2 = (x - x0)^2 + (y - y0)^2
1019 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1020 // R = Sqrt() = 1/Curvature
1021 // => y = y0 +/- Sqrt(1/Curvature^2 - (x - x0)^2)
1022 Double_t res = (x * a + b); // = (x - x0)/y0
1024 res = 1.0 - c * a + b * b - res; // = (R^2 - (x - x0)^2)/y0^2
1026 res = TMath::Sqrt(res);
1027 y = (1.0 - res) / a;
1030 // dy: R^2 = (x - x0)^2 + (y - y0)^2
1031 // => y = +/- Sqrt(R^2 - (x - x0)^2) + y0
1032 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1033 // Curvature: cr = 1/R = a/Sqrt(1 + b^2 - c*a)
1034 // => dy/dx = (x - x0)/(1/(cr^2) - (x - x0)^2)
1035 Double_t x0 = -b / a;
1036 if (-c * a + b * b + 1 > 0) {
1037 if (1.0/(curvature * curvature) - (x - x0) * (x - x0) > 0.0) {
1038 Double_t yderiv = (x - x0) / TMath::Sqrt(1.0/(curvature * curvature) - (x - x0) * (x - x0));
1039 if (a < 0) yderiv *= -1.0;
1043 z = offset + slope * (x - xref);
1045 tracklets[iLayer].SetYref(0, y);
1046 tracklets[iLayer].SetYref(1, dy);
1047 tracklets[iLayer].SetZref(0, z);
1048 tracklets[iLayer].SetZref(1, dz);
1049 tracklets[iLayer].SetC(curvature);
1050 tracklets[iLayer].SetChi2(chi2track);
1053 /* if(fReconstructor->GetStreamLevel() >=5){
1054 TTreeSRedirector &cstreamer = *fgDebugStreamer;
1055 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1056 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1057 Double_t chi2z = CalculateChi2Z(tracklets, offset, slope, xref);
1058 cstreamer << "FitTiltedRieman0"
1059 << "EventNumber=" << eventNumber
1060 << "CandidateNumber=" << candidateNumber
1062 << "Chi2Z=" << chi2z
1069 //____________________________________________________________________
1070 Double_t AliTRDtrackerV1::FitLine(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t err, Int_t np, AliTrackPoint *points)
1072 AliTRDLeastSquare yfitter, zfitter;
1073 AliTRDcluster *cl = 0x0;
1075 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1077 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1078 if(!(tracklet = track->GetTracklet(ipl))) continue;
1079 if(!tracklet->IsOK()) continue;
1080 new(&work[ipl]) AliTRDseedV1(*tracklet);
1082 tracklets = &work[0];
1085 Double_t xref = CalculateReferenceX(tracklets);
1086 Double_t x, y, z, dx, ye, yr, tilt;
1087 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1088 if(!tracklets[ipl].IsOK()) continue;
1089 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1090 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1091 if (!tracklets[ipl].IsUsable(itb)) continue;
1095 zfitter.AddPoint(&dx, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1099 Double_t z0 = zfitter.GetFunctionParameter(0);
1100 Double_t dzdx = zfitter.GetFunctionParameter(1);
1101 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1102 if(!tracklets[ipl].IsOK()) continue;
1103 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1104 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1105 if (!tracklets[ipl].IsUsable(itb)) continue;
1109 tilt = tracklets[ipl].GetTilt();
1111 yr = y + tilt*(z - z0 - dzdx*dx);
1112 // error definition changes for the different calls
1113 ye = tilt*TMath::Sqrt(cl->GetSigmaZ2());
1114 ye += err ? tracklets[ipl].GetSigmaY() : 0.2;
1115 yfitter.AddPoint(&dx, yr, ye);
1119 Double_t y0 = yfitter.GetFunctionParameter(0);
1120 Double_t dydx = yfitter.GetFunctionParameter(1);
1121 Double_t chi2 = 0.;//yfitter.GetChisquare()/Double_t(nPoints);
1123 //update track points array
1126 for(int ip=0; ip<np; ip++){
1127 points[ip].GetXYZ(xyz);
1128 xyz[1] = y0 + dydx * (xyz[0] - xref);
1129 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1130 points[ip].SetXYZ(xyz);
1137 //_________________________________________________________________________
1138 Double_t AliTRDtrackerV1::FitRiemanTilt(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t sigError, Int_t np, AliTrackPoint *points)
1141 // Performs a Riemann fit taking tilting pad correction into account
1142 // The equation of a Riemann circle, where the y position is substituted by the
1143 // measured y-position taking pad tilting into account, has to be transformed
1144 // into a 4-dimensional hyperplane equation
1145 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
1146 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
1147 // zc: center of the pad row
1148 // zt: z-position of the track
1149 // The z-position of the track is assumed to be linear dependent on the x-position
1150 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
1151 // Transformation: u = 2 * x * t
1152 // v = 2 * tan(phiT) * t
1153 // w = 2 * tan(phiT) * (x - xref) * t
1154 // t = 1 / (x^2 + ymeas^2)
1155 // Parameters: a = -1/y0
1157 // c = (R^2 -x0^2 - y0^2)/y0
1160 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
1161 // results from the simple riemann fit. Afterwards the fit is redone.
1162 // The curvature is calculated according to the formula:
1163 // curv = a/(1 + b^2 + c*a) = 1/R
1165 // Paramters: - Array of tracklets (connected to the track candidate)
1166 // - Flag selecting the error definition
1167 // Output: - Chi2 values of the track (in Parameter list)
1169 TLinearFitter *fitter = GetTiltedRiemanFitter();
1170 fitter->StoreData(kTRUE);
1171 fitter->ClearPoints();
1172 AliTRDLeastSquare zfitter;
1173 AliTRDcluster *cl = 0x0;
1175 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1177 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1178 if(!(tracklet = track->GetTracklet(ipl))) continue;
1179 if(!tracklet->IsOK()) continue;
1180 new(&work[ipl]) AliTRDseedV1(*tracklet);
1182 tracklets = &work[0];
1185 Double_t xref = CalculateReferenceX(tracklets);
1186 Double_t x, y, z, t, tilt, dx, w, we;
1189 // Containers for Least-square fitter
1190 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1191 if(!tracklets[ipl].IsOK()) continue;
1192 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1193 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1194 if (!tracklets[ipl].IsUsable(itb)) continue;
1198 tilt = tracklets[ipl].GetTilt();
1202 uvt[0] = 2. * x * t;
1204 uvt[2] = 2. * tilt * t;
1205 uvt[3] = 2. * tilt * dx * t;
1206 w = 2. * (y + tilt*z) * t;
1207 // error definition changes for the different calls
1209 we *= sigError ? tracklets[ipl].GetSigmaY() : 0.2;
1210 fitter->AddPoint(uvt, w, we);
1211 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1215 if(fitter->Eval()) return 1.E10;
1217 Double_t z0 = fitter->GetParameter(3);
1218 Double_t dzdx = fitter->GetParameter(4);
1221 // Linear fitter - not possible to make boundaries
1222 // Do not accept non possible z and dzdx combinations
1223 Bool_t accept = kTRUE;
1224 Double_t zref = 0.0;
1225 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
1226 if(!tracklets[iLayer].IsOK()) continue;
1227 zref = z0 + dzdx * (tracklets[iLayer].GetX0() - xref);
1228 if (TMath::Abs(tracklets[iLayer].GetZProb() - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
1233 Double_t dzmf = zfitter.GetFunctionParameter(1);
1234 Double_t zmf = zfitter.GetFunctionValue(&xref);
1235 fitter->FixParameter(3, zmf);
1236 fitter->FixParameter(4, dzmf);
1238 fitter->ReleaseParameter(3);
1239 fitter->ReleaseParameter(4);
1240 z0 = fitter->GetParameter(3); // = zmf ?
1241 dzdx = fitter->GetParameter(4); // = dzmf ?
1244 // Calculate Curvature
1245 Double_t a = fitter->GetParameter(0);
1246 Double_t b = fitter->GetParameter(1);
1247 Double_t c = fitter->GetParameter(2);
1248 Double_t y0 = 1. / a;
1249 Double_t x0 = -b * y0;
1250 Double_t R = TMath::Sqrt(y0*y0 + x0*x0 - c*y0);
1251 Double_t C = 1.0 + b*b - c*a;
1252 if (C > 0.0) C = a / TMath::Sqrt(C);
1254 // Calculate chi2 of the fit
1255 Double_t chi2 = fitter->GetChisquare()/Double_t(nPoints);
1257 // Update the tracklets
1259 for(Int_t ip = 0; ip < kNPlanes; ip++) {
1260 x = tracklets[ip].GetX0();
1261 Double_t tmp = TMath::Sqrt(R*R-(x-x0)*(x-x0));
1263 // y: R^2 = (x - x0)^2 + (y - y0)^2
1264 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1265 tracklets[ip].SetYref(0, y0 - (y0>0.?1.:-1)*tmp);
1266 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1267 tracklets[ip].SetYref(1, (x - x0) / tmp);
1268 tracklets[ip].SetZref(0, z0 + dzdx * (x - xref));
1269 tracklets[ip].SetZref(1, dzdx);
1270 tracklets[ip].SetC(C);
1271 tracklets[ip].SetChi2(chi2);
1275 //update track points array
1278 for(int ip=0; ip<np; ip++){
1279 points[ip].GetXYZ(xyz);
1280 xyz[1] = y0 - (y0>0.?1.:-1)*TMath::Sqrt(R*R-(xyz[0]-x0)*(xyz[0]-x0));
1281 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1282 points[ip].SetXYZ(xyz);
1286 /* if(fReconstructor->GetStreamLevel() >=5){
1287 TTreeSRedirector &cstreamer = *fgDebugStreamer;
1288 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1289 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1290 Double_t chi2z = CalculateChi2Z(tracklets, z0, dzdx, xref);
1291 cstreamer << "FitRiemanTilt"
1292 << "EventNumber=" << eventNumber
1293 << "CandidateNumber=" << candidateNumber
1295 << "Chi2Z=" << chi2z
1302 //____________________________________________________________________
1303 Double_t AliTRDtrackerV1::FitKalman(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t up, Int_t np, AliTrackPoint *points)
1305 // Kalman filter implementation for the TRD.
1306 // It returns the positions of the fit in the array "points"
1308 // Author : A.Bercuci@gsi.de
1310 //printf("Start track @ x[%f]\n", track->GetX());
1312 //prepare marker points along the track
1313 Int_t ip = np ? 0 : 1;
1315 if((up?-1:1) * (track->GetX() - points[ip].GetX()) > 0.) break;
1316 //printf("AliTRDtrackerV1::FitKalman() : Skip track marker x[%d] = %7.3f. Before track start ( %7.3f ).\n", ip, points[ip].GetX(), track->GetX());
1319 //if(points) printf("First marker point @ x[%d] = %f\n", ip, points[ip].GetX());
1322 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
1324 //Loop through the TRD planes
1325 for (Int_t jplane = 0; jplane < kNPlanes; jplane++) {
1326 // GET TRACKLET OR BUILT IT
1327 Int_t iplane = up ? jplane : kNPlanes - 1 - jplane;
1329 if(!(ptrTracklet = &tracklets[iplane])) continue;
1331 if(!(ptrTracklet = track->GetTracklet(iplane))){
1332 /*AliTRDtrackerV1 *tracker = 0x0;
1333 if(!(tracker = dynamic_cast<AliTRDtrackerV1*>( AliTRDReconstructor::Tracker()))) continue;
1334 ptrTracklet = new(&tracklet) AliTRDseedV1(iplane);
1335 if(!tracker->MakeTracklet(ptrTracklet, track)) */
1339 if(!ptrTracklet->IsOK()) continue;
1341 Double_t x = ptrTracklet->GetX0();
1344 //don't do anything if next marker is after next update point.
1345 if((up?-1:1) * (points[ip].GetX() - x) - fgkMaxStep < 0) break;
1347 //printf("Propagate to x[%d] = %f\n", ip, points[ip].GetX());
1349 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1351 Double_t xyz[3]; // should also get the covariance
1352 track->GetXYZ(xyz); points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1355 //printf("plane[%d] tracklet[%p] x[%f]\n", iplane, ptrTracklet, x);
1357 //Propagate closer to the next update point
1358 if(((up?-1:1) * (x - track->GetX()) + fgkMaxStep < 0) && !PropagateToX(*track, x + (up?-1:1)*fgkMaxStep, fgkMaxStep)) return -1.;
1360 if(!AdjustSector(track)) return -1;
1361 if(TMath::Abs(track->GetSnp()) > fgkMaxSnp) return -1;
1363 //load tracklet to the tracker and the track
1365 if((index = FindTracklet(ptrTracklet)) < 0){
1366 ptrTracklet = SetTracklet(&tracklet);
1367 index = fTracklets->GetEntriesFast()-1;
1369 track->SetTracklet(ptrTracklet, index);*/
1372 // register tracklet to track with tracklet creation !!
1373 // PropagateBack : loaded tracklet to the tracker and update index
1374 // RefitInward : update index
1375 // MakeTrack : loaded tracklet to the tracker and update index
1376 if(!tracklets) track->SetTracklet(ptrTracklet, -1);
1379 //Calculate the mean material budget along the path inside the chamber
1380 Double_t xyz0[3]; track->GetXYZ(xyz0);
1381 Double_t alpha = track->GetAlpha();
1382 Double_t xyz1[3], y, z;
1383 if(!track->GetProlongation(x, y, z)) return -1;
1384 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
1385 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
1388 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
1389 Double_t xrho = param[0]*param[4]; // density*length
1390 Double_t xx0 = param[1]; // radiation length
1392 //Propagate the track
1393 track->PropagateTo(x, xx0, xrho);
1394 if (!AdjustSector(track)) break;
1397 Double_t chi2 = track->GetPredictedChi2(ptrTracklet);
1398 if(chi2<1e+10) track->Update(ptrTracklet, chi2);
1402 //Reset material budget if 2 consecutive gold
1403 if(iplane>0 && track->GetTracklet(iplane-1) && ptrTracklet->GetN() + track->GetTracklet(iplane-1)->GetN() > 20) track->SetBudget(2, 0.);
1404 } // end planes loop
1408 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1410 Double_t xyz[3]; // should also get the covariance
1411 track->GetXYZ(xyz); points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1415 return track->GetChi2();
1418 //_________________________________________________________________________
1419 Float_t AliTRDtrackerV1::CalculateChi2Z(AliTRDseedV1 *tracklets, Double_t offset, Double_t slope, Double_t xref)
1422 // Calculates the chi2-value of the track in z-Direction including tilting pad correction.
1423 // A linear dependence on the x-value serves as a model.
1424 // The parameters are related to the tilted Riemann fit.
1425 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
1426 // - the offset for the reference x
1428 // - the reference x position
1429 // Output: - The Chi2 value of the track in z-Direction
1431 Float_t chi2Z = 0, nLayers = 0;
1432 for (Int_t iLayer = 0; iLayer < AliTRDgeometry::kNlayer; iLayer++) {
1433 if(!tracklets[iLayer].IsOK()) continue;
1434 Double_t z = offset + slope * (tracklets[iLayer].GetX0() - xref);
1435 chi2Z += TMath::Abs(tracklets[iLayer].GetMeanz() - z);
1438 chi2Z /= TMath::Max((nLayers - 3.0),1.0);
1442 //_____________________________________________________________________________
1443 Int_t AliTRDtrackerV1::PropagateToX(AliTRDtrackV1 &t, Double_t xToGo, Double_t maxStep)
1446 // Starting from current X-position of track <t> this function
1447 // extrapolates the track up to radial position <xToGo>.
1448 // Returns 1 if track reaches the plane, and 0 otherwise
1451 const Double_t kEpsilon = 0.00001;
1453 // Current track X-position
1454 Double_t xpos = t.GetX();
1456 // Direction: inward or outward
1457 Double_t dir = (xpos < xToGo) ? 1.0 : -1.0;
1459 while (((xToGo - xpos) * dir) > kEpsilon) {
1468 // The next step size
1469 Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
1471 // Get the global position of the starting point
1474 // X-position after next step
1477 // Get local Y and Z at the X-position of the next step
1478 if (!t.GetProlongation(x,y,z)) {
1479 return 0; // No prolongation possible
1482 // The global position of the end point of this prolongation step
1483 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1484 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1487 // Calculate the mean material budget between start and
1488 // end point of this prolongation step
1489 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
1491 // Propagate the track to the X-position after the next step
1492 if (!t.PropagateTo(x,param[1],param[0]*param[4])) {
1496 // Rotate the track if necessary
1499 // New track X-position
1509 //_____________________________________________________________________________
1510 Int_t AliTRDtrackerV1::ReadClusters(TClonesArray* &array, TTree *clusterTree) const
1513 // Reads AliTRDclusters from the file.
1514 // The names of the cluster tree and branches
1515 // should match the ones used in AliTRDclusterizer::WriteClusters()
1518 Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
1519 TObjArray *clusterArray = new TObjArray(nsize+1000);
1521 TBranch *branch = clusterTree->GetBranch("TRDcluster");
1523 AliError("Can't get the branch !");
1526 branch->SetAddress(&clusterArray);
1529 array = new TClonesArray("AliTRDcluster", nsize);
1530 array->SetOwner(kTRUE);
1533 // Loop through all entries in the tree
1534 Int_t nEntries = (Int_t) clusterTree->GetEntries();
1537 AliTRDcluster *c = 0x0;
1538 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1540 nbytes += clusterTree->GetEvent(iEntry);
1542 // Get the number of points in the detector
1543 Int_t nCluster = clusterArray->GetEntriesFast();
1544 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1545 if(!(c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster))) continue;
1547 new((*fClusters)[ncl++]) AliTRDcluster(*c);
1548 delete (clusterArray->RemoveAt(iCluster));
1552 delete clusterArray;
1557 //_____________________________________________________________________________
1558 Int_t AliTRDtrackerV1::LoadClusters(TTree *cTree)
1561 // Fills clusters into TRD tracking sectors
1564 if(!(fClusters = AliTRDReconstructor::GetClusters())){
1565 if (ReadClusters(fClusters, cTree)) {
1566 AliError("Problem with reading the clusters !");
1572 if(!fClusters->GetEntriesFast()){
1573 AliInfo("No TRD clusters");
1578 BuildTrackingContainers();
1580 //Int_t ncl = fClusters->GetEntriesFast();
1581 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1586 //_____________________________________________________________________________
1587 Int_t AliTRDtrackerV1::LoadClusters(TClonesArray *clusters)
1590 // Fills clusters into TRD tracking sectors
1591 // Function for use in the HLT
1593 if(!clusters || !clusters->GetEntriesFast()){
1594 AliInfo("No TRD clusters");
1598 fClusters = clusters;
1602 BuildTrackingContainers();
1604 //Int_t ncl = fClusters->GetEntriesFast();
1605 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1611 //____________________________________________________________________
1612 Int_t AliTRDtrackerV1::BuildTrackingContainers()
1614 // Building tracking containers for clusters
1616 Int_t nin =0, icl = fClusters->GetEntriesFast();
1618 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(icl);
1619 if(c->IsInChamber()) nin++;
1620 Int_t detector = c->GetDetector();
1621 Int_t sector = fGeom->GetSector(detector);
1622 Int_t stack = fGeom->GetStack(detector);
1623 Int_t layer = fGeom->GetLayer(detector);
1625 fTrSec[sector].GetChamber(stack, layer, kTRUE)->InsertCluster(c, icl);
1628 for(int isector =0; isector<AliTRDgeometry::kNsector; isector++){
1629 if(!fTrSec[isector].GetNChambers()) continue;
1630 fTrSec[isector].Init(fReconstructor);
1638 //____________________________________________________________________
1639 void AliTRDtrackerV1::UnloadClusters()
1642 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1645 if(fTracks) fTracks->Delete();
1646 if(fTracklets) fTracklets->Delete();
1647 if(fClusters && IsClustersOwner()) fClusters->Delete();
1649 for (int i = 0; i < AliTRDgeometry::kNsector; i++) fTrSec[i].Clear();
1651 // Increment the Event Number
1652 AliTRDtrackerDebug::SetEventNumber(AliTRDtrackerDebug::GetEventNumber() + 1);
1655 //_____________________________________________________________________________
1656 Bool_t AliTRDtrackerV1::AdjustSector(AliTRDtrackV1 *track)
1659 // Rotates the track when necessary
1662 Double_t alpha = AliTRDgeometry::GetAlpha();
1663 Double_t y = track->GetY();
1664 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
1667 if (!track->Rotate( alpha)) {
1671 else if (y < -ymax) {
1672 if (!track->Rotate(-alpha)) {
1682 //____________________________________________________________________
1683 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(AliTRDtrackV1 *track, Int_t p, Int_t &idx)
1685 // Find tracklet for TRD track <track>
1694 // Detailed description
1696 idx = track->GetTrackletIndex(p);
1697 AliTRDseedV1 *tracklet = (idx==0xffff) ? 0x0 : (AliTRDseedV1*)fTracklets->UncheckedAt(idx);
1702 //____________________________________________________________________
1703 AliTRDseedV1* AliTRDtrackerV1::SetTracklet(AliTRDseedV1 *tracklet)
1705 // Add this tracklet to the list of tracklets stored in the tracker
1708 // - tracklet : pointer to the tracklet to be added to the list
1711 // - the index of the new tracklet in the tracker tracklets list
1713 // Detailed description
1714 // Build the tracklets list if it is not yet created (late initialization)
1715 // and adds the new tracklet to the list.
1718 fTracklets = new TClonesArray("AliTRDseedV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1719 fTracklets->SetOwner(kTRUE);
1721 Int_t nentries = fTracklets->GetEntriesFast();
1722 return new ((*fTracklets)[nentries]) AliTRDseedV1(*tracklet);
1725 //____________________________________________________________________
1726 AliTRDtrackV1* AliTRDtrackerV1::SetTrack(AliTRDtrackV1 *track)
1728 // Add this track to the list of tracks stored in the tracker
1731 // - track : pointer to the track to be added to the list
1734 // - the pointer added
1736 // Detailed description
1737 // Build the tracks list if it is not yet created (late initialization)
1738 // and adds the new track to the list.
1741 fTracks = new TClonesArray("AliTRDtrackV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1742 fTracks->SetOwner(kTRUE);
1744 Int_t nentries = fTracks->GetEntriesFast();
1745 return new ((*fTracks)[nentries]) AliTRDtrackV1(*track);
1750 //____________________________________________________________________
1751 Int_t AliTRDtrackerV1::Clusters2TracksSM(Int_t sector, AliESDEvent *esd)
1754 // Steer tracking for one SM.
1757 // sector : Array of (SM) propagation layers containing clusters
1758 // esd : The current ESD event. On output it contains the also
1759 // the ESD (TRD) tracks found in this SM.
1762 // Number of tracks found in this TRD supermodule.
1764 // Detailed description
1766 // 1. Unpack AliTRDpropagationLayers objects for each stack.
1767 // 2. Launch stack tracking.
1768 // See AliTRDtrackerV1::Clusters2TracksStack() for details.
1769 // 3. Pack results in the ESD event.
1772 // allocate space for esd tracks in this SM
1773 TClonesArray esdTrackList("AliESDtrack", 2*kMaxTracksStack);
1774 esdTrackList.SetOwner();
1777 Int_t nChambers = 0;
1778 AliTRDtrackingChamber **stack = 0x0, *chamber = 0x0;
1779 for(int istack = 0; istack<AliTRDgeometry::kNstack; istack++){
1780 if(!(stack = fTrSec[sector].GetStack(istack))) continue;
1782 for(int ilayer=0; ilayer<AliTRDgeometry::kNlayer; ilayer++){
1783 if(!(chamber = stack[ilayer])) continue;
1784 if(chamber->GetNClusters() < fgNTimeBins * fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
1786 //AliInfo(Form("sector %d stack %d layer %d clusters %d", sector, istack, ilayer, chamber->GetNClusters()));
1788 if(nChambers < 4) continue;
1789 //AliInfo(Form("Doing stack %d", istack));
1790 nTracks += Clusters2TracksStack(stack, &esdTrackList);
1792 //AliInfo(Form("Found %d tracks in SM %d [%d]\n", nTracks, sector, esd->GetNumberOfTracks()));
1794 for(int itrack=0; itrack<nTracks; itrack++)
1795 esd->AddTrack((AliESDtrack*)esdTrackList[itrack]);
1797 // Reset Track and Candidate Number
1798 AliTRDtrackerDebug::SetCandidateNumber(0);
1799 AliTRDtrackerDebug::SetTrackNumber(0);
1803 //____________________________________________________________________
1804 Int_t AliTRDtrackerV1::Clusters2TracksStack(AliTRDtrackingChamber **stack, TClonesArray *esdTrackList)
1807 // Make tracks in one TRD stack.
1810 // layer : Array of stack propagation layers containing clusters
1811 // esdTrackList : Array of ESD tracks found by the stand alone tracker.
1812 // On exit the tracks found in this stack are appended.
1815 // Number of tracks found in this stack.
1817 // Detailed description
1819 // 1. Find the 3 most useful seeding chambers. See BuildSeedingConfigs() for details.
1820 // 2. Steer AliTRDtrackerV1::MakeSeeds() for 3 seeding layer configurations.
1821 // See AliTRDtrackerV1::MakeSeeds() for more details.
1822 // 3. Arrange track candidates in decreasing order of their quality
1823 // 4. Classify tracks in 5 categories according to:
1824 // a) number of layers crossed
1826 // 5. Sign clusters by tracks in decreasing order of track quality
1827 // 6. Build AliTRDtrack out of seeding tracklets
1829 // 8. Build ESD track and register it to the output list
1832 AliTRDtrackingChamber *chamber = 0x0;
1833 AliTRDseedV1 sseed[kMaxTracksStack*6]; // to be initialized
1834 Int_t pars[4]; // MakeSeeds parameters
1836 //Double_t alpha = AliTRDgeometry::GetAlpha();
1837 //Double_t shift = .5 * alpha;
1838 Int_t configs[kNConfigs];
1840 // Build initial seeding configurations
1841 Double_t quality = BuildSeedingConfigs(stack, configs);
1842 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
1843 AliInfo(Form("Plane config %d %d %d Quality %f"
1844 , configs[0], configs[1], configs[2], quality));
1847 // Initialize contors
1848 Int_t ntracks, // number of TRD track candidates
1849 ntracks1, // number of registered TRD tracks/iter
1850 ntracks2 = 0; // number of all registered TRD tracks in stack
1853 // Loop over seeding configurations
1854 ntracks = 0; ntracks1 = 0;
1855 for (Int_t iconf = 0; iconf<3; iconf++) {
1856 pars[0] = configs[iconf];
1858 ntracks = MakeSeeds(stack, &sseed[6*ntracks], pars);
1859 if(ntracks == kMaxTracksStack) break;
1861 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1) AliInfo(Form("Candidate TRD tracks %d in iteration %d.", ntracks, fSieveSeeding));
1865 // Sort the seeds according to their quality
1866 Int_t sort[kMaxTracksStack];
1867 TMath::Sort(ntracks, fTrackQuality, sort, kTRUE);
1869 // Initialize number of tracks so far and logic switches
1870 Int_t ntracks0 = esdTrackList->GetEntriesFast();
1871 Bool_t signedTrack[kMaxTracksStack];
1872 Bool_t fakeTrack[kMaxTracksStack];
1873 for (Int_t i=0; i<ntracks; i++){
1874 signedTrack[i] = kFALSE;
1875 fakeTrack[i] = kFALSE;
1877 //AliInfo("Selecting track candidates ...");
1879 // Sieve clusters in decreasing order of track quality
1880 Double_t trackParams[7];
1881 // AliTRDseedV1 *lseed = 0x0;
1882 Int_t jSieve = 0, candidates;
1884 //AliInfo(Form("\t\tITER = %i ", jSieve));
1886 // Check track candidates
1888 for (Int_t itrack = 0; itrack < ntracks; itrack++) {
1889 Int_t trackIndex = sort[itrack];
1890 if (signedTrack[trackIndex] || fakeTrack[trackIndex]) continue;
1893 // Calculate track parameters from tracklets seeds
1894 Int_t labelsall[1000];
1895 Int_t nlabelsall = 0;
1896 Int_t naccepted = 0;
1901 for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
1902 Int_t jseed = kNPlanes*trackIndex+jLayer;
1903 if(!sseed[jseed].IsOK()) continue;
1904 if (TMath::Abs(sseed[jseed].GetYref(0) / sseed[jseed].GetX0()) < 0.15) findable++;
1906 sseed[jseed].UpdateUsed();
1907 ncl += sseed[jseed].GetN2();
1908 nused += sseed[jseed].GetNUsed();
1912 // for (Int_t itime = 0; itime < fgNTimeBins; itime++) {
1913 // if(!sseed[jseed].IsUsable(itime)) continue;
1915 // Int_t tindex = 0, ilab = 0;
1916 // while(ilab<3 && (tindex = sseed[jseed].GetClusters(itime)->GetLabel(ilab)) >= 0){
1917 // labelsall[nlabelsall++] = tindex;
1923 // Filter duplicated tracks
1925 //printf("Skip %d nused %d\n", trackIndex, nused);
1926 fakeTrack[trackIndex] = kTRUE;
1929 if (Float_t(nused)/ncl >= .25){
1930 //printf("Skip %d nused/ncl >= .25\n", trackIndex);
1931 fakeTrack[trackIndex] = kTRUE;
1936 Bool_t skip = kFALSE;
1939 if(nlayers < 6) {skip = kTRUE; break;}
1940 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
1944 if(nlayers < findable){skip = kTRUE; break;}
1945 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -4.){skip = kTRUE; break;}
1949 if ((nlayers == findable) || (nlayers == 6)) { skip = kTRUE; break;}
1950 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -6.0){skip = kTRUE; break;}
1954 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
1958 if (nlayers == 3){skip = kTRUE; break;}
1959 //if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) - nused/(nlayers-3.0) < -15.0){skip = kTRUE; break;}
1964 //printf("REJECTED : %d [%d] nlayers %d trackQuality = %e nused %d\n", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused);
1967 signedTrack[trackIndex] = kTRUE;
1971 AliTRDcluster *cl = 0x0; Int_t clusterIndex = -1;
1972 for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
1973 Int_t jseed = kNPlanes*trackIndex+jLayer;
1974 if(!sseed[jseed].IsOK()) continue;
1975 if(TMath::Abs(sseed[jseed].GetYfit(1) - sseed[jseed].GetYfit(1)) >= .2) continue; // check this condition with Marian
1976 sseed[jseed].UseClusters();
1979 while(!(cl = sseed[jseed].GetClusters(ic))) ic++;
1980 clusterIndex = sseed[jseed].GetIndexes(ic);
1986 // Build track parameters
1987 AliTRDseedV1 *lseed =&sseed[trackIndex*6];
1989 while(idx<3 && !lseed->IsOK()) {
1993 Double_t x = lseed->GetX0();// - 3.5;
1994 trackParams[0] = x; //NEW AB
1995 trackParams[1] = lseed->GetYref(0); // lseed->GetYat(x);
1996 trackParams[2] = lseed->GetZref(0); // lseed->GetZat(x);
1997 trackParams[3] = TMath::Sin(TMath::ATan(lseed->GetYref(1)));
1998 trackParams[4] = lseed->GetZref(1) / TMath::Sqrt(1. + lseed->GetYref(1) * lseed->GetYref(1));
1999 trackParams[5] = lseed->GetC();
2000 Int_t ich = 0; while(!(chamber = stack[ich])) ich++;
2001 trackParams[6] = fGeom->GetSector(chamber->GetDetector());/* *alpha+shift; // Supermodule*/
2003 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2004 AliInfo(Form("Track %d [%d] nlayers %d trackQuality = %e nused %d, yref = %3.3f", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused, trackParams[1]));
2006 Int_t nclusters = 0;
2007 AliTRDseedV1 *dseed[6];
2009 // Build track label - what happens if measured data ???
2013 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
2014 Int_t jseed = kNPlanes*trackIndex+iLayer;
2015 dseed[iLayer] = new AliTRDseedV1(sseed[jseed]);
2016 dseed[iLayer]->SetOwner();
2017 nclusters += sseed[jseed].GetN2();
2018 if(!sseed[jseed].IsOK()) continue;
2019 for(int ilab=0; ilab<2; ilab++){
2020 if(sseed[jseed].GetLabels(ilab) < 0) continue;
2021 labels[nlab] = sseed[jseed].GetLabels(ilab);
2025 Freq(nlab,labels,outlab,kFALSE);
2026 Int_t label = outlab[0];
2027 Int_t frequency = outlab[1];
2028 Freq(nlabelsall,labelsall,outlab,kFALSE);
2029 Int_t label1 = outlab[0];
2030 Int_t label2 = outlab[2];
2031 Float_t fakeratio = (naccepted - outlab[1]) / Float_t(naccepted);
2033 //Int_t eventNrInFile = esd->GetEventNumberInFile();
2034 //AliInfo(Form("Number of clusters %d.", nclusters));
2035 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2036 Int_t trackNumber = AliTRDtrackerDebug::GetTrackNumber();
2037 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2038 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2039 cstreamer << "Clusters2TracksStack"
2040 << "EventNumber=" << eventNumber
2041 << "TrackNumber=" << trackNumber
2042 << "CandidateNumber=" << candidateNumber
2043 << "Iter=" << fSieveSeeding
2044 << "Like=" << fTrackQuality[trackIndex]
2045 << "S0.=" << dseed[0]
2046 << "S1.=" << dseed[1]
2047 << "S2.=" << dseed[2]
2048 << "S3.=" << dseed[3]
2049 << "S4.=" << dseed[4]
2050 << "S5.=" << dseed[5]
2051 << "p0=" << trackParams[0]
2052 << "p1=" << trackParams[1]
2053 << "p2=" << trackParams[2]
2054 << "p3=" << trackParams[3]
2055 << "p4=" << trackParams[4]
2056 << "p5=" << trackParams[5]
2057 << "p6=" << trackParams[6]
2058 << "Label=" << label
2059 << "Label1=" << label1
2060 << "Label2=" << label2
2061 << "FakeRatio=" << fakeratio
2062 << "Freq=" << frequency
2064 << "NLayers=" << nlayers
2065 << "Findable=" << findable
2066 << "NUsed=" << nused
2070 AliTRDtrackV1 *track = MakeTrack(&sseed[trackIndex*kNPlanes], trackParams);
2072 AliWarning("Fail to build a TRD Track.");
2075 //AliInfo("End of MakeTrack()");
2076 AliESDtrack *esdTrack = new ((*esdTrackList)[ntracks0++]) AliESDtrack();
2077 esdTrack->UpdateTrackParams(track, AliESDtrack::kTRDout);
2078 esdTrack->SetLabel(track->GetLabel());
2079 track->UpdateESDtrack(esdTrack);
2080 // write ESD-friends if neccessary
2081 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0){
2082 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(*track);
2083 calibTrack->SetOwner();
2084 esdTrack->AddCalibObject(calibTrack);
2087 AliTRDtrackerDebug::SetTrackNumber(AliTRDtrackerDebug::GetTrackNumber() + 1);
2091 } while(jSieve<5 && candidates); // end track candidates sieve
2092 if(!ntracks1) break;
2094 // increment counters
2095 ntracks2 += ntracks1;
2097 if(fReconstructor->IsHLT()) break;
2100 // Rebuild plane configurations and indices taking only unused clusters into account
2101 quality = BuildSeedingConfigs(stack, configs);
2102 if(quality < 1.E-7) break; //fReconstructor->GetRecoParam() ->GetPlaneQualityThreshold()) break;
2104 for(Int_t ip = 0; ip < kNPlanes; ip++){
2105 if(!(chamber = stack[ip])) continue;
2106 chamber->Build(fGeom);//Indices(fSieveSeeding);
2109 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2110 AliInfo(Form("Sieve level %d Plane config %d %d %d Quality %f", fSieveSeeding, configs[0], configs[1], configs[2], quality));
2112 } while(fSieveSeeding<10); // end stack clusters sieve
2116 //AliInfo(Form("Registered TRD tracks %d in stack %d.", ntracks2, pars[1]));
2121 //___________________________________________________________________
2122 Double_t AliTRDtrackerV1::BuildSeedingConfigs(AliTRDtrackingChamber **stack, Int_t *configs)
2125 // Assign probabilities to chambers according to their
2126 // capability of producing seeds.
2130 // layers : Array of stack propagation layers for all 6 chambers in one stack
2131 // configs : On exit array of configuration indexes (see GetSeedingConfig()
2132 // for details) in the decreasing order of their seeding probabilities.
2136 // Return top configuration quality
2138 // Detailed description:
2140 // To each chamber seeding configuration (see GetSeedingConfig() for
2141 // the list of all configurations) one defines 2 quality factors:
2142 // - an apriori topological quality (see GetSeedingConfig() for details) and
2143 // - a data quality based on the uniformity of the distribution of
2144 // clusters over the x range (time bins population). See CookChamberQA() for details.
2145 // The overall chamber quality is given by the product of this 2 contributions.
2148 Double_t chamberQ[kNPlanes];
2149 AliTRDtrackingChamber *chamber = 0x0;
2150 for(int iplane=0; iplane<kNPlanes; iplane++){
2151 if(!(chamber = stack[iplane])) continue;
2152 chamberQ[iplane] = (chamber = stack[iplane]) ? chamber->GetQuality() : 0.;
2155 Double_t tconfig[kNConfigs];
2157 for(int iconf=0; iconf<kNConfigs; iconf++){
2158 GetSeedingConfig(iconf, planes);
2159 tconfig[iconf] = fgTopologicQA[iconf];
2160 for(int iplane=0; iplane<4; iplane++) tconfig[iconf] *= chamberQ[planes[iplane]];
2163 TMath::Sort((Int_t)kNConfigs, tconfig, configs, kTRUE);
2164 // AliInfo(Form("q[%d] = %f", configs[0], tconfig[configs[0]]));
2165 // AliInfo(Form("q[%d] = %f", configs[1], tconfig[configs[1]]));
2166 // AliInfo(Form("q[%d] = %f", configs[2], tconfig[configs[2]]));
2168 return tconfig[configs[0]];
2171 //____________________________________________________________________
2172 Int_t AliTRDtrackerV1::MakeSeeds(AliTRDtrackingChamber **stack, AliTRDseedV1 *sseed, Int_t *ipar)
2175 // Make tracklet seeds in the TRD stack.
2178 // layers : Array of stack propagation layers containing clusters
2179 // sseed : Array of empty tracklet seeds. On exit they are filled.
2180 // ipar : Control parameters:
2181 // ipar[0] -> seeding chambers configuration
2182 // ipar[1] -> stack index
2183 // ipar[2] -> number of track candidates found so far
2186 // Number of tracks candidates found.
2188 // Detailed description
2190 // The following steps are performed:
2191 // 1. Select seeding layers from seeding chambers
2192 // 2. Select seeding clusters from the seeding AliTRDpropagationLayerStack.
2193 // The clusters are taken from layer 3, layer 0, layer 1 and layer 2, in
2194 // this order. The parameters controling the range of accepted clusters in
2195 // layer 0, 1, and 2 are defined in AliTRDchamberTimeBin::BuildCond().
2196 // 3. Helix fit of the cluster set. (see AliTRDtrackerFitter::FitRieman(AliTRDcluster**))
2197 // 4. Initialize seeding tracklets in the seeding chambers.
2199 // Chi2 in the Y direction less than threshold ... (1./(3. - sLayer))
2200 // Chi2 in the Z direction less than threshold ... (1./(3. - sLayer))
2201 // 6. Attach clusters to seeding tracklets and find linear approximation of
2202 // the tracklet (see AliTRDseedV1::AttachClustersIter()). The number of used
2203 // clusters used by current seeds should not exceed ... (25).
2205 // All 4 seeding tracklets should be correctly constructed (see
2206 // AliTRDseedV1::AttachClustersIter())
2207 // 8. Helix fit of the seeding tracklets
2209 // Likelihood calculation of the fit. (See AliTRDtrackerV1::CookLikelihood() for details)
2210 // 10. Extrapolation of the helix fit to the other 2 chambers:
2211 // a) Initialization of extrapolation tracklet with fit parameters
2212 // b) Helix fit of tracklets
2213 // c) Attach clusters and linear interpolation to extrapolated tracklets
2214 // d) Helix fit of tracklets
2215 // 11. Improve seeding tracklets quality by reassigning clusters.
2216 // See AliTRDtrackerV1::ImproveSeedQuality() for details.
2217 // 12. Helix fit of all 6 seeding tracklets and chi2 calculation
2218 // 13. Hyperplane fit and track quality calculation. See AliTRDtrackerFitter::FitHyperplane() for details.
2219 // 14. Cooking labels for tracklets. Should be done only for MC
2220 // 15. Register seeds.
2223 AliTRDtrackingChamber *chamber = 0x0;
2224 AliTRDcluster *c[kNSeedPlanes] = {0x0, 0x0, 0x0, 0x0}; // initilize seeding clusters
2225 AliTRDseedV1 *cseed = &sseed[0]; // initialize tracklets for first track
2226 Int_t ncl, mcl; // working variable for looping over clusters
2227 Int_t index[AliTRDchamberTimeBin::kMaxClustersLayer], jndex[AliTRDchamberTimeBin::kMaxClustersLayer];
2229 // chi2[0] = tracklet chi2 on the Z direction
2230 // chi2[1] = tracklet chi2 on the R direction
2233 // Default positions for the anode wire in all 6 Layers in case of a stack with missing clusters
2234 // Positions taken using cosmic data taken with SM3 after rebuild
2235 Double_t x_def[kNPlanes] = {300.2, 312.8, 325.4, 338, 350.6, 363.2};
2237 // this should be data member of AliTRDtrack
2238 Double_t seedQuality[kMaxTracksStack];
2240 // unpack control parameters
2241 Int_t config = ipar[0];
2242 Int_t ntracks = ipar[1];
2243 Int_t planes[kNSeedPlanes]; GetSeedingConfig(config, planes);
2244 Int_t planesExt[kNPlanes-kNSeedPlanes]; GetExtrapolationConfig(config, planesExt);
2247 // Init chambers geometry
2248 Int_t ic = 0; while(!(chamber = stack[ic])) ic++;
2249 Int_t istack = fGeom->GetStack(chamber->GetDetector());
2250 Double_t hL[kNPlanes]; // Tilting angle
2251 Float_t padlength[kNPlanes]; // pad lenghts
2252 AliTRDpadPlane *pp = 0x0;
2253 for(int iplane=0; iplane<kNPlanes; iplane++){
2254 pp = fGeom->GetPadPlane(iplane, istack);
2255 hL[iplane] = TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle());
2256 padlength[iplane] = pp->GetLengthIPad();
2259 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2260 AliInfo(Form("Making seeds Stack[%d] Config[%d] Tracks[%d]...", istack, config, ntracks));
2265 for(int isl=0; isl<kNSeedPlanes; isl++){
2266 if(!(chamber = stack[planes[isl]])) continue;
2267 if(!chamber->GetSeedingLayer(fSeedTB[isl], fGeom, fReconstructor)) continue;
2270 if(nlayers < 4) return 0;
2273 // Start finding seeds
2274 Double_t cond0[4], cond1[4], cond2[4];
2276 while((c[3] = (*fSeedTB[3])[icl++])){
2278 fSeedTB[0]->BuildCond(c[3], cond0, 0);
2279 fSeedTB[0]->GetClusters(cond0, index, ncl);
2280 //printf("Found c[3] candidates 0 %d\n", ncl);
2283 c[0] = (*fSeedTB[0])[index[jcl++]];
2285 Double_t dx = c[3]->GetX() - c[0]->GetX();
2286 Double_t theta = (c[3]->GetZ() - c[0]->GetZ())/dx;
2287 Double_t phi = (c[3]->GetY() - c[0]->GetY())/dx;
2288 fSeedTB[1]->BuildCond(c[0], cond1, 1, theta, phi);
2289 fSeedTB[1]->GetClusters(cond1, jndex, mcl);
2290 //printf("Found c[0] candidates 1 %d\n", mcl);
2294 c[1] = (*fSeedTB[1])[jndex[kcl++]];
2296 fSeedTB[2]->BuildCond(c[1], cond2, 2, theta, phi);
2297 c[2] = fSeedTB[2]->GetNearestCluster(cond2);
2298 //printf("Found c[1] candidate 2 %p\n", c[2]);
2301 // AliInfo("Seeding clusters found. Building seeds ...");
2302 // 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());
2304 for (Int_t il = 0; il < kNPlanes; il++) cseed[il].Reset();
2308 AliTRDseedV1 *tseed = 0x0;
2309 for(int iLayer=0; iLayer<kNPlanes; iLayer++){
2310 tseed = &cseed[iLayer];
2311 tseed->SetPlane(iLayer);
2312 tseed->SetTilt(hL[iLayer]);
2313 tseed->SetPadLength(padlength[iLayer]);
2314 tseed->SetReconstructor(fReconstructor);
2315 Double_t x_anode = stack[iLayer] ? stack[iLayer]->GetX() : x_def[iLayer];
2316 tseed->SetX0(x_anode);
2317 tseed->Init(GetRiemanFitter());
2320 Bool_t isFake = kFALSE;
2321 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2322 if (c[0]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2323 if (c[1]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2324 if (c[2]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2327 for(Int_t l = 0; l < kNSeedPlanes; l++) xpos[l] = fSeedTB[l]->GetX();
2329 for(int il=0; il<4; il++) yref[il] = cseed[planes[il]].GetYref(0);
2330 Int_t ll = c[3]->GetLabel(0);
2331 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2332 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2333 AliRieman *rim = GetRiemanFitter();
2334 TTreeSRedirector &cs0 = *fgDebugStreamer;
2336 <<"EventNumber=" << eventNumber
2337 <<"CandidateNumber=" << candidateNumber
2338 <<"isFake=" << isFake
2339 <<"config=" << config
2341 <<"chi2z=" << chi2[0]
2342 <<"chi2y=" << chi2[1]
2343 <<"Y2exp=" << cond2[0]
2344 <<"Z2exp=" << cond2[1]
2345 <<"X0=" << xpos[0] //layer[sLayer]->GetX()
2346 <<"X1=" << xpos[1] //layer[sLayer + 1]->GetX()
2347 <<"X2=" << xpos[2] //layer[sLayer + 2]->GetX()
2348 <<"X3=" << xpos[3] //layer[sLayer + 3]->GetX()
2349 <<"yref0=" << yref[0]
2350 <<"yref1=" << yref[1]
2351 <<"yref2=" << yref[2]
2352 <<"yref3=" << yref[3]
2357 <<"Seed0.=" << &cseed[planes[0]]
2358 <<"Seed1.=" << &cseed[planes[1]]
2359 <<"Seed2.=" << &cseed[planes[2]]
2360 <<"Seed3.=" << &cseed[planes[3]]
2361 <<"RiemanFitter.=" << rim
2365 if(chi2[0] > fReconstructor->GetRecoParam() ->GetChi2Z()/*7./(3. - sLayer)*//*iter*/){
2366 //AliInfo(Form("Failed chi2 filter on chi2Z [%f].", chi2[0]));
2367 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2370 if(chi2[1] > fReconstructor->GetRecoParam() ->GetChi2Y()/*1./(3. - sLayer)*//*iter*/){
2371 //AliInfo(Form("Failed chi2 filter on chi2Y [%f].", chi2[1]));
2372 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2375 //AliInfo("Passed chi2 filter.");
2377 // try attaching clusters to tracklets
2380 for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
2381 Int_t jLayer = planes[iLayer];
2382 if(!cseed[jLayer].AttachClustersIter(stack[jLayer], 5., kFALSE, c[iLayer])) continue;
2383 nUsedCl += cseed[jLayer].GetNUsed();
2384 if(nUsedCl > 25) break;
2388 if(mlayers < kNSeedPlanes){
2389 //AliInfo(Form("Failed updating all seeds %d [%d].", mlayers, kNSeedPlanes));
2390 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2394 // temporary exit door for the HLT
2395 if(fReconstructor->IsHLT()){
2396 // attach clusters to extrapolation chambers
2397 for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
2398 Int_t jLayer = planesExt[iLayer];
2399 if(!(chamber = stack[jLayer])) continue;
2400 cseed[jLayer].AttachClustersIter(chamber, 1000.);
2402 fTrackQuality[ntracks] = 1.; // dummy value
2404 if(ntracks == kMaxTracksStack){
2405 AliWarning(Form("Number of seeds reached maximum allowed (%d) in stack.", kMaxTracksStack));
2413 // fit tracklets and cook likelihood
2414 FitTiltedRieman(&cseed[0], kTRUE);// Update Seeds and calculate Likelihood
2415 chi2[0] = GetChi2Y(&cseed[0]);
2416 chi2[1] = GetChi2Z(&cseed[0]);
2417 //Chi2 definitions in testing stage
2418 //chi2[0] = GetChi2YTest(&cseed[0]);
2419 //chi2[1] = GetChi2ZTest(&cseed[0]);
2420 Double_t like = CookLikelihood(&cseed[0], planes, chi2); // to be checked
2422 if (TMath::Log(1.E-9 + like) < fReconstructor->GetRecoParam() ->GetTrackLikelihood()){
2423 //AliInfo(Form("Failed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2424 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2427 //AliInfo(Form("Passed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2429 // book preliminary results
2430 seedQuality[ntracks] = like;
2431 fSeedLayer[ntracks] = config;/*sLayer;*/
2433 // attach clusters to the extrapolation seeds
2434 Int_t nusedf = 0; // debug value
2435 for(int iLayer=0; iLayer<kNPlanes-kNSeedPlanes; iLayer++){
2436 Int_t jLayer = planesExt[iLayer];
2437 if(!(chamber = stack[jLayer])) continue;
2439 // fit extrapolated seed
2440 if ((jLayer == 0) && !(cseed[1].IsOK())) continue;
2441 if ((jLayer == 5) && !(cseed[4].IsOK())) continue;
2442 AliTRDseedV1 pseed = cseed[jLayer];
2443 if(!pseed.AttachClustersIter(chamber, 1000.)) continue;
2444 cseed[jLayer] = pseed;
2445 nusedf += cseed[jLayer].GetNUsed(); // debug value
2446 FitTiltedRieman(cseed, kTRUE);
2449 // AliInfo("Extrapolation done.");
2450 // Debug Stream containing all the 6 tracklets
2451 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2452 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2453 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2454 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2455 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2456 cstreamer << "MakeSeeds1"
2457 << "EventNumber=" << eventNumber
2458 << "CandidateNumber=" << candidateNumber
2459 << "S0.=" << &cseed[0]
2460 << "S1.=" << &cseed[1]
2461 << "S2.=" << &cseed[2]
2462 << "S3.=" << &cseed[3]
2463 << "S4.=" << &cseed[4]
2464 << "S5.=" << &cseed[5]
2465 << "FitterT.=" << tiltedRieman
2469 if(fReconstructor->GetRecoParam()->HasImproveTracklets() && ImproveSeedQuality(stack, cseed) < 4){
2470 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2473 //AliInfo("Improve seed quality done.");
2475 // fit full track and cook likelihoods
2476 // Double_t curv = FitRieman(&cseed[0], chi2);
2477 // Double_t chi2ZF = chi2[0] / TMath::Max((mlayers - 3.), 1.);
2478 // Double_t chi2RF = chi2[1] / TMath::Max((mlayers - 3.), 1.);
2480 // do the final track fitting (Once with vertex constraint and once without vertex constraint)
2481 Double_t chi2Vals[3];
2482 chi2Vals[0] = FitTiltedRieman(&cseed[0], kFALSE);
2483 if(fReconstructor->GetRecoParam() ->IsVertexConstrained())
2484 chi2Vals[1] = FitTiltedRiemanConstraint(&cseed[0], GetZ()); // Do Vertex Constrained fit if desired
2487 chi2Vals[2] = GetChi2Z(&cseed[0]) / TMath::Max((mlayers - 3.), 1.);
2488 // Chi2 definitions in testing stage
2489 //chi2Vals[2] = GetChi2ZTest(&cseed[0]);
2490 fTrackQuality[ntracks] = CalculateTrackLikelihood(&cseed[0], &chi2Vals[0]);
2491 //AliInfo("Hyperplane fit done\n");
2493 // finalize tracklets
2497 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2498 if (!cseed[iLayer].IsOK()) continue;
2500 if (cseed[iLayer].GetLabels(0) >= 0) {
2501 labels[nlab] = cseed[iLayer].GetLabels(0);
2505 if (cseed[iLayer].GetLabels(1) >= 0) {
2506 labels[nlab] = cseed[iLayer].GetLabels(1);
2510 Freq(nlab,labels,outlab,kFALSE);
2511 Int_t label = outlab[0];
2512 Int_t frequency = outlab[1];
2513 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2514 cseed[iLayer].SetFreq(frequency);
2515 cseed[iLayer].SetChi2Z(chi2[1]);
2518 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2519 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2520 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2521 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2522 TLinearFitter *fitterTC = GetTiltedRiemanFitterConstraint();
2523 TLinearFitter *fitterT = GetTiltedRiemanFitter();
2524 cstreamer << "MakeSeeds2"
2525 << "EventNumber=" << eventNumber
2526 << "CandidateNumber=" << candidateNumber
2527 << "Chi2TR=" << chi2Vals[0]
2528 << "Chi2TC=" << chi2Vals[1]
2529 << "Nlayers=" << mlayers
2530 << "NUsedS=" << nUsedCl
2531 << "NUsed=" << nusedf
2533 << "S0.=" << &cseed[0]
2534 << "S1.=" << &cseed[1]
2535 << "S2.=" << &cseed[2]
2536 << "S3.=" << &cseed[3]
2537 << "S4.=" << &cseed[4]
2538 << "S5.=" << &cseed[5]
2539 << "Label=" << label
2540 << "Freq=" << frequency
2541 << "FitterT.=" << fitterT
2542 << "FitterTC.=" << fitterTC
2547 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2548 if(ntracks == kMaxTracksStack){
2549 AliWarning(Form("Number of seeds reached maximum allowed (%d) in stack.", kMaxTracksStack));
2560 //_____________________________________________________________________________
2561 AliTRDtrackV1* AliTRDtrackerV1::MakeTrack(AliTRDseedV1 *seeds, Double_t *params)
2564 // Build a TRD track out of tracklet candidates
2567 // seeds : array of tracklets
2568 // params : track parameters (see MakeSeeds() function body for a detailed description)
2573 // Detailed description
2575 // To be discussed with Marian !!
2578 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
2579 if (!calibra) AliInfo("Could not get Calibra instance\n");
2581 Double_t alpha = AliTRDgeometry::GetAlpha();
2582 Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
2586 c[ 1] = 0.0; c[ 2] = 2.0;
2587 c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
2588 c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
2589 c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
2591 AliTRDtrackV1 track(seeds, ¶ms[1], c, params[0], params[6]*alpha+shift);
2592 track.PropagateTo(params[0]-5.0);
2593 track.ResetCovariance(1);
2594 Int_t nc = TMath::Abs(FollowBackProlongation(track));
2595 if (nc < 30) return 0x0;
2597 AliTRDtrackV1 *ptrTrack = SetTrack(&track);
2598 ptrTrack->CookLabel(.9);
2599 // computes PID for track
2600 ptrTrack->CookPID();
2601 // update calibration references using this track
2602 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(ptrTrack);
2608 //____________________________________________________________________
2609 Int_t AliTRDtrackerV1::ImproveSeedQuality(AliTRDtrackingChamber **stack, AliTRDseedV1 *cseed)
2612 // Sort tracklets according to "quality" and try to "improve" the first 4 worst
2615 // layers : Array of propagation layers for a stack/supermodule
2616 // cseed : Array of 6 seeding tracklets which has to be improved
2619 // cssed : Improved seeds
2621 // Detailed description
2623 // Iterative procedure in which new clusters are searched for each
2624 // tracklet seed such that the seed quality (see AliTRDseed::GetQuality())
2625 // can be maximized. If some optimization is found the old seeds are replaced.
2630 // make a local working copy
2631 AliTRDtrackingChamber *chamber = 0x0;
2632 AliTRDseedV1 bseed[6];
2634 for (Int_t jLayer = 0; jLayer < 6; jLayer++) bseed[jLayer] = cseed[jLayer];
2636 Float_t lastquality = 10000.0;
2637 Float_t lastchi2 = 10000.0;
2638 Float_t chi2 = 1000.0;
2640 for (Int_t iter = 0; iter < 4; iter++) {
2641 Float_t sumquality = 0.0;
2642 Float_t squality[6];
2643 Int_t sortindexes[6];
2645 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2646 squality[jLayer] = bseed[jLayer].IsOK() ? bseed[jLayer].GetQuality(kTRUE) : 1000.;
2647 sumquality += squality[jLayer];
2649 if ((sumquality >= lastquality) || (chi2 > lastchi2)) break;
2652 lastquality = sumquality;
2654 if (iter > 0) for (Int_t jLayer = 0; jLayer < 6; jLayer++) cseed[jLayer] = bseed[jLayer];
2656 TMath::Sort(6, squality, sortindexes, kFALSE);
2657 for (Int_t jLayer = 5; jLayer > 1; jLayer--) {
2658 Int_t bLayer = sortindexes[jLayer];
2659 if(!(chamber = stack[bLayer])) continue;
2660 bseed[bLayer].AttachClustersIter(chamber, squality[bLayer], kTRUE);
2661 if(bseed[bLayer].IsOK()) nLayers++;
2664 chi2 = FitTiltedRieman(bseed, kTRUE);
2665 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 7){
2666 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2667 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2668 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2669 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2670 cstreamer << "ImproveSeedQuality"
2671 << "EventNumber=" << eventNumber
2672 << "CandidateNumber=" << candidateNumber
2673 << "Iteration=" << iter
2674 << "S0.=" << &bseed[0]
2675 << "S1.=" << &bseed[1]
2676 << "S2.=" << &bseed[2]
2677 << "S3.=" << &bseed[3]
2678 << "S4.=" << &bseed[4]
2679 << "S5.=" << &bseed[5]
2680 << "FitterT.=" << tiltedRieman
2685 // we are sure that at least 2 tracklets are OK !
2689 //_________________________________________________________________________
2690 Double_t AliTRDtrackerV1::CalculateTrackLikelihood(AliTRDseedV1 *tracklets, Double_t *chi2){
2692 // Calculates the Track Likelihood value. This parameter serves as main quality criterion for
2693 // the track selection
2694 // The likelihood value containes:
2695 // - The chi2 values from the both fitters and the chi2 values in z-direction from a linear fit
2696 // - The Sum of the Parameter |slope_ref - slope_fit|/Sigma of the tracklets
2697 // For all Parameters an exponential dependency is used
2699 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
2700 // - Array of chi2 values:
2701 // * Non-Constrained Tilted Riemann fit
2702 // * Vertex-Constrained Tilted Riemann fit
2703 // * z-Direction from Linear fit
2704 // Output: - The calculated track likelihood
2709 Double_t sumdaf = 0, nLayers = 0;
2710 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
2711 if(!tracklets[iLayer].IsOK()) continue;
2712 sumdaf += TMath::Abs((tracklets[iLayer].GetYfit(1) - tracklets[iLayer].GetYref(1))/ tracklets[iLayer].GetSigmaY2());
2715 sumdaf /= Float_t (nLayers - 2.0);
2717 Double_t likeChi2Z = TMath::Exp(-chi2[2] * 0.14); // Chi2Z
2718 Double_t likeChi2TC = (fReconstructor->GetRecoParam() ->IsVertexConstrained()) ?
2719 TMath::Exp(-chi2[1] * 0.677) : 1; // Constrained Tilted Riemann
2720 Double_t likeChi2TR = TMath::Exp(-chi2[0] * 0.78); // Non-constrained Tilted Riemann
2721 Double_t likeAF = TMath::Exp(-sumdaf * 3.23);
2722 Double_t trackLikelihood = likeChi2Z * likeChi2TR * likeAF;
2724 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2725 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2726 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2727 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2728 cstreamer << "CalculateTrackLikelihood0"
2729 << "EventNumber=" << eventNumber
2730 << "CandidateNumber=" << candidateNumber
2731 << "LikeChi2Z=" << likeChi2Z
2732 << "LikeChi2TR=" << likeChi2TR
2733 << "LikeChi2TC=" << likeChi2TC
2734 << "LikeAF=" << likeAF
2735 << "TrackLikelihood=" << trackLikelihood
2739 return trackLikelihood;
2742 //____________________________________________________________________
2743 Double_t AliTRDtrackerV1::CookLikelihood(AliTRDseedV1 *cseed, Int_t planes[4]
2747 // Calculate the probability of this track candidate.
2750 // cseeds : array of candidate tracklets
2751 // planes : array of seeding planes (see seeding configuration)
2752 // chi2 : chi2 values (on the Z and Y direction) from the rieman fit of the track.
2757 // Detailed description
2759 // The track quality is estimated based on the following 4 criteria:
2760 // 1. precision of the rieman fit on the Y direction (likea)
2761 // 2. chi2 on the Y direction (likechi2y)
2762 // 3. chi2 on the Z direction (likechi2z)
2763 // 4. number of attached clusters compared to a reference value
2764 // (see AliTRDrecoParam::fkFindable) (likeN)
2766 // The distributions for each type of probabilities are given below as of
2767 // (date). They have to be checked to assure consistency of estimation.
2770 // ratio of the total number of clusters/track which are expected to be found by the tracker.
2771 Float_t fgFindable = fReconstructor->GetRecoParam() ->GetFindableClusters();
2774 Int_t nclusters = 0;
2775 Double_t sumda = 0.;
2776 for(UChar_t ilayer = 0; ilayer < 4; ilayer++){
2777 Int_t jlayer = planes[ilayer];
2778 nclusters += cseed[jlayer].GetN2();
2779 sumda += TMath::Abs(cseed[jlayer].GetYfitR(1) - cseed[jlayer].GetYref(1));
2781 Double_t likea = TMath::Exp(-sumda*10.6);
2782 Double_t likechi2y = 0.0000000001;
2783 if (chi2[0] < 0.5) likechi2y += TMath::Exp(-TMath::Sqrt(chi2[0]) * 7.73);
2784 Double_t likechi2z = TMath::Exp(-chi2[1] * 0.088) / TMath::Exp(-chi2[1] * 0.019);
2785 Int_t enc = Int_t(fgFindable*4.*fgNTimeBins); // Expected Number Of Clusters, normally 72
2786 Double_t likeN = TMath::Exp(-(enc - nclusters) * 0.19);
2788 Double_t like = likea * likechi2y * likechi2z * likeN;
2790 // 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));
2791 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2792 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2793 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2794 // The Debug Stream contains the seed
2795 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2796 cstreamer << "CookLikelihood"
2797 << "EventNumber=" << eventNumber
2798 << "CandidateNumber=" << candidateNumber
2799 << "tracklet0.=" << &cseed[0]
2800 << "tracklet1.=" << &cseed[1]
2801 << "tracklet2.=" << &cseed[2]
2802 << "tracklet3.=" << &cseed[3]
2803 << "tracklet4.=" << &cseed[4]
2804 << "tracklet5.=" << &cseed[5]
2805 << "sumda=" << sumda
2806 << "chi0=" << chi2[0]
2807 << "chi1=" << chi2[1]
2808 << "likea=" << likea
2809 << "likechi2y=" << likechi2y
2810 << "likechi2z=" << likechi2z
2811 << "nclusters=" << nclusters
2812 << "likeN=" << likeN
2822 //____________________________________________________________________
2823 void AliTRDtrackerV1::GetSeedingConfig(Int_t iconfig, Int_t planes[4])
2826 // Map seeding configurations to detector planes.
2829 // iconfig : configuration index
2830 // planes : member planes of this configuration. On input empty.
2833 // planes : contains the planes which are defining the configuration
2835 // Detailed description
2837 // Here is the list of seeding planes configurations together with
2838 // their topological classification:
2856 // The topologic quality is modeled as follows:
2857 // 1. The general model is define by the equation:
2858 // p(conf) = exp(-conf/2)
2859 // 2. According to the topologic classification, configurations from the same
2860 // class are assigned the agerage value over the model values.
2861 // 3. Quality values are normalized.
2863 // The topologic quality distribution as function of configuration is given below:
2865 // <img src="gif/topologicQA.gif">
2870 case 0: // 5432 TQ 0
2876 case 1: // 4321 TQ 0
2882 case 2: // 3210 TQ 0
2888 case 3: // 5321 TQ 1
2894 case 4: // 4210 TQ 1
2900 case 5: // 5431 TQ 1
2906 case 6: // 4320 TQ 1
2912 case 7: // 5430 TQ 2
2918 case 8: // 5210 TQ 2
2924 case 9: // 5421 TQ 3
2930 case 10: // 4310 TQ 3
2936 case 11: // 5410 TQ 4
2942 case 12: // 5420 TQ 5
2948 case 13: // 5320 TQ 5
2954 case 14: // 5310 TQ 5
2963 //____________________________________________________________________
2964 void AliTRDtrackerV1::GetExtrapolationConfig(Int_t iconfig, Int_t planes[2])
2967 // Returns the extrapolation planes for a seeding configuration.
2970 // iconfig : configuration index
2971 // planes : planes which are not in this configuration. On input empty.
2974 // planes : contains the planes which are not in the configuration
2976 // Detailed description
2980 case 0: // 5432 TQ 0
2984 case 1: // 4321 TQ 0
2988 case 2: // 3210 TQ 0
2992 case 3: // 5321 TQ 1
2996 case 4: // 4210 TQ 1
3000 case 5: // 5431 TQ 1
3004 case 6: // 4320 TQ 1
3008 case 7: // 5430 TQ 2
3012 case 8: // 5210 TQ 2
3016 case 9: // 5421 TQ 3
3020 case 10: // 4310 TQ 3
3024 case 11: // 5410 TQ 4
3028 case 12: // 5420 TQ 5
3032 case 13: // 5320 TQ 5
3036 case 14: // 5310 TQ 5
3043 //____________________________________________________________________
3044 AliCluster* AliTRDtrackerV1::GetCluster(Int_t idx) const
3046 Int_t ncls = fClusters->GetEntriesFast();
3047 return idx >= 0 && idx < ncls ? (AliCluster*)fClusters->UncheckedAt(idx) : 0x0;
3050 //____________________________________________________________________
3051 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(Int_t idx) const
3053 Int_t ntrklt = fTracklets->GetEntriesFast();
3054 return idx >= 0 && idx < ntrklt ? (AliTRDseedV1*)fTracklets->UncheckedAt(idx) : 0x0;
3057 //____________________________________________________________________
3058 AliKalmanTrack* AliTRDtrackerV1::GetTrack(Int_t idx) const
3060 Int_t ntrk = fTracks->GetEntriesFast();
3061 return idx >= 0 && idx < ntrk ? (AliKalmanTrack*)fTracks->UncheckedAt(idx) : 0x0;
3064 //____________________________________________________________________
3065 Float_t AliTRDtrackerV1::CalculateReferenceX(AliTRDseedV1 *tracklets){
3067 // Calculates the reference x-position for the tilted Rieman fit defined as middle
3068 // of the stack (middle between layers 2 and 3). For the calculation all the tracklets
3069 // are taken into account
3071 // Parameters: - Array of tracklets(AliTRDseedV1)
3073 // Output: - The reference x-position(Float_t)
3075 Int_t nDistances = 0;
3076 Float_t meanDistance = 0.;
3077 Int_t startIndex = 5;
3078 for(Int_t il =5; il > 0; il--){
3079 if(tracklets[il].IsOK() && tracklets[il -1].IsOK()){
3080 Float_t xdiff = tracklets[il].GetX0() - tracklets[il -1].GetX0();
3081 meanDistance += xdiff;
3084 if(tracklets[il].IsOK()) startIndex = il;
3086 if(tracklets[0].IsOK()) startIndex = 0;
3088 // We should normally never get here
3089 Float_t xpos[2]; memset(xpos, 0, sizeof(Float_t) * 2);
3090 Int_t iok = 0, idiff = 0;
3091 // This attempt is worse and should be avoided:
3092 // check for two chambers which are OK and repeat this without taking the mean value
3093 // Strategy avoids a division by 0;
3094 for(Int_t il = 5; il >= 0; il--){
3095 if(tracklets[il].IsOK()){
3096 xpos[iok] = tracklets[il].GetX0();
3100 if(iok) idiff++; // to get the right difference;
3104 meanDistance = (xpos[0] - xpos[1])/idiff;
3107 // we have do not even have 2 layers which are OK? The we do not need to fit at all
3112 meanDistance /= nDistances;
3114 return tracklets[startIndex].GetX0() + (2.5 - startIndex) * meanDistance - 0.5 * (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
3117 //_____________________________________________________________________________
3118 Int_t AliTRDtrackerV1::Freq(Int_t n, const Int_t *inlist
3119 , Int_t *outlist, Bool_t down)
3122 // Sort eleements according occurancy
3123 // The size of output array has is 2*n
3130 Int_t *sindexS = new Int_t[n]; // Temporary array for sorting
3131 Int_t *sindexF = new Int_t[2*n];
3132 for (Int_t i = 0; i < n; i++) {
3136 TMath::Sort(n,inlist,sindexS,down);
3138 Int_t last = inlist[sindexS[0]];
3141 sindexF[0+n] = last;
3145 for (Int_t i = 1; i < n; i++) {
3146 val = inlist[sindexS[i]];
3148 sindexF[countPos]++;
3152 sindexF[countPos+n] = val;
3153 sindexF[countPos]++;
3161 // Sort according frequency
3162 TMath::Sort(countPos,sindexF,sindexS,kTRUE);
3164 for (Int_t i = 0; i < countPos; i++) {
3165 outlist[2*i ] = sindexF[sindexS[i]+n];
3166 outlist[2*i+1] = sindexF[sindexS[i]];
3177 //____________________________________________________________________
3178 void AliTRDtrackerV1::SetReconstructor(const AliTRDReconstructor *rec)
3180 fReconstructor = rec;
3181 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
3182 if(!fgDebugStreamer){
3183 TDirectory *savedir = gDirectory;
3184 fgDebugStreamer = new TTreeSRedirector("TRD.TrackerDebug.root");
3190 //_____________________________________________________________________________
3191 Float_t AliTRDtrackerV1::GetChi2Y(AliTRDseedV1 *tracklets) const
3193 // Chi2 definition on y-direction
3196 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
3197 if(!tracklets[ipl].IsOK()) continue;
3198 Double_t distLayer = tracklets[ipl].GetYfit(0) - tracklets[ipl].GetYref(0);
3199 chi2 += distLayer * distLayer;
3204 //____________________________________________________________________
3205 void AliTRDtrackerV1::ResetSeedTB()
3207 // reset buffer for seeding time bin layers. If the time bin
3208 // layers are not allocated this function allocates them
3210 for(Int_t isl=0; isl<kNSeedPlanes; isl++){
3211 if(!fSeedTB[isl]) fSeedTB[isl] = new AliTRDchamberTimeBin();
3212 else fSeedTB[isl]->Clear();
3216 //_____________________________________________________________________________
3217 Float_t AliTRDtrackerV1::GetChi2Z(AliTRDseedV1 *tracklets) const
3219 // Chi2 definition on z-direction
3222 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
3223 if(!tracklets[ipl].IsOK()) continue;
3224 Double_t distLayer = tracklets[ipl].GetMeanz() - tracklets[ipl].GetZref(0);
3225 chi2 += distLayer * distLayer;
3230 ///////////////////////////////////////////////////////
3232 // Resources of class AliTRDLeastSquare //
3234 ///////////////////////////////////////////////////////
3236 //_____________________________________________________________________________
3237 AliTRDtrackerV1::AliTRDLeastSquare::AliTRDLeastSquare(){
3239 // Constructor of the nested class AliTRDtrackFitterLeastSquare
3241 memset(fParams, 0, sizeof(Double_t) * 2);
3242 memset(fSums, 0, sizeof(Double_t) * 5);
3243 memset(fCovarianceMatrix, 0, sizeof(Double_t) * 3);
3247 //_____________________________________________________________________________
3248 void AliTRDtrackerV1::AliTRDLeastSquare::AddPoint(Double_t *x, Double_t y, Double_t sigmaY){
3250 // Adding Point to the fitter
3252 Double_t weight = 1/(sigmaY * sigmaY);
3254 // printf("Adding point x = %f, y = %f, sigma = %f\n", xpt, y, sigmaY);
3256 fSums[1] += weight * xpt;
3257 fSums[2] += weight * y;
3258 fSums[3] += weight * xpt * y;
3259 fSums[4] += weight * xpt * xpt;
3260 fSums[5] += weight * y * y;
3263 //_____________________________________________________________________________
3264 void AliTRDtrackerV1::AliTRDLeastSquare::RemovePoint(Double_t *x, Double_t y, Double_t sigmaY){
3266 // Remove Point from the sample
3268 Double_t weight = 1/(sigmaY * sigmaY);
3271 fSums[1] -= weight * xpt;
3272 fSums[2] -= weight * y;
3273 fSums[3] -= weight * xpt * y;
3274 fSums[4] -= weight * xpt * xpt;
3275 fSums[5] -= weight * y * y;
3278 //_____________________________________________________________________________
3279 void AliTRDtrackerV1::AliTRDLeastSquare::Eval(){
3281 // Evaluation of the fit:
3282 // Calculation of the parameters
3283 // Calculation of the covariance matrix
3286 Double_t denominator = fSums[0] * fSums[4] - fSums[1] *fSums[1];
3287 if(denominator==0) return;
3289 // for(Int_t isum = 0; isum < 5; isum++)
3290 // printf("fSums[%d] = %f\n", isum, fSums[isum]);
3291 // printf("denominator = %f\n", denominator);
3292 fParams[0] = (fSums[2] * fSums[4] - fSums[1] * fSums[3])/ denominator;
3293 fParams[1] = (fSums[0] * fSums[3] - fSums[1] * fSums[2]) / denominator;
3294 // printf("fParams[0] = %f, fParams[1] = %f\n", fParams[0], fParams[1]);
3296 // Covariance matrix
3297 fCovarianceMatrix[0] = fSums[4] - fSums[1] * fSums[1] / fSums[0];
3298 fCovarianceMatrix[1] = fSums[5] - fSums[2] * fSums[2] / fSums[0];
3299 fCovarianceMatrix[2] = fSums[3] - fSums[1] * fSums[2] / fSums[0];
3302 //_____________________________________________________________________________
3303 Double_t AliTRDtrackerV1::AliTRDLeastSquare::GetFunctionValue(Double_t *xpos) const {
3305 // Returns the Function value of the fitted function at a given x-position
3307 return fParams[0] + fParams[1] * (*xpos);
3310 //_____________________________________________________________________________
3311 void AliTRDtrackerV1::AliTRDLeastSquare::GetCovarianceMatrix(Double_t *storage) const {
3313 // Copies the values of the covariance matrix into the storage
3315 memcpy(storage, fCovarianceMatrix, sizeof(Double_t) * 3);