]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/EMCAL/AliEmcalClusTrackMatcherTask.cxx
forgotton
[u/mrichter/AliRoot.git] / PWG / EMCAL / AliEmcalClusTrackMatcherTask.cxx
1 // $Id$
2 //
3 // Track/cluster matcher
4 // 
5 // Author: C.Loizides, S.Aiola
6
7 #include "AliEmcalClusTrackMatcherTask.h"
8
9 #include <TClonesArray.h>
10
11 #include "AliEmcalParticle.h"
12 #include "AliLog.h"
13 #include "AliPicoTrack.h"
14 #include "AliVCluster.h"
15 #include "AliVTrack.h"
16
17 ClassImp(AliEmcalClusTrackMatcherTask)
18
19 //________________________________________________________________________
20 AliEmcalClusTrackMatcherTask::AliEmcalClusTrackMatcherTask() : 
21   AliAnalysisTaskEmcal("AliEmcalClusTrackMatcherTask",kFALSE),
22   fMaxDistance(0.06)
23 {
24   // Constructor.
25 }
26
27 //________________________________________________________________________
28 AliEmcalClusTrackMatcherTask::AliEmcalClusTrackMatcherTask(const char *name) : 
29   AliAnalysisTaskEmcal(name,kFALSE),
30   fMaxDistance(0.06)
31 {
32   // Standard constructor.
33 }
34
35 //________________________________________________________________________
36 AliEmcalClusTrackMatcherTask::~AliEmcalClusTrackMatcherTask()
37 {
38   // Destructor.
39 }
40
41 //________________________________________________________________________
42 Bool_t AliEmcalClusTrackMatcherTask::Run() 
43 {
44   // Run the matching for the selected options.
45   
46   const Double_t maxd2 = fMaxDistance*fMaxDistance;
47
48   const Int_t nC = fCaloClusters->GetEntries();
49   const Int_t nT = fTracks->GetEntries();
50
51   // set the links between tracks and clusters
52   for (Int_t c = 0; c < nC; ++c) {
53     AliEmcalParticle *partC = static_cast<AliEmcalParticle*>(fCaloClusters->At(c));
54     if (!partC)
55       continue;
56     if (!AcceptEmcalPart(partC))
57       continue;
58     for (Int_t t = 0; t < nT; ++t) {
59       AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(t));
60       if (!partT)
61         continue;
62       if (!AcceptEmcalPart(partT))
63         continue;
64       AliVCluster *clust = partC->GetCluster();
65       AliVTrack   *track = partT->GetTrack()  ;
66       Double_t deta = 999;
67       Double_t dphi = 999;
68       AliPicoTrack::GetEtaPhiDiff(track, clust, dphi, deta);
69       Double_t d2 = deta * deta + dphi * dphi;
70       if (d2 > maxd2)
71         continue;
72       Double_t d = TMath::Sqrt(d2);
73       partC->AddMatchedObj(t, d);
74       partT->AddMatchedObj(c, d);
75     }
76   }
77
78   for (Int_t c = 0; c < nC; ++c) {
79     AliEmcalParticle *partC = static_cast<AliEmcalParticle*>(fCaloClusters->At(c));
80     if (!partC)
81       continue;
82     if (partC->GetNumberOfMatchedObj() <= 0)
83       continue;
84     const UInt_t matchedId = partC->GetMatchedObjId();
85     AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(matchedId));
86     AliVCluster *clust = partC->GetCluster();
87     AliVTrack   *track = partT->GetTrack()  ;
88     Double_t deta = 999;
89     Double_t dphi = 999;
90     AliPicoTrack::GetEtaPhiDiff(track, clust, dphi, deta);
91     clust->SetEmcCpvDistance(matchedId);
92     clust->SetTrackDistance(deta, dphi);
93   }
94
95   for (Int_t t = 0; t < nT; ++t) {
96     AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(t));
97     if (!partT)
98       continue;
99     if (partT->GetNumberOfMatchedObj() <= 0)
100       continue;
101     AliVTrack *track = partT->GetTrack();
102     track->SetEMCALcluster(partT->GetMatchedObjId());
103   }
104
105   return kTRUE;
106 }
107