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