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