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