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