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
108 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
109 TDirectory *savedir = gDirectory;
110 fgDebugStreamer = new TTreeSRedirector("TRD.TrackerDebug.root");
116 //____________________________________________________________________
117 AliTRDtrackerV1::~AliTRDtrackerV1()
123 if(fgDebugStreamer) delete fgDebugStreamer;
124 if(fgRieman) delete fgRieman;
125 if(fgTiltedRieman) delete fgTiltedRieman;
126 if(fgTiltedRiemanConstrained) delete fgTiltedRiemanConstrained;
127 for(Int_t isl =0; isl<kNSeedPlanes; isl++) if(fSeedTB[isl]) delete fSeedTB[isl];
128 if(fTracks) {fTracks->Delete(); delete fTracks;}
129 if(fTracklets) {fTracklets->Delete(); delete fTracklets;}
130 if(fClusters) {fClusters->Delete(); delete fClusters;}
131 if(fGeom) delete fGeom;
134 //____________________________________________________________________
135 Int_t AliTRDtrackerV1::Clusters2Tracks(AliESDEvent *esd)
138 // Steering stand alone tracking for full TRD detector
141 // esd : The ESD event. On output it contains
142 // the ESD tracks found in TRD.
145 // Number of tracks found in the TRD detector.
147 // Detailed description
148 // 1. Launch individual SM trackers.
149 // See AliTRDtrackerV1::Clusters2TracksSM() for details.
152 if(!fReconstructor->GetRecoParam() ){
153 AliError("Reconstruction configuration not initialized. Call first AliTRDReconstructor::SetRecoParam().");
157 //AliInfo("Start Track Finder ...");
159 for(int ism=0; ism<AliTRDgeometry::kNsector; ism++){
160 // for(int ism=1; ism<2; ism++){
161 //AliInfo(Form("Processing supermodule %i ...", ism));
162 ntracks += Clusters2TracksSM(ism, esd);
164 AliInfo(Form("Number of found tracks : %d", ntracks));
169 //_____________________________________________________________________________
170 Bool_t AliTRDtrackerV1::GetTrackPoint(Int_t index, AliTrackPoint &p) const
172 //AliInfo(Form("Asking for tracklet %d", index));
174 AliTRDseedV1 *tracklet = GetTracklet(index);
175 if (!tracklet) return kFALSE;
177 // get detector for this tracklet
178 AliTRDcluster *cl = 0x0;
179 Int_t ic = 0; do; while(!(cl = tracklet->GetClusters(ic++)));
180 Int_t idet = cl->GetDetector();
183 local[0] = tracklet->GetX0();
184 local[1] = tracklet->GetYfit(0);
185 local[2] = tracklet->GetZfit(0);
187 fGeom->RotateBack(idet, local, global);
188 p.SetXYZ(global[0],global[1],global[2]);
192 AliGeomManager::ELayerID iLayer = AliGeomManager::kTRD1;
193 switch (fGeom->GetLayer(idet)) {
195 iLayer = AliGeomManager::kTRD1;
198 iLayer = AliGeomManager::kTRD2;
201 iLayer = AliGeomManager::kTRD3;
204 iLayer = AliGeomManager::kTRD4;
207 iLayer = AliGeomManager::kTRD5;
210 iLayer = AliGeomManager::kTRD6;
213 Int_t modId = fGeom->GetSector(idet) * fGeom->Nstack() + fGeom->GetStack(idet);
214 UShort_t volid = AliGeomManager::LayerToVolUID(iLayer, modId);
215 p.SetVolumeID(volid);
220 //____________________________________________________________________
221 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitter()
223 if(!fgTiltedRieman) fgTiltedRieman = new TLinearFitter(4, "hyp4");
224 return fgTiltedRieman;
227 //____________________________________________________________________
228 TLinearFitter* AliTRDtrackerV1::GetTiltedRiemanFitterConstraint()
230 if(!fgTiltedRiemanConstrained) fgTiltedRiemanConstrained = new TLinearFitter(2, "hyp2");
231 return fgTiltedRiemanConstrained;
234 //____________________________________________________________________
235 AliRieman* AliTRDtrackerV1::GetRiemanFitter()
237 if(!fgRieman) fgRieman = new AliRieman(AliTRDtrackingChamber::kNTimeBins * AliTRDgeometry::kNlayer);
241 //_____________________________________________________________________________
242 Int_t AliTRDtrackerV1::PropagateBack(AliESDEvent *event)
245 // Gets seeds from ESD event. The seeds are AliTPCtrack's found and
246 // backpropagated by the TPC tracker. Each seed is first propagated
247 // to the TRD, and then its prolongation is searched in the TRD.
248 // If sufficiently long continuation of the track is found in the TRD
249 // the track is updated, otherwise it's stored as originaly defined
250 // by the TPC tracker.
253 // Calibration monitor
254 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
255 if (!calibra) AliInfo("Could not get Calibra instance\n");
257 Int_t found = 0; // number of tracks found
258 Float_t foundMin = 20.0;
260 Float_t *quality = 0x0;
262 Int_t nSeed = event->GetNumberOfTracks();
264 quality = new Float_t[nSeed];
265 index = new Int_t[nSeed];
266 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
267 AliESDtrack *seed = event->GetTrack(iSeed);
268 Double_t covariance[15];
269 seed->GetExternalCovariance(covariance);
270 quality[iSeed] = covariance[0] + covariance[2];
272 // Sort tracks according to covariance of local Y and Z
273 TMath::Sort(nSeed,quality,index,kFALSE);
276 // Backpropagate all seeds
279 for (Int_t iSeed = 0; iSeed < nSeed; iSeed++) {
281 // Get the seeds in sorted sequence
282 AliESDtrack *seed = event->GetTrack(index[iSeed]);
284 // Check the seed status
285 ULong_t status = seed->GetStatus();
286 if ((status & AliESDtrack::kTPCout) == 0) continue;
287 if ((status & AliESDtrack::kTRDout) != 0) continue;
289 // Do the back prolongation
290 new(&track) AliTRDtrackV1(*seed);
292 //Int_t lbl = seed->GetLabel();
293 //track.SetSeedLabel(lbl);
294 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup); // Make backup
295 Float_t p4 = track.GetC();
296 expectedClr = FollowBackProlongation(track);
298 if (expectedClr<0) continue; // Back prolongation failed
302 // computes PID for track
304 // update calibration references using this track
305 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(&track);
306 // save calibration object
307 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
308 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
310 track.UpdateESDtrack(seed);
312 // Add TRD track to ESDfriendTrack
313 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0 /*&& quality TODO*/){
314 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(track);
315 calibTrack->SetOwner();
316 seed->AddCalibObject(calibTrack);
321 if ((TMath::Abs(track.GetC() - p4) / TMath::Abs(p4) < 0.2) ||(track.Pt() > 0.8)) {
323 // Make backup for back propagation
325 Int_t foundClr = track.GetNumberOfClusters();
326 if (foundClr >= foundMin) {
327 //AliInfo(Form("Making backup track ncls [%d]...", foundClr));
329 //track.CookdEdxTimBin(seed->GetID());
330 track.CookLabel(1. - fgkLabelFraction);
331 if(track.GetBackupTrack()) UseClusters(track.GetBackupTrack());
334 // Sign only gold tracks
335 if (track.GetChi2() / track.GetNumberOfClusters() < 4) {
336 if ((seed->GetKinkIndex(0) == 0) && (track.Pt() < 1.5)) UseClusters(&track);
338 Bool_t isGold = kFALSE;
341 if (track.GetChi2() / track.GetNumberOfClusters() < 5) {
342 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
348 if ((!isGold) && (track.GetNCross() == 0) && (track.GetChi2() / track.GetNumberOfClusters() < 7)) {
349 //seed->UpdateTrackParams(track, AliESDtrack::kTRDbackup);
350 if (track.GetBackupTrack()) seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
355 if ((!isGold) && (track.GetBackupTrack())) {
356 if ((track.GetBackupTrack()->GetNumberOfClusters() > foundMin) && ((track.GetBackupTrack()->GetChi2()/(track.GetBackupTrack()->GetNumberOfClusters()+1)) < 7)) {
357 seed->UpdateTrackParams(track.GetBackupTrack(),AliESDtrack::kTRDbackup);
362 //if ((track->StatusForTOF() > 0) && (track->GetNCross() == 0) && (Float_t(track->GetNumberOfClusters()) / Float_t(track->GetNExpected()) > 0.4)) {
363 //seed->UpdateTrackParams(track->GetBackupTrack(), AliESDtrack::kTRDbackup);
368 // Propagation to the TOF (I.Belikov)
369 if (track.IsStopped() == kFALSE) {
370 Double_t xtof = 371.0;
371 Double_t xTOF0 = 370.0;
373 Double_t c2 = track.GetSnp() + track.GetC() * (xtof - track.GetX());
374 if (TMath::Abs(c2) >= 0.99) continue;
376 if (!PropagateToX(track, xTOF0, fgkMaxStep)) continue;
378 // Energy losses taken to the account - check one more time
379 c2 = track.GetSnp() + track.GetC() * (xtof - track.GetX());
380 if (TMath::Abs(c2) >= 0.99) continue;
382 //if (!PropagateToX(*track,xTOF0,fgkMaxStep)) {
383 // fHBackfit->Fill(7);
388 Double_t ymax = xtof * TMath::Tan(0.5 * AliTRDgeometry::GetAlpha());
390 track.GetYAt(xtof,GetBz(),y);
392 if (!track.Rotate( AliTRDgeometry::GetAlpha())) continue;
393 }else if (y < -ymax) {
394 if (!track.Rotate(-AliTRDgeometry::GetAlpha())) continue;
397 if (track.PropagateTo(xtof)) {
398 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
399 track.UpdateESDtrack(seed);
402 if ((track.GetNumberOfClusters() > 15) && (track.GetNumberOfClusters() > 0.5*expectedClr)) {
403 seed->UpdateTrackParams(&track, AliESDtrack::kTRDout);
405 track.UpdateESDtrack(seed);
409 seed->SetTRDQuality(track.StatusForTOF());
410 seed->SetTRDBudget(track.GetBudget(0));
412 if(index) delete [] index;
413 if(quality) delete [] quality;
416 AliInfo(Form("Number of seeds: %d", nSeed));
417 AliInfo(Form("Number of back propagated TRD tracks: %d", found));
419 // run stand alone tracking
420 if (fReconstructor->IsSeeding()) Clusters2Tracks(event);
426 //____________________________________________________________________
427 Int_t AliTRDtrackerV1::RefitInward(AliESDEvent *event)
430 // Refits tracks within the TRD. The ESD event is expected to contain seeds
431 // at the outer part of the TRD.
432 // The tracks are propagated to the innermost time bin
433 // of the TRD and the ESD event is updated
434 // Origin: Thomas KUHR (Thomas.Kuhr@cern.ch)
437 Int_t nseed = 0; // contor for loaded seeds
438 Int_t found = 0; // contor for updated TRD tracks
442 for (Int_t itrack = 0; itrack < event->GetNumberOfTracks(); itrack++) {
443 AliESDtrack *seed = event->GetTrack(itrack);
444 new(&track) AliTRDtrackV1(*seed);
446 if (track.GetX() < 270.0) {
447 seed->UpdateTrackParams(&track, AliESDtrack::kTRDbackup);
451 ULong_t status = seed->GetStatus();
452 if((status & AliESDtrack::kTRDout) == 0) continue;
453 if((status & AliESDtrack::kTRDin) != 0) continue;
456 track.ResetCovariance(50.0);
458 // do the propagation and processing
459 Bool_t kUPDATE = kFALSE;
460 Double_t xTPC = 250.0;
461 if(FollowProlongation(track)){
463 if (PropagateToX(track, xTPC, fgkMaxStep)) { // -with update
464 seed->UpdateTrackParams(&track, AliESDtrack::kTRDrefit);
470 // Prolongate to TPC without update
472 AliTRDtrackV1 tt(*seed);
473 if (PropagateToX(tt, xTPC, fgkMaxStep)) seed->UpdateTrackParams(&tt, AliESDtrack::kTRDrefit);
476 AliInfo(Form("Number of loaded seeds: %d",nseed));
477 AliInfo(Form("Number of found tracks from loaded seeds: %d",found));
482 //____________________________________________________________________
483 Int_t AliTRDtrackerV1::FollowProlongation(AliTRDtrackV1 &t)
485 // Extrapolates the TRD track in the TPC direction.
488 // t : the TRD track which has to be extrapolated
491 // number of clusters attached to the track
493 // Detailed description
495 // Starting from current radial position of track <t> this function
496 // extrapolates the track through the 6 TRD layers. The following steps
497 // are being performed for each plane:
499 // a. get plane limits in the local x direction
500 // b. check crossing sectors
501 // c. check track inclination
502 // 2. search tracklet in the tracker list (see GetTracklet() for details)
503 // 3. evaluate material budget using the geo manager
504 // 4. propagate and update track using the tracklet information.
509 Int_t nClustersExpected = 0;
510 Int_t lastplane = 5; //GetLastPlane(&t);
511 for (Int_t iplane = lastplane; iplane >= 0; iplane--) {
513 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
514 if(!tracklet) continue;
515 if(!tracklet->IsOK()) AliWarning("tracklet not OK");
517 Double_t x = tracklet->GetX0();
518 // reject tracklets which are not considered for inward refit
519 if(x > t.GetX()+fgkMaxStep) continue;
521 // append tracklet to track
522 t.SetTracklet(tracklet, index);
524 if (x < (t.GetX()-fgkMaxStep) && !PropagateToX(t, x+fgkMaxStep, fgkMaxStep)) break;
525 if (!AdjustSector(&t)) break;
527 // Start global position
531 // End global position
532 Double_t alpha = t.GetAlpha(), y, z;
533 if (!t.GetProlongation(x,y,z)) break;
535 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
536 xyz1[1] = x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
539 // Get material budget
541 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
542 Double_t xrho= param[0]*param[4];
543 Double_t xx0 = param[1]; // Get mean propagation parameters
545 // Propagate and update
546 t.PropagateTo(x, xx0, xrho);
547 if (!AdjustSector(&t)) break;
549 Double_t maxChi2 = t.GetPredictedChi2(tracklet);
550 if (maxChi2 < 1e+10 && t.Update(tracklet, maxChi2)){
551 nClustersExpected += tracklet->GetN();
555 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
557 for(int iplane=0; iplane<6; iplane++){
558 AliTRDseedV1 *tracklet = GetTracklet(&t, iplane, index);
559 if(!tracklet) continue;
560 t.SetTracklet(tracklet, index);
563 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
564 TTreeSRedirector &cstreamer = *fgDebugStreamer;
565 cstreamer << "FollowProlongation"
566 << "EventNumber=" << eventNumber
567 << "ncl=" << nClustersExpected
572 return nClustersExpected;
576 //_____________________________________________________________________________
577 Int_t AliTRDtrackerV1::FollowBackProlongation(AliTRDtrackV1 &t)
579 // Extrapolates the TRD track in the TOF direction.
582 // t : the TRD track which has to be extrapolated
585 // number of clusters attached to the track
587 // Detailed description
589 // Starting from current radial position of track <t> this function
590 // extrapolates the track through the 6 TRD layers. The following steps
591 // are being performed for each plane:
593 // a. get plane limits in the local x direction
594 // b. check crossing sectors
595 // c. check track inclination
596 // 2. build tracklet (see AliTRDseed::AttachClusters() for details)
597 // 3. evaluate material budget using the geo manager
598 // 4. propagate and update track using the tracklet information.
603 Int_t nClustersExpected = 0;
604 Double_t clength = AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick();
605 AliTRDtrackingChamber *chamber = 0x0;
607 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
609 // Loop through the TRD layers
610 for (Int_t ilayer = 0; ilayer < AliTRDgeometry::Nlayer(); ilayer++) {
611 // BUILD TRACKLET IF NOT ALREADY BUILT
612 Double_t x = 0., y, z, alpha;
613 ptrTracklet = t.GetTracklet(ilayer);
615 ptrTracklet = new(&tracklet) AliTRDseedV1(ilayer);
616 ptrTracklet->SetReconstructor(fReconstructor);
617 alpha = t.GetAlpha();
618 Int_t sector = Int_t(alpha/AliTRDgeometry::GetAlpha() + (alpha>0. ? 0 : AliTRDgeometry::kNsector));
620 if(!fTrSec[sector].GetNChambers()) continue;
622 if((x = fTrSec[sector].GetX(ilayer)) < 1.) continue;
624 if (!t.GetProlongation(x, y, z)) return -nClustersExpected;
625 Int_t stack = fGeom->GetStack(z, ilayer);
626 Int_t nCandidates = stack >= 0 ? 1 : 2;
627 z -= stack >= 0 ? 0. : 4.;
629 for(int icham=0; icham<nCandidates; icham++, z+=8){
630 if((stack = fGeom->GetStack(z, ilayer)) < 0) continue;
632 if(!(chamber = fTrSec[sector].GetChamber(stack, ilayer))) continue;
634 if(chamber->GetNClusters() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
638 AliTRDpadPlane *pp = fGeom->GetPadPlane(ilayer, stack);
639 tracklet.SetTilt(TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle()));
640 tracklet.SetPadLength(pp->GetLengthIPad());
641 tracklet.SetPlane(ilayer);
643 if(!tracklet.Init(&t)){
645 return nClustersExpected;
647 if(!tracklet.AttachClustersIter(chamber, 1000.)) continue;
650 if(tracklet.GetN() < fgNTimeBins*fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
655 if(!ptrTracklet->IsOK()){
656 if(x < 1.) continue; //temporary
657 if(!PropagateToX(t, x-fgkMaxStep, fgkMaxStep)) return -nClustersExpected;
658 if(!AdjustSector(&t)) return -nClustersExpected;
659 if(TMath::Abs(t.GetSnp()) > fgkMaxSnp) return -nClustersExpected;
663 // Propagate closer to the current chamber if neccessary
665 if (x > (fgkMaxStep + t.GetX()) && !PropagateToX(t, x-fgkMaxStep, fgkMaxStep)) return -nClustersExpected;
666 if (!AdjustSector(&t)) return -nClustersExpected;
667 if (TMath::Abs(t.GetSnp()) > fgkMaxSnp) return -nClustersExpected;
669 // load tracklet to the tracker and the track
670 ptrTracklet = SetTracklet(ptrTracklet);
671 t.SetTracklet(ptrTracklet, fTracklets->GetEntriesFast()-1);
674 // Calculate the mean material budget along the path inside the chamber
675 //Calculate global entry and exit positions of the track in chamber (only track prolongation)
676 Double_t xyz0[3]; // entry point
678 alpha = t.GetAlpha();
679 x = ptrTracklet->GetX0();
680 if (!t.GetProlongation(x, y, z)) return -nClustersExpected;
681 Double_t xyz1[3]; // exit point
682 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
683 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
686 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
687 // The mean propagation parameters
688 Double_t xrho = param[0]*param[4]; // density*length
689 Double_t xx0 = param[1]; // radiation length
691 // Propagate and update track
692 if (!t.PropagateTo(x, xx0, xrho)) return -nClustersExpected;
693 if (!AdjustSector(&t)) return -nClustersExpected;
694 Double_t maxChi2 = t.GetPredictedChi2(ptrTracklet);
695 if (!t.Update(ptrTracklet, maxChi2)) return -nClustersExpected;
697 nClustersExpected += ptrTracklet->GetN();
698 //t.SetTracklet(&tracklet, index);
700 // Reset material budget if 2 consecutive gold
701 if(ilayer>0 && t.GetTracklet(ilayer-1) && ptrTracklet->GetN() + t.GetTracklet(ilayer-1)->GetN() > 20) t.SetBudget(2, 0.);
703 // Make backup of the track until is gold
704 // TO DO update quality check of the track.
705 // consider comparison with fTimeBinsRange
706 Float_t ratio0 = ptrTracklet->GetN() / Float_t(fgNTimeBins);
707 //Float_t ratio1 = Float_t(t.GetNumberOfClusters()+1) / Float_t(t.GetNExpected()+1);
708 //printf("tracklet.GetChi2() %f [< 18.0]\n", tracklet.GetChi2());
709 //printf("ratio0 %f [> 0.8]\n", ratio0);
710 //printf("ratio1 %f [> 0.6]\n", ratio1);
711 //printf("ratio0+ratio1 %f [> 1.5]\n", ratio0+ratio1);
712 //printf("t.GetNCross() %d [== 0]\n", t.GetNCross());
713 //printf("TMath::Abs(t.GetSnp()) %f [< 0.85]\n", TMath::Abs(t.GetSnp()));
714 //printf("t.GetNumberOfClusters() %d [> 20]\n", t.GetNumberOfClusters());
716 if (//(tracklet.GetChi2() < 18.0) && TO DO check with FindClusters and move it to AliTRDseed::Update
719 //(ratio0+ratio1 > 1.5) &&
720 (t.GetNCross() == 0) &&
721 (TMath::Abs(t.GetSnp()) < 0.85) &&
722 (t.GetNumberOfClusters() > 20)) t.MakeBackupTrack();
726 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
727 TTreeSRedirector &cstreamer = *fgDebugStreamer;
728 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
729 //AliTRDtrackV1 *debugTrack = new AliTRDtrackV1(t);
730 //debugTrack->SetOwner();
731 cstreamer << "FollowBackProlongation"
732 << "EventNumber=" << eventNumber
733 << "ncl=" << nClustersExpected
734 //<< "track.=" << debugTrack
738 return nClustersExpected;
741 //_________________________________________________________________________
742 Float_t AliTRDtrackerV1::FitRieman(AliTRDseedV1 *tracklets, Double_t *chi2, Int_t *planes){
744 // Fits a Riemann-circle to the given points without tilting pad correction.
745 // The fit is performed using an instance of the class AliRieman (equations
746 // and transformations see documentation of this class)
747 // Afterwards all the tracklets are Updated
749 // Parameters: - Array of tracklets (AliTRDseedV1)
750 // - Storage for the chi2 values (beginning with direction z)
751 // - Seeding configuration
752 // Output: - The curvature
754 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
756 Int_t allplanes[] = {0, 1, 2, 3, 4, 5};
757 Int_t *ppl = &allplanes[0];
763 for(Int_t il = 0; il < maxLayers; il++){
764 if(!tracklets[ppl[il]].IsOK()) continue;
765 fitter->AddPoint(tracklets[ppl[il]].GetX0(), tracklets[ppl[il]].GetYfitR(0), tracklets[ppl[il]].GetZProb(),1,10);
768 // Set the reference position of the fit and calculate the chi2 values
769 memset(chi2, 0, sizeof(Double_t) * 2);
770 for(Int_t il = 0; il < maxLayers; il++){
771 // Reference positions
772 tracklets[ppl[il]].Init(fitter);
775 if((!tracklets[ppl[il]].IsOK()) && (!planes)) continue;
776 chi2[0] += tracklets[ppl[il]].GetChi2Y();
777 chi2[1] += tracklets[ppl[il]].GetChi2Z();
779 return fitter->GetC();
782 //_________________________________________________________________________
783 void AliTRDtrackerV1::FitRieman(AliTRDcluster **seedcl, Double_t chi2[2])
786 // Performs a Riemann helix fit using the seedclusters as spacepoints
787 // Afterwards the chi2 values are calculated and the seeds are updated
789 // Parameters: - The four seedclusters
790 // - The tracklet array (AliTRDseedV1)
791 // - The seeding configuration
796 AliRieman *fitter = AliTRDtrackerV1::GetRiemanFitter();
798 for(Int_t i = 0; i < 4; i++)
799 fitter->AddPoint(seedcl[i]->GetX(), seedcl[i]->GetY(), seedcl[i]->GetZ(), 1, 10);
803 // Update the seed and calculated the chi2 value
804 chi2[0] = 0; chi2[1] = 0;
805 for(Int_t ipl = 0; ipl < kNSeedPlanes; ipl++){
807 chi2[0] += (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetZ() - fitter->GetZat(seedcl[ipl]->GetX()));
808 chi2[1] += (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX())) * (seedcl[ipl]->GetY() - fitter->GetYat(seedcl[ipl]->GetX()));
813 //_________________________________________________________________________
814 Float_t AliTRDtrackerV1::FitTiltedRiemanConstraint(AliTRDseedV1 *tracklets, Double_t zVertex)
817 // Fits a helix to the clusters. Pad tilting is considered. As constraint it is
818 // assumed that the vertex position is set to 0.
819 // This method is very usefull for high-pt particles
820 // Basis for the fit: (x - x0)^2 + (y - y0)^2 - R^2 = 0
821 // x0, y0: Center of the circle
822 // Measured y-position: ymeas = y - tan(phiT)(zc - zt)
823 // zc: center of the pad row
824 // Equation which has to be fitted (after transformation):
825 // a + b * u + e * v + 2*(ymeas + tan(phiT)(z - zVertex))*t = 0
829 // v = 2 * x * tan(phiT) * t
830 // Parameters in the equation:
831 // a = -1/y0, b = x0/y0, e = dz/dx
833 // The Curvature is calculated by the following equation:
834 // - curv = a/Sqrt(b^2 + 1) = 1/R
835 // Parameters: - the 6 tracklets
836 // - the Vertex constraint
837 // Output: - the Chi2 value of the track
842 TLinearFitter *fitter = GetTiltedRiemanFitterConstraint();
843 fitter->StoreData(kTRUE);
844 fitter->ClearPoints();
845 AliTRDcluster *cl = 0x0;
847 Float_t x, y, z, w, t, error, tilt;
850 for(Int_t ilr = 0; ilr < AliTRDgeometry::kNlayer; ilr++){
851 if(!tracklets[ilr].IsOK()) continue;
852 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
853 if(!tracklets[ilr].IsUsable(itb)) continue;
854 cl = tracklets[ilr].GetClusters(itb);
858 tilt = tracklets[ilr].GetTilt();
860 t = 1./(x * x + y * y);
862 uvt[1] = 2. * x * t * tilt ;
863 w = 2. * (y + tilt * (z - zVertex)) * t;
864 error = 2. * 0.2 * t;
865 fitter->AddPoint(uvt, w, error);
871 // Calculate curvature
872 Double_t a = fitter->GetParameter(0);
873 Double_t b = fitter->GetParameter(1);
874 Double_t curvature = a/TMath::Sqrt(b*b + 1);
876 Float_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
877 for(Int_t ip = 0; ip < AliTRDtrackerV1::kNPlanes; ip++)
878 tracklets[ip].SetCC(curvature);
880 /* if(fReconstructor->GetStreamLevel() >= 5){
881 //Linear Model on z-direction
882 Double_t xref = CalculateReferenceX(tracklets); // Relative to the middle of the stack
883 Double_t slope = fitter->GetParameter(2);
884 Double_t zref = slope * xref;
885 Float_t chi2Z = CalculateChi2Z(tracklets, zref, slope, xref);
886 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
887 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
888 TTreeSRedirector &treeStreamer = *fgDebugStreamer;
889 treeStreamer << "FitTiltedRiemanConstraint"
890 << "EventNumber=" << eventNumber
891 << "CandidateNumber=" << candidateNumber
892 << "Curvature=" << curvature
893 << "Chi2Track=" << chi2track
901 //_________________________________________________________________________
902 Float_t AliTRDtrackerV1::FitTiltedRieman(AliTRDseedV1 *tracklets, Bool_t sigError)
905 // Performs a Riemann fit taking tilting pad correction into account
906 // The equation of a Riemann circle, where the y position is substituted by the
907 // measured y-position taking pad tilting into account, has to be transformed
908 // into a 4-dimensional hyperplane equation
909 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
910 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
911 // zc: center of the pad row
912 // zt: z-position of the track
913 // The z-position of the track is assumed to be linear dependent on the x-position
914 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
915 // Transformation: u = 2 * x * t
916 // v = 2 * tan(phiT) * t
917 // w = 2 * tan(phiT) * (x - xref) * t
918 // t = 1 / (x^2 + ymeas^2)
919 // Parameters: a = -1/y0
921 // c = (R^2 -x0^2 - y0^2)/y0
924 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
925 // results from the simple riemann fit. Afterwards the fit is redone.
926 // The curvature is calculated according to the formula:
927 // curv = a/(1 + b^2 + c*a) = 1/R
929 // Paramters: - Array of tracklets (connected to the track candidate)
930 // - Flag selecting the error definition
931 // Output: - Chi2 values of the track (in Parameter list)
933 TLinearFitter *fitter = GetTiltedRiemanFitter();
934 fitter->StoreData(kTRUE);
935 fitter->ClearPoints();
936 AliTRDLeastSquare zfitter;
937 AliTRDcluster *cl = 0x0;
939 Double_t xref = CalculateReferenceX(tracklets);
940 Double_t x, y, z, t, tilt, dx, w, we;
943 // Containers for Least-square fitter
944 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
945 if(!tracklets[ipl].IsOK()) continue;
946 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
947 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
948 if (!tracklets[ipl].IsUsable(itb)) continue;
952 tilt = tracklets[ipl].GetTilt();
958 uvt[2] = 2. * tilt * t;
959 uvt[3] = 2. * tilt * dx * t;
960 w = 2. * (y + tilt*z) * t;
961 // error definition changes for the different calls
963 we *= sigError ? tracklets[ipl].GetSigmaY() : 0.2;
964 fitter->AddPoint(uvt, w, we);
965 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
972 Double_t offset = fitter->GetParameter(3);
973 Double_t slope = fitter->GetParameter(4);
975 // Linear fitter - not possible to make boundaries
976 // Do not accept non possible z and dzdx combinations
977 Bool_t acceptablez = kTRUE;
979 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
980 if(!tracklets[iLayer].IsOK()) continue;
981 zref = offset + slope * (tracklets[iLayer].GetX0() - xref);
982 if (TMath::Abs(tracklets[iLayer].GetZProb() - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
983 acceptablez = kFALSE;
986 Double_t dzmf = zfitter.GetFunctionParameter(1);
987 Double_t zmf = zfitter.GetFunctionValue(&xref);
988 fgTiltedRieman->FixParameter(3, zmf);
989 fgTiltedRieman->FixParameter(4, dzmf);
991 fitter->ReleaseParameter(3);
992 fitter->ReleaseParameter(4);
993 offset = fitter->GetParameter(3);
994 slope = fitter->GetParameter(4);
997 // Calculate Curvarture
998 Double_t a = fitter->GetParameter(0);
999 Double_t b = fitter->GetParameter(1);
1000 Double_t c = fitter->GetParameter(2);
1001 Double_t curvature = 1.0 + b*b - c*a;
1002 if (curvature > 0.0)
1003 curvature = a / TMath::Sqrt(curvature);
1005 Double_t chi2track = fitter->GetChisquare()/Double_t(nPoints);
1007 // Update the tracklets
1009 for(Int_t iLayer = 0; iLayer < AliTRDtrackerV1::kNPlanes; iLayer++) {
1011 x = tracklets[iLayer].GetX0();
1017 // y: R^2 = (x - x0)^2 + (y - y0)^2
1018 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1019 // R = Sqrt() = 1/Curvature
1020 // => y = y0 +/- Sqrt(1/Curvature^2 - (x - x0)^2)
1021 Double_t res = (x * a + b); // = (x - x0)/y0
1023 res = 1.0 - c * a + b * b - res; // = (R^2 - (x - x0)^2)/y0^2
1025 res = TMath::Sqrt(res);
1026 y = (1.0 - res) / a;
1029 // dy: R^2 = (x - x0)^2 + (y - y0)^2
1030 // => y = +/- Sqrt(R^2 - (x - x0)^2) + y0
1031 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1032 // Curvature: cr = 1/R = a/Sqrt(1 + b^2 - c*a)
1033 // => dy/dx = (x - x0)/(1/(cr^2) - (x - x0)^2)
1034 Double_t x0 = -b / a;
1035 if (-c * a + b * b + 1 > 0) {
1036 if (1.0/(curvature * curvature) - (x - x0) * (x - x0) > 0.0) {
1037 Double_t yderiv = (x - x0) / TMath::Sqrt(1.0/(curvature * curvature) - (x - x0) * (x - x0));
1038 if (a < 0) yderiv *= -1.0;
1042 z = offset + slope * (x - xref);
1044 tracklets[iLayer].SetYref(0, y);
1045 tracklets[iLayer].SetYref(1, dy);
1046 tracklets[iLayer].SetZref(0, z);
1047 tracklets[iLayer].SetZref(1, dz);
1048 tracklets[iLayer].SetC(curvature);
1049 tracklets[iLayer].SetChi2(chi2track);
1052 /* if(fReconstructor->GetStreamLevel() >=5){
1053 TTreeSRedirector &cstreamer = *fgDebugStreamer;
1054 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1055 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1056 Double_t chi2z = CalculateChi2Z(tracklets, offset, slope, xref);
1057 cstreamer << "FitTiltedRieman0"
1058 << "EventNumber=" << eventNumber
1059 << "CandidateNumber=" << candidateNumber
1061 << "Chi2Z=" << chi2z
1068 //____________________________________________________________________
1069 Double_t AliTRDtrackerV1::FitLine(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t err, Int_t np, AliTrackPoint *points)
1071 AliTRDLeastSquare yfitter, zfitter;
1072 AliTRDcluster *cl = 0x0;
1074 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1076 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1077 if(!(tracklet = track->GetTracklet(ipl))) continue;
1078 if(!tracklet->IsOK()) continue;
1079 new(&work[ipl]) AliTRDseedV1(*tracklet);
1081 tracklets = &work[0];
1084 Double_t xref = CalculateReferenceX(tracklets);
1085 Double_t x, y, z, dx, ye, yr, tilt;
1086 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1087 if(!tracklets[ipl].IsOK()) continue;
1088 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1089 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1090 if (!tracklets[ipl].IsUsable(itb)) continue;
1094 zfitter.AddPoint(&dx, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1098 Double_t z0 = zfitter.GetFunctionParameter(0);
1099 Double_t dzdx = zfitter.GetFunctionParameter(1);
1100 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1101 if(!tracklets[ipl].IsOK()) continue;
1102 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1103 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1104 if (!tracklets[ipl].IsUsable(itb)) continue;
1108 tilt = tracklets[ipl].GetTilt();
1110 yr = y + tilt*(z - z0 - dzdx*dx);
1111 // error definition changes for the different calls
1112 ye = tilt*TMath::Sqrt(cl->GetSigmaZ2());
1113 ye += err ? tracklets[ipl].GetSigmaY() : 0.2;
1114 yfitter.AddPoint(&dx, yr, ye);
1118 Double_t y0 = yfitter.GetFunctionParameter(0);
1119 Double_t dydx = yfitter.GetFunctionParameter(1);
1120 Double_t chi2 = 0.;//yfitter.GetChisquare()/Double_t(nPoints);
1122 //update track points array
1125 for(int ip=0; ip<np; ip++){
1126 points[ip].GetXYZ(xyz);
1127 xyz[1] = y0 + dydx * (xyz[0] - xref);
1128 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1129 points[ip].SetXYZ(xyz);
1136 //_________________________________________________________________________
1137 Double_t AliTRDtrackerV1::FitRiemanTilt(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t sigError, Int_t np, AliTrackPoint *points)
1140 // Performs a Riemann fit taking tilting pad correction into account
1141 // The equation of a Riemann circle, where the y position is substituted by the
1142 // measured y-position taking pad tilting into account, has to be transformed
1143 // into a 4-dimensional hyperplane equation
1144 // Riemann circle: (x-x0)^2 + (y-y0)^2 -R^2 = 0
1145 // Measured y-Position: ymeas = y - tan(phiT)(zc - zt)
1146 // zc: center of the pad row
1147 // zt: z-position of the track
1148 // The z-position of the track is assumed to be linear dependent on the x-position
1149 // Transformed equation: a + b * u + c * t + d * v + e * w - 2 * (ymeas + tan(phiT) * zc) * t = 0
1150 // Transformation: u = 2 * x * t
1151 // v = 2 * tan(phiT) * t
1152 // w = 2 * tan(phiT) * (x - xref) * t
1153 // t = 1 / (x^2 + ymeas^2)
1154 // Parameters: a = -1/y0
1156 // c = (R^2 -x0^2 - y0^2)/y0
1159 // If the offset respectively the slope in z-position is impossible, the parameters are fixed using
1160 // results from the simple riemann fit. Afterwards the fit is redone.
1161 // The curvature is calculated according to the formula:
1162 // curv = a/(1 + b^2 + c*a) = 1/R
1164 // Paramters: - Array of tracklets (connected to the track candidate)
1165 // - Flag selecting the error definition
1166 // Output: - Chi2 values of the track (in Parameter list)
1168 TLinearFitter *fitter = GetTiltedRiemanFitter();
1169 fitter->StoreData(kTRUE);
1170 fitter->ClearPoints();
1171 AliTRDLeastSquare zfitter;
1172 AliTRDcluster *cl = 0x0;
1174 AliTRDseedV1 work[kNPlanes], *tracklet = 0x0;
1176 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1177 if(!(tracklet = track->GetTracklet(ipl))) continue;
1178 if(!tracklet->IsOK()) continue;
1179 new(&work[ipl]) AliTRDseedV1(*tracklet);
1181 tracklets = &work[0];
1184 Double_t xref = CalculateReferenceX(tracklets);
1185 Double_t x, y, z, t, tilt, dx, w, we;
1188 // Containers for Least-square fitter
1189 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
1190 if(!tracklets[ipl].IsOK()) continue;
1191 for(Int_t itb = 0; itb < fgNTimeBins; itb++){
1192 if(!(cl = tracklets[ipl].GetClusters(itb))) continue;
1193 if (!tracklets[ipl].IsUsable(itb)) continue;
1197 tilt = tracklets[ipl].GetTilt();
1201 uvt[0] = 2. * x * t;
1203 uvt[2] = 2. * tilt * t;
1204 uvt[3] = 2. * tilt * dx * t;
1205 w = 2. * (y + tilt*z) * t;
1206 // error definition changes for the different calls
1208 we *= sigError ? tracklets[ipl].GetSigmaY() : 0.2;
1209 fitter->AddPoint(uvt, w, we);
1210 zfitter.AddPoint(&x, z, static_cast<Double_t>(TMath::Sqrt(cl->GetSigmaZ2())));
1214 if(fitter->Eval()) return 1.E10;
1216 Double_t z0 = fitter->GetParameter(3);
1217 Double_t dzdx = fitter->GetParameter(4);
1220 // Linear fitter - not possible to make boundaries
1221 // Do not accept non possible z and dzdx combinations
1222 Bool_t accept = kTRUE;
1223 Double_t zref = 0.0;
1224 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
1225 if(!tracklets[iLayer].IsOK()) continue;
1226 zref = z0 + dzdx * (tracklets[iLayer].GetX0() - xref);
1227 if (TMath::Abs(tracklets[iLayer].GetZProb() - zref) > tracklets[iLayer].GetPadLength() * 0.5 + 1.0)
1232 Double_t dzmf = zfitter.GetFunctionParameter(1);
1233 Double_t zmf = zfitter.GetFunctionValue(&xref);
1234 fitter->FixParameter(3, zmf);
1235 fitter->FixParameter(4, dzmf);
1237 fitter->ReleaseParameter(3);
1238 fitter->ReleaseParameter(4);
1239 z0 = fitter->GetParameter(3); // = zmf ?
1240 dzdx = fitter->GetParameter(4); // = dzmf ?
1243 // Calculate Curvature
1244 Double_t a = fitter->GetParameter(0);
1245 Double_t b = fitter->GetParameter(1);
1246 Double_t c = fitter->GetParameter(2);
1247 Double_t y0 = 1. / a;
1248 Double_t x0 = -b * y0;
1249 Double_t R = TMath::Sqrt(y0*y0 + x0*x0 - c*y0);
1250 Double_t C = 1.0 + b*b - c*a;
1251 if (C > 0.0) C = a / TMath::Sqrt(C);
1253 // Calculate chi2 of the fit
1254 Double_t chi2 = fitter->GetChisquare()/Double_t(nPoints);
1256 // Update the tracklets
1258 for(Int_t ip = 0; ip < kNPlanes; ip++) {
1259 x = tracklets[ip].GetX0();
1260 Double_t tmp = TMath::Sqrt(R*R-(x-x0)*(x-x0));
1262 // y: R^2 = (x - x0)^2 + (y - y0)^2
1263 // => y = y0 +/- Sqrt(R^2 - (x - x0)^2)
1264 tracklets[ip].SetYref(0, y0 - (y0>0.?1.:-1)*tmp);
1265 // => dy/dx = (x - x0)/Sqrt(R^2 - (x - x0)^2)
1266 tracklets[ip].SetYref(1, (x - x0) / tmp);
1267 tracklets[ip].SetZref(0, z0 + dzdx * (x - xref));
1268 tracklets[ip].SetZref(1, dzdx);
1269 tracklets[ip].SetC(C);
1270 tracklets[ip].SetChi2(chi2);
1274 //update track points array
1277 for(int ip=0; ip<np; ip++){
1278 points[ip].GetXYZ(xyz);
1279 xyz[1] = y0 - (y0>0.?1.:-1)*TMath::Sqrt(R*R-(xyz[0]-x0)*(xyz[0]-x0));
1280 xyz[2] = z0 + dzdx * (xyz[0] - xref);
1281 points[ip].SetXYZ(xyz);
1285 /* if(fReconstructor->GetStreamLevel() >=5){
1286 TTreeSRedirector &cstreamer = *fgDebugStreamer;
1287 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
1288 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
1289 Double_t chi2z = CalculateChi2Z(tracklets, z0, dzdx, xref);
1290 cstreamer << "FitRiemanTilt"
1291 << "EventNumber=" << eventNumber
1292 << "CandidateNumber=" << candidateNumber
1294 << "Chi2Z=" << chi2z
1301 //____________________________________________________________________
1302 Double_t AliTRDtrackerV1::FitKalman(AliTRDtrackV1 *track, AliTRDseedV1 *tracklets, Bool_t up, Int_t np, AliTrackPoint *points)
1304 // Kalman filter implementation for the TRD.
1305 // It returns the positions of the fit in the array "points"
1307 // Author : A.Bercuci@gsi.de
1309 //printf("Start track @ x[%f]\n", track->GetX());
1311 //prepare marker points along the track
1312 Int_t ip = np ? 0 : 1;
1314 if((up?-1:1) * (track->GetX() - points[ip].GetX()) > 0.) break;
1315 //printf("AliTRDtrackerV1::FitKalman() : Skip track marker x[%d] = %7.3f. Before track start ( %7.3f ).\n", ip, points[ip].GetX(), track->GetX());
1318 //if(points) printf("First marker point @ x[%d] = %f\n", ip, points[ip].GetX());
1321 AliTRDseedV1 tracklet, *ptrTracklet = 0x0;
1323 //Loop through the TRD planes
1324 for (Int_t jplane = 0; jplane < kNPlanes; jplane++) {
1325 // GET TRACKLET OR BUILT IT
1326 Int_t iplane = up ? jplane : kNPlanes - 1 - jplane;
1328 if(!(ptrTracklet = &tracklets[iplane])) continue;
1330 if(!(ptrTracklet = track->GetTracklet(iplane))){
1331 /*AliTRDtrackerV1 *tracker = 0x0;
1332 if(!(tracker = dynamic_cast<AliTRDtrackerV1*>( AliTRDReconstructor::Tracker()))) continue;
1333 ptrTracklet = new(&tracklet) AliTRDseedV1(iplane);
1334 if(!tracker->MakeTracklet(ptrTracklet, track)) */
1338 if(!ptrTracklet->IsOK()) continue;
1340 Double_t x = ptrTracklet->GetX0();
1343 //don't do anything if next marker is after next update point.
1344 if((up?-1:1) * (points[ip].GetX() - x) - fgkMaxStep < 0) break;
1346 //printf("Propagate to x[%d] = %f\n", ip, points[ip].GetX());
1348 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1350 Double_t xyz[3]; // should also get the covariance
1351 track->GetXYZ(xyz); points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1354 //printf("plane[%d] tracklet[%p] x[%f]\n", iplane, ptrTracklet, x);
1356 //Propagate closer to the next update point
1357 if(((up?-1:1) * (x - track->GetX()) + fgkMaxStep < 0) && !PropagateToX(*track, x + (up?-1:1)*fgkMaxStep, fgkMaxStep)) return -1.;
1359 if(!AdjustSector(track)) return -1;
1360 if(TMath::Abs(track->GetSnp()) > fgkMaxSnp) return -1;
1362 //load tracklet to the tracker and the track
1364 if((index = FindTracklet(ptrTracklet)) < 0){
1365 ptrTracklet = SetTracklet(&tracklet);
1366 index = fTracklets->GetEntriesFast()-1;
1368 track->SetTracklet(ptrTracklet, index);*/
1371 // register tracklet to track with tracklet creation !!
1372 // PropagateBack : loaded tracklet to the tracker and update index
1373 // RefitInward : update index
1374 // MakeTrack : loaded tracklet to the tracker and update index
1375 if(!tracklets) track->SetTracklet(ptrTracklet, -1);
1378 //Calculate the mean material budget along the path inside the chamber
1379 Double_t xyz0[3]; track->GetXYZ(xyz0);
1380 Double_t alpha = track->GetAlpha();
1381 Double_t xyz1[3], y, z;
1382 if(!track->GetProlongation(x, y, z)) return -1;
1383 xyz1[0] = x * TMath::Cos(alpha) - y * TMath::Sin(alpha);
1384 xyz1[1] = +x * TMath::Sin(alpha) + y * TMath::Cos(alpha);
1387 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
1388 Double_t xrho = param[0]*param[4]; // density*length
1389 Double_t xx0 = param[1]; // radiation length
1391 //Propagate the track
1392 track->PropagateTo(x, xx0, xrho);
1393 if (!AdjustSector(track)) break;
1396 Double_t chi2 = track->GetPredictedChi2(ptrTracklet);
1397 if(chi2<1e+10) track->Update(ptrTracklet, chi2);
1401 //Reset material budget if 2 consecutive gold
1402 if(iplane>0 && track->GetTracklet(iplane-1) && ptrTracklet->GetN() + track->GetTracklet(iplane-1)->GetN() > 20) track->SetBudget(2, 0.);
1403 } // end planes loop
1407 if(((up?-1:1) * (points[ip].GetX() - track->GetX()) < 0) && !PropagateToX(*track, points[ip].GetX(), fgkMaxStep)) return -1.;
1409 Double_t xyz[3]; // should also get the covariance
1410 track->GetXYZ(xyz); points[ip].SetXYZ(xyz[0], xyz[1], xyz[2]);
1414 return track->GetChi2();
1417 //_________________________________________________________________________
1418 Float_t AliTRDtrackerV1::CalculateChi2Z(AliTRDseedV1 *tracklets, Double_t offset, Double_t slope, Double_t xref)
1421 // Calculates the chi2-value of the track in z-Direction including tilting pad correction.
1422 // A linear dependence on the x-value serves as a model.
1423 // The parameters are related to the tilted Riemann fit.
1424 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
1425 // - the offset for the reference x
1427 // - the reference x position
1428 // Output: - The Chi2 value of the track in z-Direction
1430 Float_t chi2Z = 0, nLayers = 0;
1431 for (Int_t iLayer = 0; iLayer < AliTRDgeometry::kNlayer; iLayer++) {
1432 if(!tracklets[iLayer].IsOK()) continue;
1433 Double_t z = offset + slope * (tracklets[iLayer].GetX0() - xref);
1434 chi2Z += TMath::Abs(tracklets[iLayer].GetMeanz() - z);
1437 chi2Z /= TMath::Max((nLayers - 3.0),1.0);
1441 //_____________________________________________________________________________
1442 Int_t AliTRDtrackerV1::PropagateToX(AliTRDtrackV1 &t, Double_t xToGo, Double_t maxStep)
1445 // Starting from current X-position of track <t> this function
1446 // extrapolates the track up to radial position <xToGo>.
1447 // Returns 1 if track reaches the plane, and 0 otherwise
1450 const Double_t kEpsilon = 0.00001;
1452 // Current track X-position
1453 Double_t xpos = t.GetX();
1455 // Direction: inward or outward
1456 Double_t dir = (xpos < xToGo) ? 1.0 : -1.0;
1458 while (((xToGo - xpos) * dir) > kEpsilon) {
1467 // The next step size
1468 Double_t step = dir * TMath::Min(TMath::Abs(xToGo-xpos),maxStep);
1470 // Get the global position of the starting point
1473 // X-position after next step
1476 // Get local Y and Z at the X-position of the next step
1477 if (!t.GetProlongation(x,y,z)) {
1478 return 0; // No prolongation possible
1481 // The global position of the end point of this prolongation step
1482 xyz1[0] = x * TMath::Cos(t.GetAlpha()) - y * TMath::Sin(t.GetAlpha());
1483 xyz1[1] = +x * TMath::Sin(t.GetAlpha()) + y * TMath::Cos(t.GetAlpha());
1486 // Calculate the mean material budget between start and
1487 // end point of this prolongation step
1488 AliTracker::MeanMaterialBudget(xyz0, xyz1, param);
1490 // Propagate the track to the X-position after the next step
1491 if (!t.PropagateTo(x,param[1],param[0]*param[4])) {
1495 // Rotate the track if necessary
1498 // New track X-position
1508 //_____________________________________________________________________________
1509 Int_t AliTRDtrackerV1::ReadClusters(TClonesArray* &array, TTree *clusterTree) const
1512 // Reads AliTRDclusters from the file.
1513 // The names of the cluster tree and branches
1514 // should match the ones used in AliTRDclusterizer::WriteClusters()
1517 Int_t nsize = Int_t(clusterTree->GetTotBytes() / (sizeof(AliTRDcluster)));
1518 TObjArray *clusterArray = new TObjArray(nsize+1000);
1520 TBranch *branch = clusterTree->GetBranch("TRDcluster");
1522 AliError("Can't get the branch !");
1525 branch->SetAddress(&clusterArray);
1528 array = new TClonesArray("AliTRDcluster", nsize);
1529 array->SetOwner(kTRUE);
1532 // Loop through all entries in the tree
1533 Int_t nEntries = (Int_t) clusterTree->GetEntries();
1536 AliTRDcluster *c = 0x0;
1537 for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {
1539 nbytes += clusterTree->GetEvent(iEntry);
1541 // Get the number of points in the detector
1542 Int_t nCluster = clusterArray->GetEntriesFast();
1543 for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) {
1544 if(!(c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster))) continue;
1546 new((*fClusters)[ncl++]) AliTRDcluster(*c);
1547 delete (clusterArray->RemoveAt(iCluster));
1551 delete clusterArray;
1556 //_____________________________________________________________________________
1557 Int_t AliTRDtrackerV1::LoadClusters(TTree *cTree)
1560 // Fills clusters into TRD tracking sectors
1563 if(!(fClusters = AliTRDReconstructor::GetClusters())){
1564 if (ReadClusters(fClusters, cTree)) {
1565 AliError("Problem with reading the clusters !");
1571 if(!fClusters->GetEntriesFast()){
1572 AliInfo("No TRD clusters");
1577 BuildTrackingContainers();
1579 //Int_t ncl = fClusters->GetEntriesFast();
1580 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1585 //_____________________________________________________________________________
1586 Int_t AliTRDtrackerV1::LoadClusters(TClonesArray *clusters)
1589 // Fills clusters into TRD tracking sectors
1590 // Function for use in the HLT
1592 if(!clusters || !clusters->GetEntriesFast()){
1593 AliInfo("No TRD clusters");
1597 fClusters = clusters;
1601 BuildTrackingContainers();
1603 //Int_t ncl = fClusters->GetEntriesFast();
1604 //AliInfo(Form("Clusters %d [%6.2f %% in the active volume]", ncl, 100.*float(nin)/ncl));
1610 //____________________________________________________________________
1611 Int_t AliTRDtrackerV1::BuildTrackingContainers()
1613 // Building tracking containers for clusters
1615 Int_t nin =0, icl = fClusters->GetEntriesFast();
1617 AliTRDcluster *c = (AliTRDcluster *) fClusters->UncheckedAt(icl);
1618 if(c->IsInChamber()) nin++;
1619 Int_t detector = c->GetDetector();
1620 Int_t sector = fGeom->GetSector(detector);
1621 Int_t stack = fGeom->GetStack(detector);
1622 Int_t layer = fGeom->GetLayer(detector);
1624 fTrSec[sector].GetChamber(stack, layer, kTRUE)->InsertCluster(c, icl);
1627 for(int isector =0; isector<AliTRDgeometry::kNsector; isector++){
1628 if(!fTrSec[isector].GetNChambers()) continue;
1629 fTrSec[isector].Init(fReconstructor);
1637 //____________________________________________________________________
1638 void AliTRDtrackerV1::UnloadClusters()
1641 // Clears the arrays of clusters and tracks. Resets sectors and timebins
1644 if(fTracks) fTracks->Delete();
1645 if(fTracklets) fTracklets->Delete();
1646 if(fClusters && IsClustersOwner()) fClusters->Delete();
1648 for (int i = 0; i < AliTRDgeometry::kNsector; i++) fTrSec[i].Clear();
1650 // Increment the Event Number
1651 AliTRDtrackerDebug::SetEventNumber(AliTRDtrackerDebug::GetEventNumber() + 1);
1654 //_____________________________________________________________________________
1655 Bool_t AliTRDtrackerV1::AdjustSector(AliTRDtrackV1 *track)
1658 // Rotates the track when necessary
1661 Double_t alpha = AliTRDgeometry::GetAlpha();
1662 Double_t y = track->GetY();
1663 Double_t ymax = track->GetX()*TMath::Tan(0.5*alpha);
1666 if (!track->Rotate( alpha)) {
1670 else if (y < -ymax) {
1671 if (!track->Rotate(-alpha)) {
1681 //____________________________________________________________________
1682 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(AliTRDtrackV1 *track, Int_t p, Int_t &idx)
1684 // Find tracklet for TRD track <track>
1693 // Detailed description
1695 idx = track->GetTrackletIndex(p);
1696 AliTRDseedV1 *tracklet = (idx==0xffff) ? 0x0 : (AliTRDseedV1*)fTracklets->UncheckedAt(idx);
1701 //____________________________________________________________________
1702 AliTRDseedV1* AliTRDtrackerV1::SetTracklet(AliTRDseedV1 *tracklet)
1704 // Add this tracklet to the list of tracklets stored in the tracker
1707 // - tracklet : pointer to the tracklet to be added to the list
1710 // - the index of the new tracklet in the tracker tracklets list
1712 // Detailed description
1713 // Build the tracklets list if it is not yet created (late initialization)
1714 // and adds the new tracklet to the list.
1717 fTracklets = new TClonesArray("AliTRDseedV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1718 fTracklets->SetOwner(kTRUE);
1720 Int_t nentries = fTracklets->GetEntriesFast();
1721 return new ((*fTracklets)[nentries]) AliTRDseedV1(*tracklet);
1724 //____________________________________________________________________
1725 AliTRDtrackV1* AliTRDtrackerV1::SetTrack(AliTRDtrackV1 *track)
1727 // Add this track to the list of tracks stored in the tracker
1730 // - track : pointer to the track to be added to the list
1733 // - the pointer added
1735 // Detailed description
1736 // Build the tracks list if it is not yet created (late initialization)
1737 // and adds the new track to the list.
1740 fTracks = new TClonesArray("AliTRDtrackV1", AliTRDgeometry::Nsector()*kMaxTracksStack);
1741 fTracks->SetOwner(kTRUE);
1743 Int_t nentries = fTracks->GetEntriesFast();
1744 return new ((*fTracks)[nentries]) AliTRDtrackV1(*track);
1749 //____________________________________________________________________
1750 Int_t AliTRDtrackerV1::Clusters2TracksSM(Int_t sector, AliESDEvent *esd)
1753 // Steer tracking for one SM.
1756 // sector : Array of (SM) propagation layers containing clusters
1757 // esd : The current ESD event. On output it contains the also
1758 // the ESD (TRD) tracks found in this SM.
1761 // Number of tracks found in this TRD supermodule.
1763 // Detailed description
1765 // 1. Unpack AliTRDpropagationLayers objects for each stack.
1766 // 2. Launch stack tracking.
1767 // See AliTRDtrackerV1::Clusters2TracksStack() for details.
1768 // 3. Pack results in the ESD event.
1771 // allocate space for esd tracks in this SM
1772 TClonesArray esdTrackList("AliESDtrack", 2*kMaxTracksStack);
1773 esdTrackList.SetOwner();
1776 Int_t nChambers = 0;
1777 AliTRDtrackingChamber **stack = 0x0, *chamber = 0x0;
1778 for(int istack = 0; istack<AliTRDgeometry::kNstack; istack++){
1779 if(!(stack = fTrSec[sector].GetStack(istack))) continue;
1781 for(int ilayer=0; ilayer<AliTRDgeometry::kNlayer; ilayer++){
1782 if(!(chamber = stack[ilayer])) continue;
1783 if(chamber->GetNClusters() < fgNTimeBins * fReconstructor->GetRecoParam() ->GetFindableClusters()) continue;
1785 //AliInfo(Form("sector %d stack %d layer %d clusters %d", sector, istack, ilayer, chamber->GetNClusters()));
1787 if(nChambers < 4) continue;
1788 //AliInfo(Form("Doing stack %d", istack));
1789 nTracks += Clusters2TracksStack(stack, &esdTrackList);
1791 //AliInfo(Form("Found %d tracks in SM %d [%d]\n", nTracks, sector, esd->GetNumberOfTracks()));
1793 for(int itrack=0; itrack<nTracks; itrack++)
1794 esd->AddTrack((AliESDtrack*)esdTrackList[itrack]);
1796 // Reset Track and Candidate Number
1797 AliTRDtrackerDebug::SetCandidateNumber(0);
1798 AliTRDtrackerDebug::SetTrackNumber(0);
1802 //____________________________________________________________________
1803 Int_t AliTRDtrackerV1::Clusters2TracksStack(AliTRDtrackingChamber **stack, TClonesArray *esdTrackList)
1806 // Make tracks in one TRD stack.
1809 // layer : Array of stack propagation layers containing clusters
1810 // esdTrackList : Array of ESD tracks found by the stand alone tracker.
1811 // On exit the tracks found in this stack are appended.
1814 // Number of tracks found in this stack.
1816 // Detailed description
1818 // 1. Find the 3 most useful seeding chambers. See BuildSeedingConfigs() for details.
1819 // 2. Steer AliTRDtrackerV1::MakeSeeds() for 3 seeding layer configurations.
1820 // See AliTRDtrackerV1::MakeSeeds() for more details.
1821 // 3. Arrange track candidates in decreasing order of their quality
1822 // 4. Classify tracks in 5 categories according to:
1823 // a) number of layers crossed
1825 // 5. Sign clusters by tracks in decreasing order of track quality
1826 // 6. Build AliTRDtrack out of seeding tracklets
1828 // 8. Build ESD track and register it to the output list
1831 AliTRDtrackingChamber *chamber = 0x0;
1832 AliTRDseedV1 sseed[kMaxTracksStack*6]; // to be initialized
1833 Int_t pars[4]; // MakeSeeds parameters
1835 //Double_t alpha = AliTRDgeometry::GetAlpha();
1836 //Double_t shift = .5 * alpha;
1837 Int_t configs[kNConfigs];
1839 // Build initial seeding configurations
1840 Double_t quality = BuildSeedingConfigs(stack, configs);
1841 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
1842 AliInfo(Form("Plane config %d %d %d Quality %f"
1843 , configs[0], configs[1], configs[2], quality));
1846 // Initialize contors
1847 Int_t ntracks, // number of TRD track candidates
1848 ntracks1, // number of registered TRD tracks/iter
1849 ntracks2 = 0; // number of all registered TRD tracks in stack
1852 // Loop over seeding configurations
1853 ntracks = 0; ntracks1 = 0;
1854 for (Int_t iconf = 0; iconf<3; iconf++) {
1855 pars[0] = configs[iconf];
1857 ntracks = MakeSeeds(stack, &sseed[6*ntracks], pars);
1858 if(ntracks == kMaxTracksStack) break;
1860 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1) AliInfo(Form("Candidate TRD tracks %d in iteration %d.", ntracks, fSieveSeeding));
1864 // Sort the seeds according to their quality
1865 Int_t sort[kMaxTracksStack];
1866 TMath::Sort(ntracks, fTrackQuality, sort, kTRUE);
1868 // Initialize number of tracks so far and logic switches
1869 Int_t ntracks0 = esdTrackList->GetEntriesFast();
1870 Bool_t signedTrack[kMaxTracksStack];
1871 Bool_t fakeTrack[kMaxTracksStack];
1872 for (Int_t i=0; i<ntracks; i++){
1873 signedTrack[i] = kFALSE;
1874 fakeTrack[i] = kFALSE;
1876 //AliInfo("Selecting track candidates ...");
1878 // Sieve clusters in decreasing order of track quality
1879 Double_t trackParams[7];
1880 // AliTRDseedV1 *lseed = 0x0;
1881 Int_t jSieve = 0, candidates;
1883 //AliInfo(Form("\t\tITER = %i ", jSieve));
1885 // Check track candidates
1887 for (Int_t itrack = 0; itrack < ntracks; itrack++) {
1888 Int_t trackIndex = sort[itrack];
1889 if (signedTrack[trackIndex] || fakeTrack[trackIndex]) continue;
1892 // Calculate track parameters from tracklets seeds
1893 Int_t labelsall[1000];
1894 Int_t nlabelsall = 0;
1895 Int_t naccepted = 0;
1900 for (Int_t jLayer = 0; jLayer < kNPlanes; jLayer++) {
1901 Int_t jseed = kNPlanes*trackIndex+jLayer;
1902 if(!sseed[jseed].IsOK()) continue;
1903 if (TMath::Abs(sseed[jseed].GetYref(0) / sseed[jseed].GetX0()) < 0.15) findable++;
1905 sseed[jseed].UpdateUsed();
1906 ncl += sseed[jseed].GetN2();
1907 nused += sseed[jseed].GetNUsed();
1911 for (Int_t itime = 0; itime < fgNTimeBins; itime++) {
1912 if(!sseed[jseed].IsUsable(itime)) continue;
1914 Int_t tindex = 0, ilab = 0;
1915 while(ilab<3 && (tindex = sseed[jseed].GetClusters(itime)->GetLabel(ilab)) >= 0){
1916 labelsall[nlabelsall++] = tindex;
1921 // Filter duplicated tracks
1923 //printf("Skip %d nused %d\n", trackIndex, nused);
1924 fakeTrack[trackIndex] = kTRUE;
1927 if (Float_t(nused)/ncl >= .25){
1928 //printf("Skip %d nused/ncl >= .25\n", trackIndex);
1929 fakeTrack[trackIndex] = kTRUE;
1934 Bool_t skip = kFALSE;
1937 if(nlayers < 6) {skip = kTRUE; break;}
1938 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
1942 if(nlayers < findable){skip = kTRUE; break;}
1943 if(TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -4.){skip = kTRUE; break;}
1947 if ((nlayers == findable) || (nlayers == 6)) { skip = kTRUE; break;}
1948 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -6.0){skip = kTRUE; break;}
1952 if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) < -5.){skip = kTRUE; break;}
1956 if (nlayers == 3){skip = kTRUE; break;}
1957 //if (TMath::Log(1.E-9+fTrackQuality[trackIndex]) - nused/(nlayers-3.0) < -15.0){skip = kTRUE; break;}
1962 //printf("REJECTED : %d [%d] nlayers %d trackQuality = %e nused %d\n", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused);
1965 signedTrack[trackIndex] = kTRUE;
1968 // Build track label - what happens if measured data ???
1972 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
1973 Int_t jseed = kNPlanes*trackIndex+iLayer;
1974 if(!sseed[jseed].IsOK()) continue;
1975 for(int ilab=0; ilab<2; ilab++){
1976 if(sseed[jseed].GetLabels(ilab) < 0) continue;
1977 labels[nlab] = sseed[jseed].GetLabels(ilab);
1981 Freq(nlab,labels,outlab,kFALSE);
1982 Int_t label = outlab[0];
1983 Int_t frequency = outlab[1];
1984 Freq(nlabelsall,labelsall,outlab,kFALSE);
1985 Int_t label1 = outlab[0];
1986 Int_t label2 = outlab[2];
1987 Float_t fakeratio = (naccepted - outlab[1]) / Float_t(naccepted);
1991 AliTRDcluster *cl = 0x0; Int_t clusterIndex = -1;
1992 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
1993 Int_t jseed = kNPlanes*trackIndex+jLayer;
1994 if(!sseed[jseed].IsOK()) continue;
1995 if(TMath::Abs(sseed[jseed].GetYfit(1) - sseed[jseed].GetYfit(1)) >= .2) continue; // check this condition with Marian
1996 sseed[jseed].UseClusters();
1999 while(!(cl = sseed[jseed].GetClusters(ic))) ic++;
2000 clusterIndex = sseed[jseed].GetIndexes(ic);
2006 // Build track parameters
2007 AliTRDseedV1 *lseed =&sseed[trackIndex*6];
2009 while(idx<3 && !lseed->IsOK()) {
2013 Double_t cR = lseed->GetC();
2014 Double_t x = lseed->GetX0();// - 3.5;
2015 trackParams[0] = x; //NEW AB
2016 trackParams[1] = lseed->GetYref(0); // lseed->GetYat(x);
2017 trackParams[2] = lseed->GetZref(0); // lseed->GetZat(x);
2018 trackParams[3] = TMath::Sin(TMath::ATan(lseed->GetYref(1)));
2019 trackParams[4] = lseed->GetZref(1) / TMath::Sqrt(1. + lseed->GetYref(1) * lseed->GetYref(1));
2020 trackParams[5] = cR;
2021 Int_t ich = 0; while(!(chamber = stack[ich])) ich++;
2022 trackParams[6] = fGeom->GetSector(chamber->GetDetector());/* *alpha+shift; // Supermodule*/
2024 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2025 AliInfo(Form("Track %d [%d] nlayers %d trackQuality = %e nused %d, yref = %3.3f", itrack, trackIndex, nlayers, fTrackQuality[trackIndex], nused, trackParams[1]));
2027 Int_t nclusters = 0;
2028 AliTRDseedV1 *dseed[6];
2029 for(int is=0; is<6; is++){
2030 dseed[is] = new AliTRDseedV1(sseed[trackIndex*6+is]);
2031 dseed[is]->SetOwner();
2032 nclusters += sseed[is].GetN2();
2034 //Int_t eventNrInFile = esd->GetEventNumberInFile();
2035 //AliInfo(Form("Number of clusters %d.", nclusters));
2036 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2037 Int_t trackNumber = AliTRDtrackerDebug::GetTrackNumber();
2038 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2039 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2040 cstreamer << "Clusters2TracksStack"
2041 << "EventNumber=" << eventNumber
2042 << "TrackNumber=" << trackNumber
2043 << "CandidateNumber=" << candidateNumber
2044 << "Iter=" << fSieveSeeding
2045 << "Like=" << fTrackQuality[trackIndex]
2046 << "S0.=" << dseed[0]
2047 << "S1.=" << dseed[1]
2048 << "S2.=" << dseed[2]
2049 << "S3.=" << dseed[3]
2050 << "S4.=" << dseed[4]
2051 << "S5.=" << dseed[5]
2052 << "p0=" << trackParams[0]
2053 << "p1=" << trackParams[1]
2054 << "p2=" << trackParams[2]
2055 << "p3=" << trackParams[3]
2056 << "p4=" << trackParams[4]
2057 << "p5=" << trackParams[5]
2058 << "p6=" << trackParams[6]
2059 << "Label=" << label
2060 << "Label1=" << label1
2061 << "Label2=" << label2
2062 << "FakeRatio=" << fakeratio
2063 << "Freq=" << frequency
2065 << "NLayers=" << nlayers
2066 << "Findable=" << findable
2067 << "NUsed=" << nused
2071 AliTRDtrackV1 *track = MakeTrack(&sseed[trackIndex*kNPlanes], trackParams);
2073 AliWarning("Fail to build a TRD Track.");
2076 //AliInfo("End of MakeTrack()");
2077 AliESDtrack *esdTrack = new ((*esdTrackList)[ntracks0++]) AliESDtrack();
2078 esdTrack->UpdateTrackParams(track, AliESDtrack::kTRDout);
2079 esdTrack->SetLabel(track->GetLabel());
2080 track->UpdateESDtrack(esdTrack);
2081 // write ESD-friends if neccessary
2082 if (fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 0){
2083 AliTRDtrackV1 *calibTrack = new AliTRDtrackV1(*track);
2084 calibTrack->SetOwner();
2085 esdTrack->AddCalibObject(calibTrack);
2088 AliTRDtrackerDebug::SetTrackNumber(AliTRDtrackerDebug::GetTrackNumber() + 1);
2092 } while(jSieve<5 && candidates); // end track candidates sieve
2093 if(!ntracks1) break;
2095 // increment counters
2096 ntracks2 += ntracks1;
2099 // Rebuild plane configurations and indices taking only unused clusters into account
2100 quality = BuildSeedingConfigs(stack, configs);
2101 if(quality < 1.E-7) break; //fReconstructor->GetRecoParam() ->GetPlaneQualityThreshold()) break;
2103 for(Int_t ip = 0; ip < kNPlanes; ip++){
2104 if(!(chamber = stack[ip])) continue;
2105 chamber->Build(fGeom);//Indices(fSieveSeeding);
2108 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2109 AliInfo(Form("Sieve level %d Plane config %d %d %d Quality %f", fSieveSeeding, configs[0], configs[1], configs[2], quality));
2111 } while(fSieveSeeding<10); // end stack clusters sieve
2115 //AliInfo(Form("Registered TRD tracks %d in stack %d.", ntracks2, pars[1]));
2120 //___________________________________________________________________
2121 Double_t AliTRDtrackerV1::BuildSeedingConfigs(AliTRDtrackingChamber **stack, Int_t *configs)
2124 // Assign probabilities to chambers according to their
2125 // capability of producing seeds.
2129 // layers : Array of stack propagation layers for all 6 chambers in one stack
2130 // configs : On exit array of configuration indexes (see GetSeedingConfig()
2131 // for details) in the decreasing order of their seeding probabilities.
2135 // Return top configuration quality
2137 // Detailed description:
2139 // To each chamber seeding configuration (see GetSeedingConfig() for
2140 // the list of all configurations) one defines 2 quality factors:
2141 // - an apriori topological quality (see GetSeedingConfig() for details) and
2142 // - a data quality based on the uniformity of the distribution of
2143 // clusters over the x range (time bins population). See CookChamberQA() for details.
2144 // The overall chamber quality is given by the product of this 2 contributions.
2147 Double_t chamberQ[kNPlanes];
2148 AliTRDtrackingChamber *chamber = 0x0;
2149 for(int iplane=0; iplane<kNPlanes; iplane++){
2150 if(!(chamber = stack[iplane])) continue;
2151 chamberQ[iplane] = (chamber = stack[iplane]) ? chamber->GetQuality() : 0.;
2154 Double_t tconfig[kNConfigs];
2156 for(int iconf=0; iconf<kNConfigs; iconf++){
2157 GetSeedingConfig(iconf, planes);
2158 tconfig[iconf] = fgTopologicQA[iconf];
2159 for(int iplane=0; iplane<4; iplane++) tconfig[iconf] *= chamberQ[planes[iplane]];
2162 TMath::Sort((Int_t)kNConfigs, tconfig, configs, kTRUE);
2163 // AliInfo(Form("q[%d] = %f", configs[0], tconfig[configs[0]]));
2164 // AliInfo(Form("q[%d] = %f", configs[1], tconfig[configs[1]]));
2165 // AliInfo(Form("q[%d] = %f", configs[2], tconfig[configs[2]]));
2167 return tconfig[configs[0]];
2170 //____________________________________________________________________
2171 Int_t AliTRDtrackerV1::MakeSeeds(AliTRDtrackingChamber **stack, AliTRDseedV1 *sseed, Int_t *ipar)
2174 // Make tracklet seeds in the TRD stack.
2177 // layers : Array of stack propagation layers containing clusters
2178 // sseed : Array of empty tracklet seeds. On exit they are filled.
2179 // ipar : Control parameters:
2180 // ipar[0] -> seeding chambers configuration
2181 // ipar[1] -> stack index
2182 // ipar[2] -> number of track candidates found so far
2185 // Number of tracks candidates found.
2187 // Detailed description
2189 // The following steps are performed:
2190 // 1. Select seeding layers from seeding chambers
2191 // 2. Select seeding clusters from the seeding AliTRDpropagationLayerStack.
2192 // The clusters are taken from layer 3, layer 0, layer 1 and layer 2, in
2193 // this order. The parameters controling the range of accepted clusters in
2194 // layer 0, 1, and 2 are defined in AliTRDchamberTimeBin::BuildCond().
2195 // 3. Helix fit of the cluster set. (see AliTRDtrackerFitter::FitRieman(AliTRDcluster**))
2196 // 4. Initialize seeding tracklets in the seeding chambers.
2198 // Chi2 in the Y direction less than threshold ... (1./(3. - sLayer))
2199 // Chi2 in the Z direction less than threshold ... (1./(3. - sLayer))
2200 // 6. Attach clusters to seeding tracklets and find linear approximation of
2201 // the tracklet (see AliTRDseedV1::AttachClustersIter()). The number of used
2202 // clusters used by current seeds should not exceed ... (25).
2204 // All 4 seeding tracklets should be correctly constructed (see
2205 // AliTRDseedV1::AttachClustersIter())
2206 // 8. Helix fit of the seeding tracklets
2208 // Likelihood calculation of the fit. (See AliTRDtrackerV1::CookLikelihood() for details)
2209 // 10. Extrapolation of the helix fit to the other 2 chambers:
2210 // a) Initialization of extrapolation tracklet with fit parameters
2211 // b) Helix fit of tracklets
2212 // c) Attach clusters and linear interpolation to extrapolated tracklets
2213 // d) Helix fit of tracklets
2214 // 11. Improve seeding tracklets quality by reassigning clusters.
2215 // See AliTRDtrackerV1::ImproveSeedQuality() for details.
2216 // 12. Helix fit of all 6 seeding tracklets and chi2 calculation
2217 // 13. Hyperplane fit and track quality calculation. See AliTRDtrackerFitter::FitHyperplane() for details.
2218 // 14. Cooking labels for tracklets. Should be done only for MC
2219 // 15. Register seeds.
2222 AliTRDtrackingChamber *chamber = 0x0;
2223 AliTRDcluster *c[4] = {0x0, 0x0, 0x0, 0x0}; // initilize seeding clusters
2224 AliTRDseedV1 *cseed = &sseed[0]; // initialize tracklets for first track
2225 Int_t ncl, mcl; // working variable for looping over clusters
2226 Int_t index[AliTRDchamberTimeBin::kMaxClustersLayer], jndex[AliTRDchamberTimeBin::kMaxClustersLayer];
2228 // chi2[0] = tracklet chi2 on the Z direction
2229 // chi2[1] = tracklet chi2 on the R direction
2232 // Default positions for the anode wire in all 6 Layers in case of a stack with missing clusters
2233 // Positions taken using cosmic data taken with SM3 after rebuild
2234 Double_t x_def[kNPlanes] = {300.2, 312.8, 325.4, 338, 350.6, 363.2};
2236 // this should be data member of AliTRDtrack
2237 Double_t seedQuality[kMaxTracksStack];
2239 // unpack control parameters
2240 Int_t config = ipar[0];
2241 Int_t ntracks = ipar[1];
2242 Int_t planes[kNSeedPlanes]; GetSeedingConfig(config, planes);
2244 // Init chambers geometry
2245 Int_t ic = 0; while(!(chamber = stack[ic])) ic++;
2246 Int_t istack = fGeom->GetStack(chamber->GetDetector());
2247 Double_t hL[kNPlanes]; // Tilting angle
2248 Float_t padlength[kNPlanes]; // pad lenghts
2249 AliTRDpadPlane *pp = 0x0;
2250 for(int iplane=0; iplane<kNPlanes; iplane++){
2251 pp = fGeom->GetPadPlane(iplane, istack);
2252 hL[iplane] = TMath::Tan(TMath::DegToRad()*pp->GetTiltingAngle());
2253 padlength[iplane] = pp->GetLengthIPad();
2256 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
2257 AliInfo(Form("Making seeds Stack[%d] Config[%d] Tracks[%d]...", istack, config, ntracks));
2262 for(int isl=0; isl<kNSeedPlanes; isl++){
2263 if(!(chamber = stack[planes[isl]])) continue;
2264 if(!chamber->GetSeedingLayer(fSeedTB[isl], fGeom, fReconstructor)) continue;
2267 if(nlayers < 4) return 0;
2270 // Start finding seeds
2271 Double_t cond0[4], cond1[4], cond2[4];
2273 while((c[3] = (*fSeedTB[3])[icl++])){
2275 fSeedTB[0]->BuildCond(c[3], cond0, 0);
2276 fSeedTB[0]->GetClusters(cond0, index, ncl);
2277 //printf("Found c[3] candidates 0 %d\n", ncl);
2280 c[0] = (*fSeedTB[0])[index[jcl++]];
2282 Double_t dx = c[3]->GetX() - c[0]->GetX();
2283 Double_t theta = (c[3]->GetZ() - c[0]->GetZ())/dx;
2284 Double_t phi = (c[3]->GetY() - c[0]->GetY())/dx;
2285 fSeedTB[1]->BuildCond(c[0], cond1, 1, theta, phi);
2286 fSeedTB[1]->GetClusters(cond1, jndex, mcl);
2287 //printf("Found c[0] candidates 1 %d\n", mcl);
2291 c[1] = (*fSeedTB[1])[jndex[kcl++]];
2293 fSeedTB[2]->BuildCond(c[1], cond2, 2, theta, phi);
2294 c[2] = fSeedTB[2]->GetNearestCluster(cond2);
2295 //printf("Found c[1] candidate 2 %p\n", c[2]);
2298 // AliInfo("Seeding clusters found. Building seeds ...");
2299 // 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());
2301 for (Int_t il = 0; il < 6; il++) cseed[il].Reset();
2305 AliTRDseedV1 *tseed = 0x0;
2306 for(int iLayer=0; iLayer<kNPlanes; iLayer++){
2307 tseed = &cseed[iLayer];
2308 tseed->SetPlane(iLayer);
2309 tseed->SetTilt(hL[iLayer]);
2310 tseed->SetPadLength(padlength[iLayer]);
2311 tseed->SetReconstructor(fReconstructor);
2312 Double_t x_anode = stack[iLayer] ? stack[iLayer]->GetX() : x_def[iLayer];
2313 tseed->SetX0(x_anode);
2314 //if(stack[iLayer]) tseed->SetX0(stack[iLayer]->GetX());
2315 tseed->Init(GetRiemanFitter());
2318 Bool_t isFake = kFALSE;
2319 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2320 if (c[0]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2321 if (c[1]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2322 if (c[2]->GetLabel(0) != c[3]->GetLabel(0)) isFake = kTRUE;
2325 for(Int_t l = 0; l < kNSeedPlanes; l++) xpos[l] = fSeedTB[l]->GetX();
2327 for(int il=0; il<4; il++) yref[il] = cseed[planes[il]].GetYref(0);
2328 Int_t ll = c[3]->GetLabel(0);
2329 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2330 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2331 AliRieman *rim = GetRiemanFitter();
2332 TTreeSRedirector &cs0 = *fgDebugStreamer;
2334 <<"EventNumber=" << eventNumber
2335 <<"CandidateNumber=" << candidateNumber
2336 <<"isFake=" << isFake
2337 <<"config=" << config
2339 <<"chi2z=" << chi2[0]
2340 <<"chi2y=" << chi2[1]
2341 <<"Y2exp=" << cond2[0]
2342 <<"Z2exp=" << cond2[1]
2343 <<"X0=" << xpos[0] //layer[sLayer]->GetX()
2344 <<"X1=" << xpos[1] //layer[sLayer + 1]->GetX()
2345 <<"X2=" << xpos[2] //layer[sLayer + 2]->GetX()
2346 <<"X3=" << xpos[3] //layer[sLayer + 3]->GetX()
2347 <<"yref0=" << yref[0]
2348 <<"yref1=" << yref[1]
2349 <<"yref2=" << yref[2]
2350 <<"yref3=" << yref[3]
2355 <<"Seed0.=" << &cseed[planes[0]]
2356 <<"Seed1.=" << &cseed[planes[1]]
2357 <<"Seed2.=" << &cseed[planes[2]]
2358 <<"Seed3.=" << &cseed[planes[3]]
2359 <<"RiemanFitter.=" << rim
2363 if(chi2[0] > fReconstructor->GetRecoParam() ->GetChi2Z()/*7./(3. - sLayer)*//*iter*/){
2364 //AliInfo(Form("Failed chi2 filter on chi2Z [%f].", chi2[0]));
2365 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2368 if(chi2[1] > fReconstructor->GetRecoParam() ->GetChi2Y()/*1./(3. - sLayer)*//*iter*/){
2369 //AliInfo(Form("Failed chi2 filter on chi2Y [%f].", chi2[1]));
2370 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2373 //AliInfo("Passed chi2 filter.");
2375 // try attaching clusters to tracklets
2378 for(int iLayer=0; iLayer<kNSeedPlanes; iLayer++){
2379 Int_t jLayer = planes[iLayer];
2380 if(!cseed[jLayer].AttachClustersIter(stack[jLayer], 5., kFALSE, c[iLayer])) continue;
2381 nUsedCl += cseed[jLayer].GetNUsed();
2382 if(nUsedCl > 25) break;
2385 if(mlayers < kNSeedPlanes){
2386 //AliInfo(Form("Failed updating all seeds %d [%d].", mlayers, kNSeedPlanes));
2387 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2390 // fit tracklets and cook likelihood
2391 FitTiltedRieman(&cseed[0], kTRUE);// Update Seeds and calculate Likelihood
2392 chi2[0] = GetChi2Y(&cseed[0]);
2393 chi2[1] = GetChi2Z(&cseed[0]);
2394 //Chi2 definitions in testing stage
2395 //chi2[0] = GetChi2YTest(&cseed[0]);
2396 //chi2[1] = GetChi2ZTest(&cseed[0]);
2397 Double_t like = CookLikelihood(&cseed[0], planes, chi2); // to be checked
2399 if (TMath::Log(1.E-9 + like) < fReconstructor->GetRecoParam() ->GetTrackLikelihood()){
2400 //AliInfo(Form("Failed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2401 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2404 //AliInfo(Form("Passed likelihood %f[%e].", TMath::Log(1.E-9 + like), like));
2406 // book preliminary results
2407 seedQuality[ntracks] = like;
2408 fSeedLayer[ntracks] = config;/*sLayer;*/
2410 // attach clusters to the extrapolation seeds
2412 GetExtrapolationConfig(config, lextrap);
2413 Int_t nusedf = 0; // debug value
2414 for(int iLayer=0; iLayer<2; iLayer++){
2415 Int_t jLayer = lextrap[iLayer];
2416 if(!(chamber = stack[jLayer])) continue;
2418 // fit extrapolated seed
2419 if ((jLayer == 0) && !(cseed[1].IsOK())) continue;
2420 if ((jLayer == 5) && !(cseed[4].IsOK())) continue;
2421 AliTRDseedV1 pseed = cseed[jLayer];
2422 if(!pseed.AttachClustersIter(chamber, 1000.)) continue;
2423 cseed[jLayer] = pseed;
2424 nusedf += cseed[jLayer].GetNUsed(); // debug value
2425 FitTiltedRieman(cseed, kTRUE);
2428 // AliInfo("Extrapolation done.");
2429 // Debug Stream containing all the 6 tracklets
2430 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2431 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2432 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2433 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2434 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2435 cstreamer << "MakeSeeds1"
2436 << "EventNumber=" << eventNumber
2437 << "CandidateNumber=" << candidateNumber
2438 << "S0.=" << &cseed[0]
2439 << "S1.=" << &cseed[1]
2440 << "S2.=" << &cseed[2]
2441 << "S3.=" << &cseed[3]
2442 << "S4.=" << &cseed[4]
2443 << "S5.=" << &cseed[5]
2444 << "FitterT.=" << tiltedRieman
2448 if(fReconstructor->GetRecoParam()->HasImproveTracklets() && ImproveSeedQuality(stack, cseed) < 4){
2449 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2452 //AliInfo("Improve seed quality done.");
2454 // fit full track and cook likelihoods
2455 // Double_t curv = FitRieman(&cseed[0], chi2);
2456 // Double_t chi2ZF = chi2[0] / TMath::Max((mlayers - 3.), 1.);
2457 // Double_t chi2RF = chi2[1] / TMath::Max((mlayers - 3.), 1.);
2459 // do the final track fitting (Once with vertex constraint and once without vertex constraint)
2460 Double_t chi2Vals[3];
2461 chi2Vals[0] = FitTiltedRieman(&cseed[0], kFALSE);
2462 if(fReconstructor->GetRecoParam() ->IsVertexConstrained())
2463 chi2Vals[1] = FitTiltedRiemanConstraint(&cseed[0], GetZ()); // Do Vertex Constrained fit if desired
2466 chi2Vals[2] = GetChi2Z(&cseed[0]) / TMath::Max((mlayers - 3.), 1.);
2467 // Chi2 definitions in testing stage
2468 //chi2Vals[2] = GetChi2ZTest(&cseed[0]);
2469 fTrackQuality[ntracks] = CalculateTrackLikelihood(&cseed[0], &chi2Vals[0]);
2470 //AliInfo("Hyperplane fit done\n");
2472 // finalize tracklets
2476 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2477 if (!cseed[iLayer].IsOK()) continue;
2479 if (cseed[iLayer].GetLabels(0) >= 0) {
2480 labels[nlab] = cseed[iLayer].GetLabels(0);
2484 if (cseed[iLayer].GetLabels(1) >= 0) {
2485 labels[nlab] = cseed[iLayer].GetLabels(1);
2489 Freq(nlab,labels,outlab,kFALSE);
2490 Int_t label = outlab[0];
2491 Int_t frequency = outlab[1];
2492 for (Int_t iLayer = 0; iLayer < 6; iLayer++) {
2493 cseed[iLayer].SetFreq(frequency);
2494 cseed[iLayer].SetChi2Z(chi2[1]);
2497 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2498 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2499 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2500 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2501 TLinearFitter *fitterTC = GetTiltedRiemanFitterConstraint();
2502 TLinearFitter *fitterT = GetTiltedRiemanFitter();
2503 cstreamer << "MakeSeeds2"
2504 << "EventNumber=" << eventNumber
2505 << "CandidateNumber=" << candidateNumber
2506 << "Chi2TR=" << chi2Vals[0]
2507 << "Chi2TC=" << chi2Vals[1]
2508 << "Nlayers=" << mlayers
2509 << "NUsedS=" << nUsedCl
2510 << "NUsed=" << nusedf
2512 << "S0.=" << &cseed[0]
2513 << "S1.=" << &cseed[1]
2514 << "S2.=" << &cseed[2]
2515 << "S3.=" << &cseed[3]
2516 << "S4.=" << &cseed[4]
2517 << "S5.=" << &cseed[5]
2518 << "Label=" << label
2519 << "Freq=" << frequency
2520 << "FitterT.=" << fitterT
2521 << "FitterTC.=" << fitterTC
2526 AliTRDtrackerDebug::SetCandidateNumber(AliTRDtrackerDebug::GetCandidateNumber() + 1);
2527 if(ntracks == kMaxTracksStack){
2528 AliWarning(Form("Number of seeds reached maximum allowed (%d) in stack.", kMaxTracksStack));
2539 //_____________________________________________________________________________
2540 AliTRDtrackV1* AliTRDtrackerV1::MakeTrack(AliTRDseedV1 *seeds, Double_t *params)
2543 // Build a TRD track out of tracklet candidates
2546 // seeds : array of tracklets
2547 // params : track parameters (see MakeSeeds() function body for a detailed description)
2552 // Detailed description
2554 // To be discussed with Marian !!
2557 AliTRDCalibraFillHisto *calibra = AliTRDCalibraFillHisto::Instance();
2558 if (!calibra) AliInfo("Could not get Calibra instance\n");
2560 Double_t alpha = AliTRDgeometry::GetAlpha();
2561 Double_t shift = AliTRDgeometry::GetAlpha()/2.0;
2565 c[ 1] = 0.0; c[ 2] = 2.0;
2566 c[ 3] = 0.0; c[ 4] = 0.0; c[ 5] = 0.02;
2567 c[ 6] = 0.0; c[ 7] = 0.0; c[ 8] = 0.0; c[ 9] = 0.1;
2568 c[10] = 0.0; c[11] = 0.0; c[12] = 0.0; c[13] = 0.0; c[14] = params[5]*params[5]*0.01;
2570 AliTRDtrackV1 track(seeds, ¶ms[1], c, params[0], params[6]*alpha+shift);
2571 track.PropagateTo(params[0]-5.0);
2572 track.ResetCovariance(1);
2573 Int_t nc = TMath::Abs(FollowBackProlongation(track));
2574 if (nc < 30) return 0x0;
2576 AliTRDtrackV1 *ptrTrack = SetTrack(&track);
2577 ptrTrack->CookLabel(.9);
2578 // computes PID for track
2579 ptrTrack->CookPID();
2580 // update calibration references using this track
2581 if(calibra->GetHisto2d()) calibra->UpdateHistogramsV1(ptrTrack);
2587 //____________________________________________________________________
2588 Int_t AliTRDtrackerV1::ImproveSeedQuality(AliTRDtrackingChamber **stack, AliTRDseedV1 *cseed)
2591 // Sort tracklets according to "quality" and try to "improve" the first 4 worst
2594 // layers : Array of propagation layers for a stack/supermodule
2595 // cseed : Array of 6 seeding tracklets which has to be improved
2598 // cssed : Improved seeds
2600 // Detailed description
2602 // Iterative procedure in which new clusters are searched for each
2603 // tracklet seed such that the seed quality (see AliTRDseed::GetQuality())
2604 // can be maximized. If some optimization is found the old seeds are replaced.
2609 // make a local working copy
2610 AliTRDtrackingChamber *chamber = 0x0;
2611 AliTRDseedV1 bseed[6];
2613 for (Int_t jLayer = 0; jLayer < 6; jLayer++) bseed[jLayer] = cseed[jLayer];
2615 Float_t lastquality = 10000.0;
2616 Float_t lastchi2 = 10000.0;
2617 Float_t chi2 = 1000.0;
2619 for (Int_t iter = 0; iter < 4; iter++) {
2620 Float_t sumquality = 0.0;
2621 Float_t squality[6];
2622 Int_t sortindexes[6];
2624 for (Int_t jLayer = 0; jLayer < 6; jLayer++) {
2625 squality[jLayer] = bseed[jLayer].IsOK() ? bseed[jLayer].GetQuality(kTRUE) : 1000.;
2626 sumquality += squality[jLayer];
2628 if ((sumquality >= lastquality) || (chi2 > lastchi2)) break;
2631 lastquality = sumquality;
2633 if (iter > 0) for (Int_t jLayer = 0; jLayer < 6; jLayer++) cseed[jLayer] = bseed[jLayer];
2635 TMath::Sort(6, squality, sortindexes, kFALSE);
2636 for (Int_t jLayer = 5; jLayer > 1; jLayer--) {
2637 Int_t bLayer = sortindexes[jLayer];
2638 if(!(chamber = stack[bLayer])) continue;
2639 bseed[bLayer].AttachClustersIter(chamber, squality[bLayer], kTRUE);
2640 if(bseed[bLayer].IsOK()) nLayers++;
2643 chi2 = FitTiltedRieman(bseed, kTRUE);
2644 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 7){
2645 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2646 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2647 TLinearFitter *tiltedRieman = GetTiltedRiemanFitter();
2648 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2649 cstreamer << "ImproveSeedQuality"
2650 << "EventNumber=" << eventNumber
2651 << "CandidateNumber=" << candidateNumber
2652 << "Iteration=" << iter
2653 << "S0.=" << &bseed[0]
2654 << "S1.=" << &bseed[1]
2655 << "S2.=" << &bseed[2]
2656 << "S3.=" << &bseed[3]
2657 << "S4.=" << &bseed[4]
2658 << "S5.=" << &bseed[5]
2659 << "FitterT.=" << tiltedRieman
2664 // we are sure that at least 2 tracklets are OK !
2668 //_________________________________________________________________________
2669 Double_t AliTRDtrackerV1::CalculateTrackLikelihood(AliTRDseedV1 *tracklets, Double_t *chi2){
2671 // Calculates the Track Likelihood value. This parameter serves as main quality criterion for
2672 // the track selection
2673 // The likelihood value containes:
2674 // - The chi2 values from the both fitters and the chi2 values in z-direction from a linear fit
2675 // - The Sum of the Parameter |slope_ref - slope_fit|/Sigma of the tracklets
2676 // For all Parameters an exponential dependency is used
2678 // Parameters: - Array of tracklets (AliTRDseedV1) related to the track candidate
2679 // - Array of chi2 values:
2680 // * Non-Constrained Tilted Riemann fit
2681 // * Vertex-Constrained Tilted Riemann fit
2682 // * z-Direction from Linear fit
2683 // Output: - The calculated track likelihood
2688 Double_t sumdaf = 0, nLayers = 0;
2689 for (Int_t iLayer = 0; iLayer < kNPlanes; iLayer++) {
2690 if(!tracklets[iLayer].IsOK()) continue;
2691 sumdaf += TMath::Abs((tracklets[iLayer].GetYfit(1) - tracklets[iLayer].GetYref(1))/ tracklets[iLayer].GetSigmaY2());
2694 sumdaf /= Float_t (nLayers - 2.0);
2696 Double_t likeChi2Z = TMath::Exp(-chi2[2] * 0.14); // Chi2Z
2697 Double_t likeChi2TC = (fReconstructor->GetRecoParam() ->IsVertexConstrained()) ?
2698 TMath::Exp(-chi2[1] * 0.677) : 1; // Constrained Tilted Riemann
2699 Double_t likeChi2TR = TMath::Exp(-chi2[0] * 0.78); // Non-constrained Tilted Riemann
2700 Double_t likeAF = TMath::Exp(-sumdaf * 3.23);
2701 Double_t trackLikelihood = likeChi2Z * likeChi2TR * likeAF;
2703 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2704 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2705 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2706 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2707 cstreamer << "CalculateTrackLikelihood0"
2708 << "EventNumber=" << eventNumber
2709 << "CandidateNumber=" << candidateNumber
2710 << "LikeChi2Z=" << likeChi2Z
2711 << "LikeChi2TR=" << likeChi2TR
2712 << "LikeChi2TC=" << likeChi2TC
2713 << "LikeAF=" << likeAF
2714 << "TrackLikelihood=" << trackLikelihood
2718 return trackLikelihood;
2721 //____________________________________________________________________
2722 Double_t AliTRDtrackerV1::CookLikelihood(AliTRDseedV1 *cseed, Int_t planes[4]
2726 // Calculate the probability of this track candidate.
2729 // cseeds : array of candidate tracklets
2730 // planes : array of seeding planes (see seeding configuration)
2731 // chi2 : chi2 values (on the Z and Y direction) from the rieman fit of the track.
2736 // Detailed description
2738 // The track quality is estimated based on the following 4 criteria:
2739 // 1. precision of the rieman fit on the Y direction (likea)
2740 // 2. chi2 on the Y direction (likechi2y)
2741 // 3. chi2 on the Z direction (likechi2z)
2742 // 4. number of attached clusters compared to a reference value
2743 // (see AliTRDrecoParam::fkFindable) (likeN)
2745 // The distributions for each type of probabilities are given below as of
2746 // (date). They have to be checked to assure consistency of estimation.
2749 // ratio of the total number of clusters/track which are expected to be found by the tracker.
2750 Float_t fgFindable = fReconstructor->GetRecoParam() ->GetFindableClusters();
2753 Int_t nclusters = 0;
2754 Double_t sumda = 0.;
2755 for(UChar_t ilayer = 0; ilayer < 4; ilayer++){
2756 Int_t jlayer = planes[ilayer];
2757 nclusters += cseed[jlayer].GetN2();
2758 sumda += TMath::Abs(cseed[jlayer].GetYfitR(1) - cseed[jlayer].GetYref(1));
2760 Double_t likea = TMath::Exp(-sumda*10.6);
2761 Double_t likechi2y = 0.0000000001;
2762 if (chi2[0] < 0.5) likechi2y += TMath::Exp(-TMath::Sqrt(chi2[0]) * 7.73);
2763 Double_t likechi2z = TMath::Exp(-chi2[1] * 0.088) / TMath::Exp(-chi2[1] * 0.019);
2764 Int_t enc = Int_t(fgFindable*4.*fgNTimeBins); // Expected Number Of Clusters, normally 72
2765 Double_t likeN = TMath::Exp(-(enc - nclusters) * 0.19);
2767 Double_t like = likea * likechi2y * likechi2z * likeN;
2769 // 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));
2770 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) >= 2){
2771 Int_t eventNumber = AliTRDtrackerDebug::GetEventNumber();
2772 Int_t candidateNumber = AliTRDtrackerDebug::GetCandidateNumber();
2773 // The Debug Stream contains the seed
2774 TTreeSRedirector &cstreamer = *fgDebugStreamer;
2775 cstreamer << "CookLikelihood"
2776 << "EventNumber=" << eventNumber
2777 << "CandidateNumber=" << candidateNumber
2778 << "tracklet0.=" << &cseed[0]
2779 << "tracklet1.=" << &cseed[1]
2780 << "tracklet2.=" << &cseed[2]
2781 << "tracklet3.=" << &cseed[3]
2782 << "tracklet4.=" << &cseed[4]
2783 << "tracklet5.=" << &cseed[5]
2784 << "sumda=" << sumda
2785 << "chi0=" << chi2[0]
2786 << "chi1=" << chi2[1]
2787 << "likea=" << likea
2788 << "likechi2y=" << likechi2y
2789 << "likechi2z=" << likechi2z
2790 << "nclusters=" << nclusters
2791 << "likeN=" << likeN
2801 //____________________________________________________________________
2802 void AliTRDtrackerV1::GetSeedingConfig(Int_t iconfig, Int_t planes[4])
2805 // Map seeding configurations to detector planes.
2808 // iconfig : configuration index
2809 // planes : member planes of this configuration. On input empty.
2812 // planes : contains the planes which are defining the configuration
2814 // Detailed description
2816 // Here is the list of seeding planes configurations together with
2817 // their topological classification:
2835 // The topologic quality is modeled as follows:
2836 // 1. The general model is define by the equation:
2837 // p(conf) = exp(-conf/2)
2838 // 2. According to the topologic classification, configurations from the same
2839 // class are assigned the agerage value over the model values.
2840 // 3. Quality values are normalized.
2842 // The topologic quality distribution as function of configuration is given below:
2844 // <img src="gif/topologicQA.gif">
2849 case 0: // 5432 TQ 0
2855 case 1: // 4321 TQ 0
2861 case 2: // 3210 TQ 0
2867 case 3: // 5321 TQ 1
2873 case 4: // 4210 TQ 1
2879 case 5: // 5431 TQ 1
2885 case 6: // 4320 TQ 1
2891 case 7: // 5430 TQ 2
2897 case 8: // 5210 TQ 2
2903 case 9: // 5421 TQ 3
2909 case 10: // 4310 TQ 3
2915 case 11: // 5410 TQ 4
2921 case 12: // 5420 TQ 5
2927 case 13: // 5320 TQ 5
2933 case 14: // 5310 TQ 5
2942 //____________________________________________________________________
2943 void AliTRDtrackerV1::GetExtrapolationConfig(Int_t iconfig, Int_t planes[2])
2946 // Returns the extrapolation planes for a seeding configuration.
2949 // iconfig : configuration index
2950 // planes : planes which are not in this configuration. On input empty.
2953 // planes : contains the planes which are not in the configuration
2955 // Detailed description
2959 case 0: // 5432 TQ 0
2963 case 1: // 4321 TQ 0
2967 case 2: // 3210 TQ 0
2971 case 3: // 5321 TQ 1
2975 case 4: // 4210 TQ 1
2979 case 5: // 5431 TQ 1
2983 case 6: // 4320 TQ 1
2987 case 7: // 5430 TQ 2
2991 case 8: // 5210 TQ 2
2995 case 9: // 5421 TQ 3
2999 case 10: // 4310 TQ 3
3003 case 11: // 5410 TQ 4
3007 case 12: // 5420 TQ 5
3011 case 13: // 5320 TQ 5
3015 case 14: // 5310 TQ 5
3022 //____________________________________________________________________
3023 AliCluster* AliTRDtrackerV1::GetCluster(Int_t idx) const
3025 Int_t ncls = fClusters->GetEntriesFast();
3026 return idx >= 0 && idx < ncls ? (AliCluster*)fClusters->UncheckedAt(idx) : 0x0;
3029 //____________________________________________________________________
3030 AliTRDseedV1* AliTRDtrackerV1::GetTracklet(Int_t idx) const
3032 Int_t ntrklt = fTracklets->GetEntriesFast();
3033 return idx >= 0 && idx < ntrklt ? (AliTRDseedV1*)fTracklets->UncheckedAt(idx) : 0x0;
3036 //____________________________________________________________________
3037 AliKalmanTrack* AliTRDtrackerV1::GetTrack(Int_t idx) const
3039 Int_t ntrk = fTracks->GetEntriesFast();
3040 return idx >= 0 && idx < ntrk ? (AliKalmanTrack*)fTracks->UncheckedAt(idx) : 0x0;
3043 //____________________________________________________________________
3044 Float_t AliTRDtrackerV1::CalculateReferenceX(AliTRDseedV1 *tracklets){
3046 // Calculates the reference x-position for the tilted Rieman fit defined as middle
3047 // of the stack (middle between layers 2 and 3). For the calculation all the tracklets
3048 // are taken into account
3050 // Parameters: - Array of tracklets(AliTRDseedV1)
3052 // Output: - The reference x-position(Float_t)
3054 Int_t nDistances = 0;
3055 Float_t meanDistance = 0.;
3056 Int_t startIndex = 5;
3057 for(Int_t il =5; il > 0; il--){
3058 if(tracklets[il].IsOK() && tracklets[il -1].IsOK()){
3059 Float_t xdiff = tracklets[il].GetX0() - tracklets[il -1].GetX0();
3060 meanDistance += xdiff;
3063 if(tracklets[il].IsOK()) startIndex = il;
3065 if(tracklets[0].IsOK()) startIndex = 0;
3067 // We should normally never get here
3068 Float_t xpos[2]; memset(xpos, 0, sizeof(Float_t) * 2);
3069 Int_t iok = 0, idiff = 0;
3070 // This attempt is worse and should be avoided:
3071 // check for two chambers which are OK and repeat this without taking the mean value
3072 // Strategy avoids a division by 0;
3073 for(Int_t il = 5; il >= 0; il--){
3074 if(tracklets[il].IsOK()){
3075 xpos[iok] = tracklets[il].GetX0();
3079 if(iok) idiff++; // to get the right difference;
3083 meanDistance = (xpos[0] - xpos[1])/idiff;
3086 // we have do not even have 2 layers which are OK? The we do not need to fit at all
3091 meanDistance /= nDistances;
3093 return tracklets[startIndex].GetX0() + (2.5 - startIndex) * meanDistance - 0.5 * (AliTRDgeometry::AmThick() + AliTRDgeometry::DrThick());
3096 //_____________________________________________________________________________
3097 Int_t AliTRDtrackerV1::Freq(Int_t n, const Int_t *inlist
3098 , Int_t *outlist, Bool_t down)
3101 // Sort eleements according occurancy
3102 // The size of output array has is 2*n
3109 Int_t *sindexS = new Int_t[n]; // Temporary array for sorting
3110 Int_t *sindexF = new Int_t[2*n];
3111 for (Int_t i = 0; i < n; i++) {
3115 TMath::Sort(n,inlist,sindexS,down);
3117 Int_t last = inlist[sindexS[0]];
3120 sindexF[0+n] = last;
3124 for (Int_t i = 1; i < n; i++) {
3125 val = inlist[sindexS[i]];
3127 sindexF[countPos]++;
3131 sindexF[countPos+n] = val;
3132 sindexF[countPos]++;
3140 // Sort according frequency
3141 TMath::Sort(countPos,sindexF,sindexS,kTRUE);
3143 for (Int_t i = 0; i < countPos; i++) {
3144 outlist[2*i ] = sindexF[sindexS[i]+n];
3145 outlist[2*i+1] = sindexF[sindexS[i]];
3155 void AliTRDtrackerV1::SetReconstructor(const AliTRDReconstructor *rec){
3156 fReconstructor = rec;
3157 if(fReconstructor->GetStreamLevel(AliTRDReconstructor::kTracker) > 1){
3158 if(!fgDebugStreamer)
3159 fgDebugStreamer = new TTreeSRedirector("TRD.TrackerDebug.root");
3163 //_____________________________________________________________________________
3164 Float_t AliTRDtrackerV1::GetChi2Y(AliTRDseedV1 *tracklets) const
3166 // Chi2 definition on y-direction
3169 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
3170 if(!tracklets[ipl].IsOK()) continue;
3171 Double_t distLayer = tracklets[ipl].GetYfit(0) - tracklets[ipl].GetYref(0);
3172 chi2 += distLayer * distLayer;
3177 //____________________________________________________________________
3178 void AliTRDtrackerV1::ResetSeedTB()
3180 // reset buffer for seeding time bin layers. If the time bin
3181 // layers are not allocated this function allocates them
3183 for(Int_t isl=0; isl<kNSeedPlanes; isl++){
3184 if(!fSeedTB[isl]) fSeedTB[isl] = new AliTRDchamberTimeBin();
3185 else fSeedTB[isl]->Clear();
3189 //_____________________________________________________________________________
3190 Float_t AliTRDtrackerV1::GetChi2Z(AliTRDseedV1 *tracklets) const
3192 // Chi2 definition on z-direction
3195 for(Int_t ipl = 0; ipl < kNPlanes; ipl++){
3196 if(!tracklets[ipl].IsOK()) continue;
3197 Double_t distLayer = tracklets[ipl].GetMeanz() - tracklets[ipl].GetZref(0);
3198 chi2 += distLayer * distLayer;
3203 ///////////////////////////////////////////////////////
3205 // Resources of class AliTRDLeastSquare //
3207 ///////////////////////////////////////////////////////
3209 //_____________________________________________________________________________
3210 AliTRDtrackerV1::AliTRDLeastSquare::AliTRDLeastSquare(){
3212 // Constructor of the nested class AliTRDtrackFitterLeastSquare
3214 memset(fParams, 0, sizeof(Double_t) * 2);
3215 memset(fSums, 0, sizeof(Double_t) * 5);
3216 memset(fCovarianceMatrix, 0, sizeof(Double_t) * 3);
3220 //_____________________________________________________________________________
3221 void AliTRDtrackerV1::AliTRDLeastSquare::AddPoint(Double_t *x, Double_t y, Double_t sigmaY){
3223 // Adding Point to the fitter
3225 Double_t weight = 1/(sigmaY * sigmaY);
3227 // printf("Adding point x = %f, y = %f, sigma = %f\n", xpt, y, sigmaY);
3229 fSums[1] += weight * xpt;
3230 fSums[2] += weight * y;
3231 fSums[3] += weight * xpt * y;
3232 fSums[4] += weight * xpt * xpt;
3233 fSums[5] += weight * y * y;
3236 //_____________________________________________________________________________
3237 void AliTRDtrackerV1::AliTRDLeastSquare::RemovePoint(Double_t *x, Double_t y, Double_t sigmaY){
3239 // Remove Point from the sample
3241 Double_t weight = 1/(sigmaY * sigmaY);
3244 fSums[1] -= weight * xpt;
3245 fSums[2] -= weight * y;
3246 fSums[3] -= weight * xpt * y;
3247 fSums[4] -= weight * xpt * xpt;
3248 fSums[5] -= weight * y * y;
3251 //_____________________________________________________________________________
3252 void AliTRDtrackerV1::AliTRDLeastSquare::Eval(){
3254 // Evaluation of the fit:
3255 // Calculation of the parameters
3256 // Calculation of the covariance matrix
3259 Double_t denominator = fSums[0] * fSums[4] - fSums[1] *fSums[1];
3260 if(denominator==0) return;
3262 // for(Int_t isum = 0; isum < 5; isum++)
3263 // printf("fSums[%d] = %f\n", isum, fSums[isum]);
3264 // printf("denominator = %f\n", denominator);
3265 fParams[0] = (fSums[2] * fSums[4] - fSums[1] * fSums[3])/ denominator;
3266 fParams[1] = (fSums[0] * fSums[3] - fSums[1] * fSums[2]) / denominator;
3267 // printf("fParams[0] = %f, fParams[1] = %f\n", fParams[0], fParams[1]);
3269 // Covariance matrix
3270 fCovarianceMatrix[0] = fSums[4] - fSums[1] * fSums[1] / fSums[0];
3271 fCovarianceMatrix[1] = fSums[5] - fSums[2] * fSums[2] / fSums[0];
3272 fCovarianceMatrix[2] = fSums[3] - fSums[1] * fSums[2] / fSums[0];
3275 //_____________________________________________________________________________
3276 Double_t AliTRDtrackerV1::AliTRDLeastSquare::GetFunctionValue(Double_t *xpos) const {
3278 // Returns the Function value of the fitted function at a given x-position
3280 return fParams[0] + fParams[1] * (*xpos);
3283 //_____________________________________________________________________________
3284 void AliTRDtrackerV1::AliTRDLeastSquare::GetCovarianceMatrix(Double_t *storage) const {
3286 // Copies the values of the covariance matrix into the storage
3288 memcpy(storage, fCovarianceMatrix, sizeof(Double_t) * 3);