]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG/EMCAL/AliEmcalClusTrackMatcherTask.cxx
fix for mem leak
[u/mrichter/AliRoot.git] / PWG / EMCAL / AliEmcalClusTrackMatcherTask.cxx
CommitLineData
aee86258 1// $Id$
980821ba 2//
cd231d42 3// Track/cluster matcher
4//
d258e37e 5// Author: C.Loizides, S.Aiola
aee86258 6
c79b2854 7#include "AliEmcalClusTrackMatcherTask.h"
8
64b4d2f2 9#include <TClonesArray.h>
35789a2d 10
d258e37e 11#include "AliEmcalParticle.h"
35789a2d 12#include "AliLog.h"
c79b2854 13#include "AliPicoTrack.h"
14#include "AliVCluster.h"
15#include "AliVTrack.h"
aee86258 16
17ClassImp(AliEmcalClusTrackMatcherTask)
18
d258e37e 19//________________________________________________________________________
20AliEmcalClusTrackMatcherTask::AliEmcalClusTrackMatcherTask() :
a18c300f 21 AliAnalysisTaskEmcal("AliEmcalClusTrackMatcherTask",kFALSE),
374d55cc 22 fMaxDistance(0.06)
d258e37e 23{
24 // Constructor.
25}
26
aee86258 27//________________________________________________________________________
28AliEmcalClusTrackMatcherTask::AliEmcalClusTrackMatcherTask(const char *name) :
a18c300f 29 AliAnalysisTaskEmcal(name,kFALSE),
374d55cc 30 fMaxDistance(0.06)
aee86258 31{
32 // Standard constructor.
aee86258 33}
34
35//________________________________________________________________________
36AliEmcalClusTrackMatcherTask::~AliEmcalClusTrackMatcherTask()
37{
8bec746d 38 // Destructor.
aee86258 39}
40
41//________________________________________________________________________
d258e37e 42Bool_t AliEmcalClusTrackMatcherTask::Run()
aee86258 43{
d258e37e 44 // Run the matching for the selected options.
45
e3d22acf 46 const Double_t maxd2 = fMaxDistance*fMaxDistance;
47
c79b2854 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)
d258e37e 55 continue;
c79b2854 56 if (!AcceptEmcalPart(partC))
e3d22acf 57 continue;
c79b2854 58 for (Int_t t = 0; t < nT; ++t) {
59 AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(t));
60 if (!partT)
d258e37e 61 continue;
c79b2854 62 if (!AcceptEmcalPart(partT))
e3d22acf 63 continue;
c79b2854 64 AliVCluster *clust = partC->GetCluster();
65 AliVTrack *track = partT->GetTrack() ;
d258e37e 66 Double_t deta = 999;
67 Double_t dphi = 999;
c79b2854 68 AliPicoTrack::GetEtaPhiDiff(track, clust, dphi, deta);
e3d22acf 69 Double_t d2 = deta * deta + dphi * dphi;
c79b2854 70 if (d2 > maxd2)
71 continue;
72 Double_t d = TMath::Sqrt(d2);
73 partC->AddMatchedObj(t, d);
74 partT->AddMatchedObj(c, d);
aee86258 75 }
aee86258 76 }
c79b2854 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 }
a18c300f 104
105 return kTRUE;
aee86258 106}
c79b2854 107