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