]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALTracker.cxx
1b3ca8d00479d6e8663ca0ffadc4e2a5acf2a10a
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALTracker.cxx
1 //========================================================================
2 // Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. 
3 //                                                                        
4 // Author: The ALICE Off-line Project.                                    
5 // Contributors are mentioned in the code where appropriate.              
6 //                                                                        
7 // Permission to use, copy, modify and distribute this software and its   
8 // documentation strictly for non-commercial purposes is hereby granted   
9 // without fee, provided that the above copyright notice appears in all   
10 // copies and that both the copyright notice and this permission notice   
11 // appear in the supporting documentation. The authors make no claims     
12 // about the suitability of this software for any purpose. It is          
13 // provided "as is" without express or implied warranty.                  
14 //======================================================================== 
15 //                       
16 //                       Class AliEMCALTracker 
17 //                      -----------------------
18 // Implementation of the track matching method between barrel tracks and
19 // EMCAL clusters.
20 // Besides algorithm implementation, some cuts are required to be set
21 // in order to define, for each track, an acceptance window where clusters
22 // are searched to find best match (if any).
23 // The class accepts as input an ESD container, and works directly on it,
24 // simply setting, for each of its tracks, the fEMCALindex flag, for each
25 // track which is matched to a cluster.
26 // In order to use method, one must launch PropagateBack().
27 //
28 // ------------------------------------------------------------------------
29 // author: A. Pulvirenti (alberto.pulvirenti@ct.infn.it)
30 //=========================================================================
31
32 #include <Riostream.h>
33 #include <iomanip>
34
35 #include <TFile.h>
36 #include <TTree.h>
37 #include <TList.h>
38 #include <TString.h>
39 #include <TVector3.h>
40 #include <TClonesArray.h>
41 #include <TGeoMatrix.h>
42
43 #include "AliLog.h"
44 #include "AliESDEvent.h"
45 #include "AliESDtrack.h"
46 #include "AliESDCaloCluster.h"
47 #include "AliKalmanTrack.h"
48 #include "AliEMCALRecPoint.h"
49 #include "AliRunLoader.h"
50 #include "AliEMCALTrack.h"
51 #include "AliEMCALLoader.h"
52 #include "AliEMCALGeometry.h"
53
54 #include "AliEMCALTracker.h"
55
56 ClassImp(AliEMCALTracker)
57 //
58 //------------------------------------------------------------------------------
59 //
60 AliEMCALTracker::AliEMCALTracker() 
61   : AliTracker(),
62     fNPropSteps(0),
63     fTrackCorrMode(kTrackCorrNone),
64     fCutX(50.0),
65     fCutY(50.0),
66     fCutZ(50.0),
67     fCutAlphaMin(-200.0),
68     fCutAlphaMax(200.0),
69     fCutAngle(100.0),
70     fMaxDist(100.0),
71     fRho(1.0),
72     fX0(1.0),
73     fTracks(0),
74     fClusters(0),
75     fMatches(0),
76     fGeom(0)
77 {
78         //
79         // Default constructor.
80         // Initializes al simple data members to default values,
81         // and all collections to NULL.
82         // Output file name is set to a default value.
83         //
84 }
85 //
86 //------------------------------------------------------------------------------
87 //
88 AliEMCALTracker::AliEMCALTracker(const AliEMCALTracker& copy) 
89   : AliTracker(),
90     fNPropSteps(copy.fNPropSteps),
91     fTrackCorrMode(copy.fTrackCorrMode),
92     fCutX(copy.fCutX),
93     fCutY(copy.fCutY),
94     fCutZ(copy.fCutZ),
95     fCutAlphaMin(copy.fCutAlphaMin),
96     fCutAlphaMax(copy.fCutAlphaMax),
97     fCutAngle(copy.fCutAngle),
98     fMaxDist(copy.fMaxDist),
99     fRho(copy.fRho),
100     fX0(copy.fX0),
101     fTracks((TObjArray*)copy.fTracks->Clone()),
102     fClusters((TObjArray*)copy.fClusters->Clone()),
103     fMatches((TList*)copy.fMatches->Clone()),
104     fGeom(copy.fGeom)
105 {
106         //
107         // Copy constructor
108         // Besides copying all parameters, duplicates all collections.
109         //
110 }
111 //
112 //------------------------------------------------------------------------------
113 //
114 AliEMCALTracker& AliEMCALTracker::operator=(const AliEMCALTracker& copy)
115 {
116         //
117         // Assignment operator.
118         // Besides copying all parameters, duplicates all collections.  
119         //
120         
121         fCutX = copy.fCutX;
122         fCutY = copy.fCutY;
123         fCutZ = copy.fCutZ;
124         fCutAlphaMin = copy.fCutAlphaMin;
125         fCutAlphaMax = copy.fCutAlphaMax;
126         fCutAngle = copy.fCutAngle;
127         fMaxDist = copy.fMaxDist;
128         
129         fTracks = (TObjArray*)copy.fTracks->Clone();
130         fClusters = (TObjArray*)copy.fClusters->Clone();
131         fMatches = (TList*)copy.fMatches->Clone();
132         
133         fGeom = copy.fGeom;
134         
135         return (*this);
136 }
137 //
138 //------------------------------------------------------------------------------
139 //
140 TTree* AliEMCALTracker::SearchTrueMatches()
141 {
142         if (!fClusters) return 0;
143         if (fClusters->IsEmpty()) return 0;
144         if (!fTracks) return 0;
145         if (fTracks->IsEmpty()) return 0;
146         
147         TTree *outTree = new TTree("tree", "True matches from event");
148         Int_t indexT, indexC, label;
149         outTree->Branch("indexC", &indexC, "indexC/I");
150         outTree->Branch("indexT", &indexT, "indexT/I");
151         outTree->Branch("label",  &label , "label/I");
152         
153         Double_t dist;
154         Int_t ic, nClusters = (Int_t)fClusters->GetEntries();
155         Int_t it, nTracks = fTracks->GetEntries();
156         
157         for (ic = 0; ic < nClusters; ic++) {
158                 AliEMCALMatchCluster *cluster = (AliEMCALMatchCluster*)fClusters->At(ic);
159                 label = cluster->Label();
160                 indexC = cluster->Index();
161                 for (it = 0; it < nTracks; it++) {
162                         AliEMCALTrack *track = (AliEMCALTrack*)fTracks->At(it);
163                         if (TMath::Abs(track->GetSeedLabel()) != label) continue;
164                         dist = CheckPair(track, cluster);
165                         if (dist <= fMaxDist) {
166                                 indexT = track->GetSeedIndex();
167                                 outTree->Fill();
168                         }
169                 }
170         }
171         
172         return outTree;
173 }
174 //
175 //------------------------------------------------------------------------------
176 //
177 void AliEMCALTracker::Clear(Option_t* option)
178 {
179         //
180         // Clearing method
181         // Clears all specified arrays and the containers themselves.
182         //
183
184         TString opt(option);
185         Bool_t doDelete = opt.Contains("DELETE");
186         Bool_t clearTracks = opt.Contains("TRACKS");
187         Bool_t clearClusters = opt.Contains("CLUSTERS");
188         Bool_t clearMatches = opt.Contains("MATCHES");
189         if (opt.Contains("ALL")) {
190                 clearTracks = kTRUE;
191                 clearClusters = kTRUE;
192                 clearMatches = kTRUE;
193         }
194         
195         if (fTracks != 0x0 && clearTracks) {
196                 if (!doDelete) {
197                         fTracks->Clear();
198                 }
199                 else {
200                         if (!fTracks->IsEmpty()) fTracks->Delete();
201                         delete fTracks;
202                         fTracks = 0;
203                 }
204         }
205         if (fClusters != 0x0 && clearClusters) {
206                 if (!doDelete) {
207                         fClusters->Clear();
208                 }
209                 else {
210                         if (!fClusters->IsEmpty()) fClusters->Delete();
211                         delete fClusters;
212                         fClusters = 0;
213                 }
214         }
215         if (fMatches != 0x0 && clearMatches) {
216                 if (!doDelete) {
217                         fMatches->Clear();
218                 }
219                 else {
220                         if (!fMatches->IsEmpty()) fMatches->Delete();
221                         delete fMatches;
222                         fMatches = 0;
223                 }
224         }
225 }
226 //
227 //------------------------------------------------------------------------------
228 //
229 Int_t AliEMCALTracker::LoadClusters(TTree *cTree) 
230 {
231         //
232         // Load EMCAL clusters in the form of AliEMCALRecPoint,
233         // from simulation temporary files.
234         // (When included in reconstruction chain, this method is used automatically)
235         //
236         
237         Clear("CLUSTERS");
238
239         TBranch *branch = cTree->GetBranch("EMCALECARP");
240         if (!branch) {
241                 AliError("Can't get the branch with the EMCAL clusters");
242                 return 1;
243         }
244         
245         TClonesArray *clusters = new TClonesArray("AliEMCALRecPoint", 1000);
246         branch->SetAddress(&clusters);
247         
248         cTree->GetEvent(0);
249         Int_t nClusters = (Int_t)clusters->GetEntries();
250         fClusters = new TObjArray(0);
251         for (Int_t i = 0; i < nClusters; i++) {
252                 AliEMCALRecPoint *cluster = (AliEMCALRecPoint*)clusters->At(i);
253                 if (!cluster) continue;
254                 if (cluster->GetClusterType() != AliESDCaloCluster::kClusterv1) continue;
255                 AliEMCALMatchCluster *matchCluster = new AliEMCALMatchCluster(i, cluster);
256                 fClusters->AddLast(matchCluster);
257         }
258         if (fClusters->IsEmpty()) {
259                 AliError("No clusters collected");
260                 return 1;
261         }
262         
263         branch->SetAddress(0);
264         clusters->Delete();
265         delete clusters;
266
267         AliInfo(Form("Collected %d clusters", fClusters->GetEntries()));
268
269         return 0;
270 }
271 //
272 //------------------------------------------------------------------------------
273 //
274 Int_t AliEMCALTracker::LoadClusters(AliESDEvent *esd) 
275 {
276         //
277         // Load EMCAL clusters in the form of AliESDCaloClusters,
278         // from an AliESD object.
279         //
280
281         // make sure that tracks/clusters collections are empty
282         Clear("CLUSTERS");
283         
284         Int_t start = esd->GetFirstEMCALCluster();
285         Int_t nClustersEMC = esd->GetNumberOfEMCALClusters();
286         Int_t end = start + nClustersEMC;
287         
288         fClusters = new TObjArray(0);
289                 
290         Int_t i;
291         for (i = start; i < end; i++) {
292                 AliESDCaloCluster *cluster = esd->GetCaloCluster(i);
293                 if (!cluster) continue;
294                 if (cluster->GetClusterType() != AliESDCaloCluster::kClusterv1) continue;
295                 AliEMCALMatchCluster *matchCluster = new AliEMCALMatchCluster(i, cluster);
296                 fClusters->AddLast(matchCluster);
297         }
298         if (fClusters->IsEmpty()) {
299                 AliError("No clusters collected");
300                 return 1;
301         }
302         
303         AliInfo(Form("Collected %d clusters", fClusters->GetEntries()));
304
305         return 0;
306 }
307 //
308 //------------------------------------------------------------------------------
309 //
310 Int_t AliEMCALTracker::LoadTracks(AliESDEvent *esd)
311 {
312         //
313         // Load ESD tracks.
314         //
315         
316         Clear("TRACKS");
317         
318         Int_t nTracks = esd->GetNumberOfTracks();
319         fTracks = new TObjArray(0);
320         
321         Int_t i, j;
322         Bool_t isKink;
323         Double_t alpha; 
324         for (i = 0; i < nTracks; i++) {
325                 AliESDtrack *esdTrack = esd->GetTrack(i);
326                 // set by default the value corresponding to "no match"
327                 esdTrack->SetEMCALcluster(kUnmatched);
328 //              if (esdTrack->GetLabel() < 0) continue;
329 //              if (!(esdTrack->GetStatus() & AliESDtrack::kTOFout)) continue;
330                 isKink = kFALSE;
331                 for (j = 0; j < 3; j++) {
332                         if (esdTrack->GetKinkIndex(j) != 0) isKink = kTRUE;
333                 }
334                 if (isKink) continue;
335                 AliEMCALTrack *track = new AliEMCALTrack(*esdTrack);
336                 track->SetMass(0.13957018);
337                 // check alpha and reject the tracks which fall outside EMCAL acceptance
338                 alpha = track->GetAlpha() * TMath::RadToDeg();
339                 if (alpha >  -155.0 && alpha < 67.0) {
340                         delete track;
341                         continue;
342                 }
343 //              if (!PropagateToEMCAL(track)) {
344 //                      delete track;
345 //                      continue;
346 //              }
347                 track->SetSeedIndex(i);
348                 track->SetSeedLabel(esdTrack->GetLabel());
349                 fTracks->AddLast(track);
350         }
351         if (fTracks->IsEmpty()) {
352                 AliError("No tracks collected");
353                 return 1;
354         }
355         
356         AliInfo(Form("Collected %d tracks", fTracks->GetEntries()));
357
358         return 0;
359 }
360 //
361 //------------------------------------------------------------------------------
362 //
363 Int_t AliEMCALTracker::PropagateBack(AliESDEvent* esd)
364 {
365         //
366         // Main operation method.
367         // Gets external AliESD containing tracks to be matched.
368         // After executing match finding, stores in the same ESD object all infos
369         // and releases the object for further reconstruction steps.
370         //
371         
372         if (!esd) {
373                 AliError("NULL ESD passed");
374                 return 1;
375         }
376         
377         // step 1: 
378         // if cluster array is empty, cluster are collected
379         // from the passed ESD, and work is done with ESDCaloClusters
380         Int_t okLoadClusters, nClusters;
381         if (!fClusters || (fClusters && fClusters->IsEmpty())) {
382                 okLoadClusters = LoadClusters(esd);
383                 if (okLoadClusters) return 2;
384         }
385         nClusters = fClusters->GetEntries();
386         
387         // step 2:
388         // collect ESD tracks
389         Int_t okLoadTracks = LoadTracks(esd), nTracks;
390         if (okLoadTracks) return 3;
391         nTracks = fTracks->GetEntries();
392         
393         // step 3:
394         // each track is propagated to the "R" position of each cluster.
395         // The closest cluster is assigned as match.
396         // IF no clusters lie within the maximum allowed distance, no matches are assigned.
397         Int_t nMatches = CreateMatches();
398         if (!nMatches) {
399                 AliInfo(Form("#clusters = %d -- #tracks = %d --> No good matches found.", nClusters, nTracks));
400                 return 4;
401         }
402         else {
403                 AliInfo(Form("#clusters = %d -- #tracks = %d --> Found %d matches.", nClusters, nTracks, nMatches));
404         }
405         
406         // step 4:
407         // when more than 1 track share the same matched cluster, only the closest one is kept.
408         Int_t nRemoved = SolveCompetitions();
409         AliInfo(Form("Removed %d duplicate matches", nRemoved));
410         if (nRemoved >= nMatches) {
411                 AliError("Removed ALL matches! Check the algorithm or data. Nothing to save");
412                 return 5;
413         }
414         
415         // step 5:
416         // save obtained information setting the 'fEMCALindex' field of AliESDtrack object
417         Int_t nSaved = 0, trackID, nGood = 0, nFake = 0;
418         TListIter iter(fMatches);
419         AliEMCALMatch *match = 0;
420         while ( (match = (AliEMCALMatch*)iter.Next()) ) {
421                 if (!match->CanBeSaved()) continue;
422                 AliEMCALTrack *track = (AliEMCALTrack*)fTracks->At(match->GetIndexT());
423                 AliEMCALMatchCluster *cluster = (AliEMCALMatchCluster*)fClusters->At(match->GetIndexC());
424                 trackID = track->GetSeedIndex();
425                 AliESDtrack *esdTrack = esd->GetTrack(trackID);
426                 if (!esdTrack) continue;
427                 if (TMath::Abs(esdTrack->GetLabel()) == cluster->Label()) {
428                         esdTrack->SetEMCALcluster(cluster->Index());
429                         nGood++;
430                 }
431                 else {
432                         esdTrack->SetEMCALcluster(-cluster->Index());
433                         nFake++;
434                 }
435                 nSaved++;
436         }
437         /*
438         AliEMCALTrack *track = 0;
439         TObjArrayIter tracks(fTracks);
440         while ( (track = (AliEMCALTrack*)tracks.Next()) ) {
441                 trackID = track->GetSeedIndex();
442                 clusterID = track->GetMatchedClusterIndex();
443                 AliESDtrack *esdTrack = esd->GetTrack(trackID);
444                 if (!esdTrack) continue;
445                 if (clusterID < 0) {
446                         esdTrack->SetEMCALcluster(kUnmatched);
447                 }
448                 else {
449                         AliEMCALMatchCluster *cluster = (AliEMCALMatchCluster*)fClusters->At(clusterID);
450                         if (!cluster) continue;
451                         if (esdTrack->GetLabel() == cluster->Label()) {
452                                 nGood++;
453                                 esdTrack->SetEMCALcluster(cluster->Index());
454                         }
455                         else {
456                                 esdTrack->SetEMCALcluster(-cluster->Index());
457                         }
458                         nSaved++;
459                 }
460         }
461         */
462         AliInfo(Form("Saved %d matches [%d good + %d fake]", nSaved, nGood, nFake));
463
464         return 0;
465 }
466 //
467 //------------------------------------------------------------------------------
468 //
469 void AliEMCALTracker::SetTrackCorrectionMode(Option_t *option)
470 {
471         //
472         // Set track correction mode
473         // gest the choice in string format and converts into 
474         // internal enum
475         //
476         
477         TString opt(option);
478         opt.ToUpper();
479         
480         if (!opt.CompareTo("NONE")) {
481                 fTrackCorrMode = kTrackCorrNone;
482         }
483         else if (!opt.CompareTo("MMB")) {
484                 fTrackCorrMode = kTrackCorrMMB;
485         }
486         else if (!opt.CompareTo("FIXED")) {
487                 fTrackCorrMode = kTrackCorrFixed;
488         }
489         else {
490                 cerr << "E-AliEMCALTracker::SetTrackCorrectionMode '" << option << "': Unrecognized option" << endl;
491         }
492 }
493
494 //
495 //------------------------------------------------------------------------------
496 //
497 Double_t AliEMCALTracker::AngleDiff(Double_t angle1, Double_t angle2)
498 {
499         // 
500         // [PRIVATE]
501         // Given two angles in radiants, it converts them in the range 0-2pi
502         // then computes their true difference, i.e. if the difference a1-a2
503         // results to be larger than 180 degrees, it returns 360 - diff.
504         //
505         
506         if (angle1 < 0.0) angle1 += TMath::TwoPi();
507         if (angle1 > TMath::TwoPi()) angle1 -= TMath::TwoPi();
508         if (angle2 < 0.0) angle2 += TMath::TwoPi();
509         if (angle2 > TMath::TwoPi()) angle2 -= TMath::TwoPi();
510         
511         Double_t diff = TMath::Abs(angle1 - angle2);
512         if (diff > TMath::Pi()) diff = TMath::TwoPi() - diff;
513         
514         if (angle2 > angle1) diff = -diff;
515         
516         return diff;
517 }
518 //
519 //------------------------------------------------------------------------------
520 //
521 Double_t AliEMCALTracker::CheckPair
522 (AliEMCALTrack *track, AliEMCALMatchCluster *cl)
523 {
524         //
525         // Given a track and a cluster,
526         // propagates the first to the radius of the second.
527         // Then, checks the propagation point against all cuts.
528         // If at least a cut is not passed, a valuer equal to 
529         // twice the maximum allowed distance is passed (so the value returned
530         // will not be taken into account when creating matches)
531         //
532         
533         // TEMP
534         Bool_t isTrue = kFALSE;
535 //      if (tr->GetSeedLabel() == cl->Label()) {
536 //              isTrue = kTRUE;
537 //      }
538         
539         // copy track into temporary variable
540         AliEMCALTrack *tr = new AliEMCALTrack(*track);
541         
542         Double_t distance = 2.0 * fMaxDist;
543         
544         // check against cut on difference 'alpha - phi'
545         Double_t phi = TMath::ATan2(cl->Y(), cl->X());
546         phi = AngleDiff(phi, tr->GetAlpha());
547         if (phi < fCutAlphaMin || phi > fCutAlphaMax) return distance;
548         
549         // try to propagate to cluster radius
550         // (return the 'distance' value if it fails)
551         Double_t pos[3], &x = pos[0], &y = pos[1], &z = pos[2];
552         Double_t x0, rho;
553         tr->GetXYZ(pos);
554         Double_t rt = TMath::Sqrt(x*x + y*y);
555         Double_t rc = TMath::Sqrt(cl->X()*cl->X() + cl->Y()*cl->Y());
556         
557         if (fTrackCorrMode == kTrackCorrMMB) {
558                 Double_t pos1[3], pos2[3], param[6];
559                 pos1[0] = x;
560                 pos1[1] = y;
561                 pos1[2] = z;
562                 pos2[0] = cl->X();
563                 pos2[1] = cl->Y();
564                 pos2[2] = cl->Z();
565                 AliKalmanTrack::MeanMaterialBudget(pos1, pos2, param);
566                 rho = param[0];
567                 x0 = param[1];
568         }
569         else if (fTrackCorrMode == kTrackCorrFixed) {
570                 rho = fRho;
571                 x0 = fX0;
572         }
573         else {
574                 rho = 0.0;
575                 x0 = 0.0;
576         }
577         if (fNPropSteps) {
578                 Int_t i;
579                 Double_t r;
580                 cout.setf(ios::fixed);
581                 cout.precision(5);
582                 if (isTrue) cout << "Init : " << rt << ' ' << x << ' ' << y << ' ' << z << endl;
583                 for (i = 0; i < fNPropSteps; i++) {
584                         r = rt + (rc - rt) * ((Double_t)(i+1)/(Double_t)fNPropSteps);
585                         if (!tr->PropagateTo(r, x0, rho)) return distance;
586                         tr->GetXYZ(pos);
587                         if (isTrue) cout << "Step : " << r << ' ' << x << ' ' << y << ' ' << z << endl;
588                 }
589                 if (isTrue) cout << "Clstr: " << rc << ' ' << cl->X() << ' ' << cl->Y() << ' ' << cl->Z() << endl;
590         }
591         else {
592                 // when no steps are used, no correction makes sense
593                 //if (!tr->PropagateTo(rc, 0.0, 0.0)) return distance;
594                 if (!tr->PropagateToGlobal(cl->X(), cl->Y(), cl->Z(), 0.0, 0.0)) return distance;
595                 /*
596                 Bool_t propOK = kFALSE;
597                 cout << "START" << endl;
598                 Double_t dist, rCHK, bestDist = 10000000.0;
599                 for (Double_t rTMP = rc; rTMP> rc*0.95; rTMP -= 0.1) {
600                         if (!tr->PropagateTo(rTMP)) continue;
601                         propOK = kTRUE;
602                         tr->GetXYZ(pos);
603                         rCHK = TMath::Sqrt(x*x + y*y);
604                         dist = TMath::Abs(rCHK - rc);
605                         cout << rCHK << " vs. " << rc << endl;
606                         
607                         if (TMath::Abs(rCHK - rc) < 0.01) break;
608                 }
609                 cout << "STOP" << endl;
610                 if (!propOK) return distance;
611                 */
612         }
613         
614         // get global propagation of track at end of propagation
615         tr->GetXYZ(pos);
616         
617         // check angle cut
618         TVector3 vc(cl->X(), cl->Y(), cl->Z());
619         TVector3 vt(x, y, z);
620         Double_t angle = TMath::Abs(vc.Angle(vt)) * TMath::RadToDeg();
621         // check: where is the track?
622         Double_t r, phiT, phiC;
623         r = TMath::Sqrt(pos[0]*pos[0] + pos[1]*pos[1]);
624         phiT = TMath::ATan2(pos[1], pos[0]) * TMath::RadToDeg();
625         phiC = vc.Phi() * TMath::RadToDeg();
626         //cout << "Propagated R, phiT, phiC = " << r << ' ' << phiT << ' ' << phiC << endl;
627         
628         if (angle > fCutAngle) {
629                 //cout << "angle" << endl;
630                 return distance;
631         }
632                 
633         // compute differences wr to each coordinate
634         x -= cl->X();
635         if (TMath::Abs(x) > fCutX) {
636                 //cout << "cut X" << endl;
637                 return distance;
638         }
639         y -= cl->Y();
640         if (TMath::Abs(y) > fCutY) {
641                 //cout << "cut Y" << endl;
642                 return distance;
643         }
644         z -= cl->Z();
645         if (TMath::Abs(z) > fCutZ) {
646                 //cout << "cut Z" << endl;
647                 return distance;
648         }
649         
650         // compute true distance
651         distance = TMath::Sqrt(x*x + y*y + z*z);
652         //Double_t temp = CheckPairV2(tr, cl);
653         //if (temp < distance) return temp; else 
654         
655         // delete temporary object
656         delete tr;
657         
658         return distance;
659 }
660 //
661 //------------------------------------------------------------------------------
662 //
663 Double_t AliEMCALTracker::CheckPairV2
664 (AliEMCALTrack *tr, AliEMCALMatchCluster *cl)
665 {
666         //
667         // Given a track and a cluster,
668         // propagates the first to the radius of the second.
669         // Then, checks the propagation point against all cuts.
670         // If at least a cut is not passed, a valuer equal to 
671         // twice the maximum allowed distance is passed (so the value returned
672         // will not be taken into account when creating matches)
673         //
674         
675         // TEMP
676 //      Bool_t isTrue = kFALSE;
677 //      if (tr->GetSeedLabel() == cl->Label()) {
678 //              isTrue = kTRUE;
679 //              cout << "TRUE MATCH!!!" << endl;
680 //      }
681         
682         Double_t distance = 2.0 * fMaxDist;
683         
684         Double_t x0, rho;
685         if (fTrackCorrMode == kTrackCorrMMB) {
686                 Double_t pos1[3], pos2[3], param[6];
687                 tr->GetXYZ(pos1);
688 //              pos1[0] = x;
689 //              pos1[1] = y;
690 //              pos1[2] = z;
691                 pos2[0] = cl->X();
692                 pos2[1] = cl->Y();
693                 pos2[2] = cl->Z();
694                 AliKalmanTrack::MeanMaterialBudget(pos1, pos2, param);
695                 rho = param[0];
696                 x0 = param[1];
697         }
698         else if (fTrackCorrMode == kTrackCorrFixed) {
699                 rho = fRho;
700                 x0 = fX0;
701         }
702         else {
703                 rho = 0.0;
704                 x0 = 0.0;
705         }
706         
707         // check against cut on difference 'alpha - phi'
708         Double_t phi = TMath::ATan2(cl->Y(), cl->X());
709         phi = AngleDiff(phi, tr->GetAlpha());
710         if (phi < fCutAlphaMin || phi > fCutAlphaMax) return distance;
711         
712         // get cluster position and put them into a vector
713         TVector3 vc(cl->X(), cl->Y(), cl->Z());
714         // rotate the vector in order to put all clusters on a plane intersecting 
715         // vertically the X axis; the angle depends on the sector
716         Double_t clusterRot, clusterPhi = vc.Phi() * TMath::RadToDeg();
717         if (clusterPhi < 0.0) clusterPhi += 360.0;
718         if (clusterPhi < 100.0) {
719                 clusterRot = -90.0;
720         }
721         else if (clusterPhi < 120.0) {
722                 clusterRot = -110.0;
723         }
724         else if (clusterPhi < 140.0) {
725                 clusterRot = -130.0;
726         }
727         else if (clusterPhi < 160.0) {
728                 clusterRot = -150.0;
729         }
730         else if (clusterPhi < 180.0) {
731                 clusterRot = -170.0;
732         }
733         else {
734                 clusterRot = -190.0;
735         }
736         vc.RotateZ(clusterRot * TMath::DegToRad());
737         // generate a track from the ESD track selected
738         AliEMCALTrack *track = new AliEMCALTrack(*tr);
739         // compute the 'phi' coordinate of the intersection point to 
740         // the EMCAL surface
741         Double_t x = vc.X();
742         Double_t y;
743         track->GetYAt(vc.X(), track->GetBz(), y);
744         Double_t tmp = x*TMath::Cos(track->GetAlpha()) - y*TMath::Sin(track->GetAlpha());
745         y = x*TMath::Sin(track->GetAlpha()) + y*TMath::Cos(track->GetAlpha());
746         x = tmp;
747         Double_t trackPhi = TMath::ATan2(y, x) * TMath::RadToDeg();
748         // compute phi difference
749         Double_t dphi = trackPhi - clusterPhi;
750         if (TMath::Abs(dphi) > 180.0) {
751                 dphi = 360.0 - TMath::Abs(dphi);
752                 if (clusterPhi > trackPhi) dphi = -dphi;
753         }
754         // propagate track to the X position of rotated cluster
755         // and get the vector of X, Y, Z in the local ref. frame of the track
756         track->PropagateTo(vc.X(), x0, rho);
757         TVector3 vt(track->GetX(), track->GetY(), track->GetZ());
758         vt.RotateZ((clusterPhi - trackPhi) * TMath::DegToRad());
759         TVector3 vdiff = vt-vc;
760                 
761         // compute differences wr to each coordinate
762         if (vdiff.X() > fCutX) return distance;
763         if (vdiff.Y() > fCutY) return distance;
764         if (vdiff.Z() > fCutZ) return distance;
765         
766         // compute true distance
767         distance = vdiff.Mag();
768         return distance;
769 }
770 //
771 //------------------------------------------------------------------------------
772 //
773 Double_t AliEMCALTracker::CheckPairV3
774 (AliEMCALTrack *track, AliEMCALMatchCluster *cl)
775 {
776         //
777         // Given a track and a cluster,
778         // propagates the first to the radius of the second.
779         // Then, checks the propagation point against all cuts.
780         // If at least a cut is not passed, a valuer equal to 
781         // twice the maximum allowed distance is passed (so the value returned
782         // will not be taken into account when creating matches)
783         //
784         
785         AliEMCALTrack tr(*track);
786         
787         Int_t    sector;
788         Double_t distance = 2.0 * fMaxDist;
789         Double_t dx, dy, dz;
790         Double_t phi, alpha, slope, tgtXnum, tgtXden, sectorWidth = 20.0 * TMath::DegToRad();
791         Double_t xcurr, xprop, param[6] = {0., 0., 0., 0., 0., 0.}, d, x0, rho, bz;
792         Double_t x[3], x1[3], x2[3];
793         
794         // get initial track position
795         xcurr = tr.GetX();
796         
797         // evaluate the EMCAL sector number
798         phi = cl->Phi();
799         if (phi < 0.0) phi += TMath::TwoPi();
800         sector = (Int_t)(phi / sectorWidth);
801         alpha = ((Double_t)sector + 0.5) * sectorWidth;
802         // evaluate the corresponding X for track propagation
803         slope = TMath::Tan(alpha - 0.5*TMath::Pi());
804         tgtXnum = cl->Y() - slope * cl->X();
805         tgtXden = TMath::Sqrt(1.0 + slope*slope);
806         xprop = TMath::Abs(tgtXnum / tgtXden);
807         
808         // propagate by small steps
809         tr.GetXYZ(x1);
810         bz = tr.GetBz();
811         if (!tr.GetXYZAt(xprop, bz, x2)) return distance;
812         //AliKalmanTrack::MeanMaterialBudget(x1, x2, param);
813         d = param[4];
814         rho = param[0];
815         x0 = param[1];
816         if (!tr.PropagateTo(xprop, d*rho/x0, x0)) return distance;
817         //if (!tr.PropagateTo(xprop, 0.0, 0.0)) return distance;
818         
819         // get propagated position at the end
820         tr.GetXYZ(x);
821         dx = TMath::Abs(x[0] - cl->X());
822         dy = TMath::Abs(x[1] - cl->Y());
823         dz = TMath::Abs(x[2] - cl->Z());
824         if (dx > fCutX || dy > fCutY || dz > fCutZ) return distance;
825         
826         distance = TMath::Sqrt(dx*dx + dy*dy + dz*dz);
827         
828         return distance;
829 }
830 //
831 //------------------------------------------------------------------------------
832 //
833 Bool_t AliEMCALTracker::PropagateToEMCAL(AliEMCALTrack *tr)
834 {
835         //
836         // Propagates the track to the proximity of the EMCAL surface
837         //
838         
839         Double_t xcurr, xtemp, xprop = 438.0, step = 10.0, param[6], d, x0, rho, bz;
840         Double_t x1[3], x2[3];
841         
842         // get initial track position
843         xcurr = tr->GetX();
844         
845         // propagate by small steps
846         for (xtemp = xcurr + step; xtemp < xprop; xtemp += step) {
847                 // to compute material budget, take current position and 
848                 // propagated hypothesis without energy loss
849                 tr->GetXYZ(x1);
850                 bz = tr->GetBz();
851                 if (!tr->GetXYZAt(xtemp, bz, x2)) return kFALSE;
852                 AliKalmanTrack::MeanMaterialBudget(x1, x2, param);
853                 d = param[4];
854                 rho = param[0];
855                 x0 = param[1];
856                 if (!tr->PropagateTo(xtemp, d*rho/x0, x0)) return kFALSE;
857         }
858         
859         return kTRUE;
860 }
861 //
862 //------------------------------------------------------------------------------
863 //
864 Int_t AliEMCALTracker::CreateMatches()
865 {
866         //
867         // Creation of matches between tracks and clusters.
868         // For each ESD track collected by ReadESD(), an AliEMCALTrack is made.
869         // If it finds a cluster close enough to its propagation to EMCAL,
870         // which passes all cuts, its index is stored.
871         // If many clusters are found which satisfy the criteria described above, 
872         // only the closest one is stored.
873         // At this level, it is possible that two tracks share the same cluster.
874         //
875         
876         // if matches collection is already present, it is deleted
877         if (fMatches) {
878                 fMatches->Delete();
879                 delete fMatches;
880         }
881         fMatches = new TList;
882         
883         // initialize counters and indexes
884         Int_t count = 0;
885         Int_t ic, nClusters = (Int_t)fClusters->GetEntries();
886         Int_t it, nTracks = fTracks->GetEntries();
887         
888         // external loop on clusters, internal loop on tracks
889         Double_t dist;
890         for (ic = 0; ic < nClusters; ic++) {
891                 AliEMCALMatchCluster *cluster = (AliEMCALMatchCluster*)fClusters->At(ic);
892                 for (it = 0; it < nTracks; it++) {
893                         AliEMCALTrack *track = (AliEMCALTrack*)fTracks->At(it);
894                         dist = CheckPair(track, cluster);
895                         //cout << dist << endl;
896                         if (dist <= fMaxDist) {
897                                 AliEMCALMatch *candidate = new AliEMCALMatch;
898                                 candidate->SetIndexT(it);
899                                 candidate->SetIndexC(ic);
900                                 candidate->SetDistance(dist);
901                                 candidate->SetPt(track->GetPt());
902                                 fMatches->Add(candidate);
903                                 count++;
904                         }
905                 }
906         }
907                 
908         /*
909         // loop on clusters and tracks
910         Int_t icBest;
911         Double_t dist, distBest;
912         for (it = 0; it < nTracks; it++) {
913                 AliEMCALTrack *track = (AliEMCALTrack*)fTracks->At(it);
914                 if (!track) continue;
915                 icBest = -1;
916                 distBest = fMaxDist;
917                 for (ic = 0; ic < nClusters; ic++) {
918                         AliEMCALMatchCluster *cluster = (AliEMCALMatchCluster*)fClusters->At(ic);
919                         if (!cluster) continue;
920                         dist = CheckPair(track, cluster);
921                         if (dist < distBest) {
922                                 distBest = dist;
923                                 icBest = ic;
924                         }
925                 }
926                 if (icBest >= 0) {
927                         track->SetMatchedClusterIndex(icBest);
928                         track->SetMatchedClusterDist(distBest);
929                         count++;
930                 }
931                 else {
932                         track->SetMatchedClusterIndex(-1);
933                 }
934         }
935         */
936         
937         return count;
938 }
939 //
940 //------------------------------------------------------------------------------
941 //
942 Int_t AliEMCALTracker::SolveCompetitions()
943 {
944         //
945         // Match selector.
946         // The match list is sorted from the best to the worst match, w.r. to the 
947         // distance between track prolongation and cluster position.
948         // Based on this criterion, starting from the first (best) match, a flag
949         // is set to both the involved track and cluster, and all matches containing
950         // an already used track or cluster are removed, leaving only the best match
951         // for each cluster.
952         //
953         
954         // sort matches with respect to track-cluster distance
955         fMatches->Sort(kSortAscending);
956         
957         // keep track of eliminated matches
958         Int_t count = 0;
959         
960         // initialize flags to check repetitions
961         Int_t ic, nClusters = (Int_t)fClusters->GetEntries();
962         Int_t it, nTracks = fTracks->GetEntries();
963         Bool_t *usedC = new Bool_t[nClusters];
964         Bool_t *usedT = new Bool_t[nTracks];
965         for (ic = 0; ic < nClusters; ic++) usedC[ic] = kFALSE;
966         for (it = 0; it < nTracks; it++) usedT[it] = kFALSE;
967         
968         // loop on matches
969         TListIter iter(fMatches);
970         AliEMCALMatch *match = 0;
971         while ( (match = (AliEMCALMatch*)iter.Next()) ) {
972                 ic = match->GetIndexC();
973                 it = match->GetIndexT();
974                 if (!usedT[it] && !usedC[ic]) {
975                         usedT[it] = kTRUE;
976                         usedC[ic] = kTRUE;
977                         match->CanBeSaved() = kTRUE;
978                 }
979                 else {
980                         count++;
981                 }
982         }
983         
984         /*
985         Int_t it1, it2, nTracks = (Int_t)fTracks->GetEntries();
986         AliEMCALTrack *track1 = 0, *track2 = 0;
987         for (it1 = 0; it1 < nTracks; it1++) {
988                 track1 = (AliEMCALTrack*)fTracks->At(it1);
989                 if (!track1) continue;
990                 if (track1->GetMatchedClusterIndex() < 0) continue;
991                 for (it2 = it1+1; it2 < nTracks; it2++) {
992                         track2 = (AliEMCALTrack*)fTracks->At(it2);
993                         if (!track2) continue;
994                         if (track2->GetMatchedClusterIndex() < 0) continue;
995                         if (track1->GetMatchedClusterIndex() != track2->GetMatchedClusterIndex()) continue;
996                         count++;
997                         if (track1->GetMatchedClusterDist() < track2->GetMatchedClusterDist()) {
998                                 track2->SetMatchedClusterIndex(-1);
999                         }
1000                         else if (track2->GetMatchedClusterDist() < track1->GetMatchedClusterDist()) {
1001                                 track1->SetMatchedClusterIndex(-1);
1002                         }
1003                 }
1004         }
1005         */
1006         
1007         return count;
1008 }
1009 //
1010 //------------------------------------------------------------------------------
1011 //
1012 void AliEMCALTracker::UnloadClusters() 
1013 {
1014         //
1015         // Free memory from clusters
1016         //
1017         
1018         Clear("CLUSTERS");
1019 }
1020 //
1021 //------------------------------------------------------------------------------
1022 //
1023 AliEMCALTracker::AliEMCALMatchCluster::AliEMCALMatchCluster(Int_t index, AliEMCALRecPoint *recPoint)
1024   : fIndex(index),
1025     fLabel(recPoint->GetPrimaryIndex()),
1026     fX(0.),
1027     fY(0.),
1028     fZ(0.)
1029 {
1030         //
1031         // Translates an AliEMCALRecPoint object into the internal format.
1032         // Index of passed cluster in its native array must be specified.
1033         //
1034         TVector3 clpos;
1035         recPoint->GetGlobalPosition(clpos);
1036         
1037         fX = clpos.X();
1038         fY = clpos.Y();
1039         fZ = clpos.Z();
1040 }
1041 //
1042 //------------------------------------------------------------------------------
1043 //
1044 AliEMCALTracker::AliEMCALMatchCluster::AliEMCALMatchCluster(Int_t index, AliESDCaloCluster *caloCluster)
1045   : fIndex(index),
1046     fLabel(caloCluster->GetLabel()),
1047     fX(0.),
1048     fY(0.),
1049     fZ(0.)
1050 {
1051         //
1052         // Translates an AliESDCaloCluster object into the internal format.
1053         // Index of passed cluster in its native array must be specified.
1054         //
1055         Float_t clpos[3];
1056         caloCluster->GetPosition(clpos);
1057         
1058         fX = (Double_t)clpos[0];
1059         fY = (Double_t)clpos[1];
1060         fZ = (Double_t)clpos[2];
1061 }
1062 //
1063 //------------------------------------------------------------------------------
1064 //
1065 Int_t AliEMCALTracker::AliEMCALMatch::Compare(const TObject *obj) const 
1066 {
1067         //
1068         // Tracks compared wrt their distance from matched point
1069         //
1070         
1071         AliEMCALTracker::AliEMCALMatch *that = (AliEMCALTracker::AliEMCALMatch*)obj;
1072         
1073         Double_t thisDist = fPt;//fDistance;
1074         Double_t thatDist = that->fPt;//that->GetDistance();
1075         
1076         if (thisDist > thatDist) return 1;
1077         else if (thisDist < thatDist) return -1;
1078         return 0;
1079 }
1080
1081 AliEMCALTracker::AliEMCALMatch::AliEMCALMatch() 
1082   : TObject(),
1083     fCanBeSaved(kFALSE), 
1084     fIndexC(0), 
1085     fIndexT(0), 
1086     fDistance(0.),
1087     fPt(0.)
1088 {
1089   //default constructor
1090
1091 }
1092
1093 AliEMCALTracker::AliEMCALMatch::AliEMCALMatch(const AliEMCALMatch& copy)
1094   : TObject(),
1095     fCanBeSaved(copy.fCanBeSaved),
1096     fIndexC(copy.fIndexC),
1097     fIndexT(copy.fIndexT),
1098     fDistance(copy.fDistance),
1099     fPt(copy.fPt)
1100 {
1101   //copy ctor
1102 }