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