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