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