AliPhysics  0bb4a45 (0bb4a45)
AliCaloTrackMatcher.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: Daniel Mühlheim *
5  * Version 1.0 *
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 
17 //---------------------------------------------
18 // Basic Track Matching Class
19 //---------------------------------------------
21 
22 
23 #include "AliAnalysisManager.h"
24 #include "AliAODEvent.h"
25 #include "AliCaloTrackMatcher.h"
26 #include "AliEMCALRecoUtils.h"
27 #include "AliESDEvent.h"
28 #include "AliESDtrack.h"
29 #include "AliESDtrackCuts.h"
30 #include "AliInputEventHandler.h"
31 #include "AliTrackerBase.h"
32 #include "AliV0ReaderV1.h"
33 
34 #include "TAxis.h"
35 #include "TChain.h"
36 #include "TH1F.h"
37 #include "TF1.h"
38 
39 #include <vector>
40 #include <map>
41 #include <utility>
42 
43 class iostream;
44 
45 using namespace std;
46 
47 
48 ClassImp(AliCaloTrackMatcher)
49 
50 //________________________________________________________________________
51 AliCaloTrackMatcher::AliCaloTrackMatcher(const char *name, Int_t clusterType) : AliAnalysisTaskSE(name),
52  fClusterType(clusterType),
53  fV0ReaderName(""),
54  fCorrTaskSetting(""),
55  fAnalysisTrainMode("Grid"),
56  fMatchingWindow(200),
57  fMatchingResidual(0.2),
58  fRunNumber(-1),
59  fGeomEMCAL(NULL),
60  fGeomPHOS(NULL),
61  fMapTrackToCluster(),
62  fMapClusterToTrack(),
63  fNEntries(1),
64  fVectorDeltaEtaDeltaPhi(0),
65  fMap_TrID_ClID_ToIndex(),
66  fSecMapTrackToCluster(),
67  fSecMapClusterToTrack(),
68  fSecNEntries(1),
69  fSecVectorDeltaEtaDeltaPhi(0),
70  fSecMap_TrID_ClID_ToIndex(),
71  fSecMap_TrID_ClID_AlreadyTried(),
72  fListHistos(NULL),
73  fHistControlMatches(NULL),
74  fSecHistControlMatches(NULL)
75 {
76  // Default constructor
77  DefineInput(0, TChain::Class());
78 }
79 
80 //________________________________________________________________________
82  // default deconstructor
83  fMapTrackToCluster.clear();
84  fMapClusterToTrack.clear();
85  fVectorDeltaEtaDeltaPhi.clear();
86  fMap_TrID_ClID_ToIndex.clear();
87 
88  fSecMapTrackToCluster.clear();
89  fSecMapClusterToTrack.clear();
90  fSecVectorDeltaEtaDeltaPhi.clear();
91  fSecMap_TrID_ClID_ToIndex.clear();
92  fSecMap_TrID_ClID_AlreadyTried.clear();
93 
94  if(fHistControlMatches) delete fHistControlMatches;
95  if(fSecHistControlMatches) delete fSecHistControlMatches;
96  if(fAnalysisTrainMode.EqualTo("Grid")){
97  if(fListHistos != NULL){
98  delete fListHistos;
99  }
100  }
101 }
102 
103 //________________________________________________________________________
105  fMapTrackToCluster.clear();
106  fMapClusterToTrack.clear();
107  fVectorDeltaEtaDeltaPhi.clear();
108  fMap_TrID_ClID_ToIndex.clear();
109 
110  fSecMapTrackToCluster.clear();
111  fSecMapClusterToTrack.clear();
112  fSecVectorDeltaEtaDeltaPhi.clear();
113  fSecMap_TrID_ClID_ToIndex.clear();
114  fSecMap_TrID_ClID_AlreadyTried.clear();
115 }
116 
117 //________________________________________________________________________
119  if(fListHistos != NULL){
120  delete fListHistos;
121  fListHistos = NULL;
122  }
123  if(fListHistos == NULL){
124  fListHistos = new TList();
125  fListHistos->SetOwner(kTRUE);
126  fListHistos->SetName(Form("CaloTrackMatcher_%i",fClusterType));
127  }
128 
129  // Create User Output Objects
130  fHistControlMatches = new TH2F(Form("ControlMatches_%i",fClusterType),Form("ControlMatches_%i",fClusterType),7,-0.5,6.5,50,0.15,200);
131  SetLogBinningYTH2(fHistControlMatches);
132  fHistControlMatches->GetYaxis()->SetTitle("track pT (GeV/c)");
133  fHistControlMatches->GetXaxis()->SetBinLabel(1,"nTr in");
134  fHistControlMatches->GetXaxis()->SetBinLabel(2,"no inner Tr params || track not in right direction");
135  fHistControlMatches->GetXaxis()->SetBinLabel(3,"failed 1st pro-step");
136  fHistControlMatches->GetXaxis()->SetBinLabel(4,"Tr not in EMCal acc");
137  fHistControlMatches->GetXaxis()->SetBinLabel(5,"failed 2nd pro-step");
138  fHistControlMatches->GetXaxis()->SetBinLabel(6,"w/o match to cluster");
139  fHistControlMatches->GetXaxis()->SetBinLabel(7,"nTr out, w/ match");
140  fListHistos->Add(fHistControlMatches);
141 
142  fSecHistControlMatches = new TH2F(Form("ControlSecMatches_%i",fClusterType),Form("ControlSecMatches_%i",fClusterType),7,-0.5,6.5,50,0.15,200);
143  SetLogBinningYTH2(fSecHistControlMatches);
144  fSecHistControlMatches->GetYaxis()->SetTitle("track pT (GeV/c)");
145  fSecHistControlMatches->GetXaxis()->SetBinLabel(1,"nTr in");
146  fSecHistControlMatches->GetXaxis()->SetBinLabel(2,"no inner Tr params");
147  fSecHistControlMatches->GetXaxis()->SetBinLabel(3,"failed 1st pro-step");
148  fSecHistControlMatches->GetXaxis()->SetBinLabel(4,"Tr not in EMCal acc");
149  fSecHistControlMatches->GetXaxis()->SetBinLabel(5,"failed 2nd pro-step");
150  fSecHistControlMatches->GetXaxis()->SetBinLabel(6,"w/o match to cluster");
151  fSecHistControlMatches->GetXaxis()->SetBinLabel(7,"nTr out, w/ match");
152  fListHistos->Add(fSecHistControlMatches);
153 }
154 
155 //________________________________________________________________________
157  // Initialize function to be called once before analysis
158  fMapTrackToCluster.clear();
159  fMapClusterToTrack.clear();
160  fNEntries = 1;
161  fVectorDeltaEtaDeltaPhi.clear();
162  fMap_TrID_ClID_ToIndex.clear();
163 
164  fSecMapTrackToCluster.clear();
165  fSecMapClusterToTrack.clear();
166  fSecNEntries = 1;
167  fSecVectorDeltaEtaDeltaPhi.clear();
168  fSecMap_TrID_ClID_ToIndex.clear();
169  fSecMap_TrID_ClID_AlreadyTried.clear();
170 
171  if(fRunNumber == -1 || fRunNumber != runNumber){
172  if(fClusterType == 1 || fClusterType == 3){
173  fGeomEMCAL = AliEMCALGeometry::GetInstanceFromRunNumber(runNumber);
174  if(!fGeomEMCAL){ AliFatal("EMCal geometry not initialized!");}
175  }else if(fClusterType == 2){
176  fGeomPHOS = AliPHOSGeometry::GetInstance();
177  if(!fGeomPHOS) AliFatal("PHOS geometry not initialized!");
178  }
179  fRunNumber = runNumber;
180  }
181 }
182 
183 //________________________________________________________________________
185 
186  // main method of AliCaloTrackMatcher, first initialize and then process event
187 
188  //DebugV0Matching();
189 
190  // do processing only for EMCal (1), DCal (3) or PHOS (2) clusters, otherwise do nothing
191  if(fClusterType == 1 || fClusterType == 2 || fClusterType == 3){
192  Initialize(fInputEvent->GetRunNumber());
193  ProcessEvent(fInputEvent);
194  }
195 
196  //DebugMatching();
197 
198  return;
199 }
200 
201 //________________________________________________________________________
202 void AliCaloTrackMatcher::ProcessEvent(AliVEvent *event){
203  Int_t nClus = 0;
204  TClonesArray * arrClusters = NULL;
205  if(!fCorrTaskSetting.CompareTo("")){
206  nClus = event->GetNumberOfCaloClusters();
207  } else {
208  arrClusters = dynamic_cast<TClonesArray*>(event->FindListObject(Form("%sClustersBranch",fCorrTaskSetting.Data())));
209  nClus = arrClusters->GetEntries();
210  }
211  Int_t nModules = 0;
212  if(fClusterType == 1 || fClusterType == 3) nModules = fGeomEMCAL->GetNumberOfSuperModules();
213  else if(fClusterType == 2) nModules = fGeomPHOS->GetNModules();
214 
215  AliESDEvent *esdev = dynamic_cast<AliESDEvent*>(event);
216  AliAODEvent *aodev = 0;
217  if (!esdev) {
218  aodev = dynamic_cast<AliAODEvent*>(event);
219  if (!aodev) {
220  AliError("Task needs AOD or ESD event, returning");
221  return;
222  }
223  }
224  static AliESDtrackCuts *EsdTrackCuts = 0x0;
225  static int prevRun = -1;
226  // Using standard function for setting Cuts
227  Int_t runNumber = fInputEvent->GetRunNumber();
228  if (esdev){
229  if (prevRun!=runNumber) {
230  delete EsdTrackCuts;
231  EsdTrackCuts = 0;
232  prevRun = runNumber;
233  }
234  if (!EsdTrackCuts) {
235  // if LHC11a or earlier or if LHC13g or if LHC12a-i -> use 2010 cuts
236  if( (runNumber<=146860) || (runNumber>=197470 && runNumber<=197692) || (runNumber>=172440 && runNumber<=193766) ){
237  EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010();
238  // else if run2 data use 2015 PbPb cuts
239  } else if (runNumber>=209122){
240  //EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2015PbPb();
241  // hard coded track cuts for the moment, because AliESDtrackCuts::GetStandardITSTPCTrackCuts2015PbPb() gives spams warnings
242  EsdTrackCuts = new AliESDtrackCuts();
243  // TPC; clusterCut = 1, cutAcceptanceEdges = kTRUE, removeDistortedRegions = kFALSE
244  EsdTrackCuts->AliESDtrackCuts::SetMinNCrossedRowsTPC(70);
245  EsdTrackCuts->AliESDtrackCuts::SetMinRatioCrossedRowsOverFindableClustersTPC(0.8);
246  EsdTrackCuts->SetCutGeoNcrNcl(2., 130., 1.5, 0.0, 0.0); // only dead zone and not clusters per length
247  //EsdTrackCuts->AliESDtrackCuts::SetCutOutDistortedRegionsTPC(kTRUE);
248  EsdTrackCuts->AliESDtrackCuts::SetMaxChi2PerClusterTPC(4);
249  EsdTrackCuts->AliESDtrackCuts::SetAcceptKinkDaughters(kFALSE);
250  EsdTrackCuts->AliESDtrackCuts::SetRequireTPCRefit(kTRUE);
251  // ITS; selPrimaries = 1
252  EsdTrackCuts->AliESDtrackCuts::SetRequireITSRefit(kTRUE);
253  EsdTrackCuts->AliESDtrackCuts::SetClusterRequirementITS(AliESDtrackCuts::kSPD,
254  AliESDtrackCuts::kAny);
255  EsdTrackCuts->AliESDtrackCuts::SetMaxDCAToVertexXYPtDep("0.0105+0.0350/pt^1.1");
256  EsdTrackCuts->AliESDtrackCuts::SetMaxChi2TPCConstrainedGlobal(36);
257  EsdTrackCuts->AliESDtrackCuts::SetMaxDCAToVertexZ(2);
258  EsdTrackCuts->AliESDtrackCuts::SetDCAToVertex2D(kFALSE);
259  EsdTrackCuts->AliESDtrackCuts::SetRequireSigmaToVertex(kFALSE);
260  EsdTrackCuts->AliESDtrackCuts::SetMaxChi2PerClusterITS(36);
261  // else use 2011 version of track cuts
262  }else{
263  EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2011();
264  }
265  EsdTrackCuts->SetMaxDCAToVertexZ(2);
266  }
267  }
268 
269  for (Int_t itr=0;itr<event->GetNumberOfTracks();itr++){
270  AliExternalTrackParam *trackParam = 0;
271  AliVTrack *inTrack = 0x0;
272  if(esdev){
273  inTrack = esdev->GetTrack(itr);
274  if(!inTrack) continue;
275  if(TMath::Abs(inTrack->Eta())>0.8 && (fClusterType==1 || fClusterType==3 )) continue;
276  if(TMath::Abs(inTrack->Eta())>0.3 && fClusterType==2 ) continue;
277  if(inTrack->Pt()<0.5) continue;
278  AliESDtrack *esdt = dynamic_cast<AliESDtrack*>(inTrack);
279  if(!EsdTrackCuts->AcceptTrack(esdt)) continue;
280  fHistControlMatches->Fill(0.,inTrack->Pt());
281 
282  const AliExternalTrackParam *in = esdt->GetInnerParam();
283  if (!in){AliDebug(2, "Could not get InnerParam of Track, continue"); fHistControlMatches->Fill(1.,inTrack->Pt()); continue;}
284  trackParam = new AliExternalTrackParam(*in);
285  } else if(aodev) {
286  inTrack = dynamic_cast<AliVTrack*>(aodev->GetTrack(itr));
287  if(!inTrack) continue;
288  if(inTrack->GetID()<0) continue; // Avoid double counting of tracks
289 
290  if(TMath::Abs(inTrack->Eta())>0.8 && (fClusterType==1 || fClusterType==3 )) continue;
291  if(TMath::Abs(inTrack->Eta())>0.3 && fClusterType==2 ) continue;
292  if(inTrack->Pt()<0.5) continue;
293 
294 
295  AliAODTrack *aodt = dynamic_cast<AliAODTrack*>(inTrack);
296  if(!aodt->IsHybridGlobalConstrainedGlobal()) continue;
297  Double_t xyz[3] = {0}, pxpypz[3] = {0}, cv[21] = {0};
298  aodt->GetPxPyPz(pxpypz);
299  aodt->GetXYZ(xyz);
300  aodt->GetCovarianceXYZPxPyPz(cv);
301  fHistControlMatches->Fill(0.,inTrack->Pt());
302 
303  // check for EMC tracks already propagated tracks are out of bounds
304  if (fClusterType == 1 || fClusterType == 3){
305  if( TMath::Abs(aodt->GetTrackEtaOnEMCal()) > 0.75 )
306  continue;
307 
308  // conditions for run1
309  if( fClusterType == 1 && nModules < 13 && ( aodt->GetTrackPhiOnEMCal() < 70*TMath::DegToRad() || aodt->GetTrackPhiOnEMCal() > 190*TMath::DegToRad()))
310  continue;
311 
312  // conditions for run2
313  if( nModules > 12 ){
314  if (fClusterType == 3 && ( aodt->GetTrackPhiOnEMCal() < 250*TMath::DegToRad() || aodt->GetTrackPhiOnEMCal() > 340*TMath::DegToRad()))
315  continue;
316  if( fClusterType == 1 && ( aodt->GetTrackPhiOnEMCal() < 70*TMath::DegToRad() || aodt->GetTrackPhiOnEMCal() > 190*TMath::DegToRad()))
317  continue;
318  }
319  }
320  trackParam = new AliExternalTrackParam(xyz,pxpypz,cv,aodt->Charge());
321  }
322 
323  if (!trackParam) {AliError("Could not get TrackParameters, continue"); fHistControlMatches->Fill(1.,inTrack->Pt()); continue;}
324  AliExternalTrackParam emcParam(*trackParam);
325  Float_t eta, phi, pt;
326 
327  //propagate tracks to emc surfaces
328  if(fClusterType == 1 || fClusterType == 3){
329  if (!AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(&emcParam, 440., 0.139, 20., eta, phi, pt)) {
330  delete trackParam;
331  fHistControlMatches->Fill(2.,inTrack->Pt());
332  continue;
333  }
334  if( TMath::Abs(eta) > 0.75 ) {
335  delete trackParam;
336  fHistControlMatches->Fill(3.,inTrack->Pt());
337  continue;
338  }
339  // Save some time and memory in case of no DCal present
340  if( fClusterType == 1 && nModules < 13 && ( phi < 70*TMath::DegToRad() || phi > 190*TMath::DegToRad())){
341  delete trackParam;
342  fHistControlMatches->Fill(3.,inTrack->Pt());
343  continue;
344  }
345  // Save some time and memory in case of run2
346  if( nModules > 12 ){
347  if (fClusterType == 3 && ( phi < 250*TMath::DegToRad() || phi > 340*TMath::DegToRad())){
348  delete trackParam;
349  fHistControlMatches->Fill(3.,inTrack->Pt());
350  continue;
351  }
352  if( fClusterType == 1 && ( phi < 70*TMath::DegToRad() || phi > 190*TMath::DegToRad())){
353  delete trackParam;
354  fHistControlMatches->Fill(3.,inTrack->Pt());
355  continue;
356  }
357  }
358  }else if(fClusterType == 2){
359  if( !AliTrackerBase::PropagateTrackToBxByBz(&emcParam, 460., 0.139, 20, kTRUE, 0.8, -1)){
360  delete trackParam;
361  fHistControlMatches->Fill(2.,inTrack->Pt());
362  continue;
363  }
364 
365  if (TMath::Abs(eta) > 0.25 && (fClusterType == 2)){
366  delete trackParam;
367  fHistControlMatches->Fill(3.,inTrack->Pt());
368  continue;
369  }
370  }
371 
372  Float_t dEta=-999, dPhi=-999;
373  Float_t clsPos[3] = {0.,0.,0.};
374  Double_t exPos[3] = {0.,0.,0.};
375  if (!emcParam.GetXYZ(exPos)){
376  fHistControlMatches->Fill(2.,inTrack->Pt());
377  delete trackParam;
378  continue;
379  }
380 
381  // cout << inTrack->GetID() << " - " << trackParam << endl;
382  // cout << "eta/phi: " << eta << ", " << phi << endl;
383  // cout << "nClus: " << nClus << endl;
384  Int_t nClusterMatchesToTrack = 0;
385  for(Int_t iclus=0;iclus < nClus;iclus++){
386  AliVCluster* cluster = NULL;
387  if(arrClusters){
388  if(esdev){
389  if(arrClusters)
390  cluster = new AliESDCaloCluster(*(AliESDCaloCluster*)arrClusters->At(iclus));
391  } else if(aodev){
392  if(arrClusters)
393  cluster = new AliAODCaloCluster(*(AliAODCaloCluster*)arrClusters->At(iclus));
394  }
395  }
396  else
397  cluster = event->GetCaloCluster(iclus);
398  if (!cluster){
399  if(arrClusters) delete cluster;
400  continue;
401  }
402  // cout << "-------------------------LOOPING: " << iclus << ", " << cluster->GetID() << endl;
403  cluster->GetPosition(clsPos);
404  Double_t dR = TMath::Sqrt(TMath::Power(exPos[0]-clsPos[0],2)+TMath::Power(exPos[1]-clsPos[1],2)+TMath::Power(exPos[2]-clsPos[2],2));
405  //cout << "dR: " << dR << endl;
406  if (dR > fMatchingWindow){
407  if(arrClusters) delete cluster;
408  continue;
409  }
410  Double_t clusterR = TMath::Sqrt( clsPos[0]*clsPos[0] + clsPos[1]*clsPos[1] );
411  AliExternalTrackParam trackParamTmp(emcParam);//Retrieve the starting point every time before the extrapolation
412  if(fClusterType == 1 || fClusterType == 3){
413  if (!cluster->IsEMCAL()){
414  if(arrClusters) delete cluster;
415  continue;
416  }
417  if(!AliEMCALRecoUtils::ExtrapolateTrackToCluster(&trackParamTmp, cluster, 0.139, 5., dEta, dPhi)){
418  fHistControlMatches->Fill(4.,inTrack->Pt());
419  if(arrClusters) delete cluster;
420  continue;
421  }
422  }else if(fClusterType == 2){
423  if (!cluster->IsPHOS()){
424  if(arrClusters) delete cluster;
425  continue;
426  }
427  if(!AliTrackerBase::PropagateTrackToBxByBz(&trackParamTmp, clusterR, 0.139, 5., kTRUE, 0.8, -1)){
428  fHistControlMatches->Fill(4.,inTrack->Pt());
429  if(arrClusters) delete cluster;
430  continue;
431  }
432  Double_t trkPos[3] = {0,0,0};
433  trackParamTmp.GetXYZ(trkPos);
434  TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
435  TVector3 clsPosVec(clsPos);
436  dPhi = clsPosVec.DeltaPhi(trkPosVec);
437  dEta = clsPosVec.Eta()-trkPosVec.Eta();
438  }
439 
440  Float_t dR2 = dPhi*dPhi + dEta*dEta;
441 
442  //cout << dEta << " - " << dPhi << " - " << dR2 << endl;
443  if(dR2 > fMatchingResidual){
444  if(arrClusters) delete cluster;
445  continue;
446  }
447  nClusterMatchesToTrack++;
448  if(aodev){
449  fMapTrackToCluster.insert(make_pair(itr,cluster->GetID()));
450  fMapClusterToTrack.insert(make_pair(cluster->GetID(),itr));
451  }else{
452  fMapTrackToCluster.insert(make_pair(inTrack->GetID(),cluster->GetID()));
453  fMapClusterToTrack.insert(make_pair(cluster->GetID(),inTrack->GetID()));
454  }
455  fVectorDeltaEtaDeltaPhi.push_back(make_pair(dEta,dPhi));
456  fMap_TrID_ClID_ToIndex[make_pair(inTrack->GetID(),cluster->GetID())] = fNEntries++;
457  if( (Int_t)fVectorDeltaEtaDeltaPhi.size() != (fNEntries-1)) AliFatal("Fatal error in AliCaloTrackMatcher, vector and map are not in sync!");
458  if(arrClusters) delete cluster;
459 
460  }
461  if(nClusterMatchesToTrack == 0) fHistControlMatches->Fill(5.,inTrack->Pt());
462  else fHistControlMatches->Fill(6.,inTrack->Pt());
463  delete trackParam;
464  }
465 
466  return;
467 }
468 
469 //________________________________________________________________________
470 Bool_t AliCaloTrackMatcher::PropagateV0TrackToClusterAndGetMatchingResidual(AliVTrack* inSecTrack, AliVCluster* cluster, AliVEvent* event, Float_t &dEta, Float_t &dPhi){
471 
472  //if V0-track to cluster match is already available return stored residuals
473  if(GetSecTrackClusterMatchingResidual(inSecTrack->GetID(),cluster->GetID(), dEta, dPhi)){
474  //cout << "RESIDUAL ALREADY AVAILABLE! - " << dEta << "/" << dPhi << endl;
475  return kTRUE;
476  }
477 
478  if(IsSecTrackClusterAlreadyTried(inSecTrack->GetID(),cluster->GetID())){
479  //cout << "PROPAGATION ALREADY FAILED! - " << inSecTrack->GetID() << "/" << cluster->GetID() << endl;
480  return kFALSE;
481  }
482 
483  //cout << "running matching! - " << inSecTrack->GetID() << "/" << cluster->GetID() << endl;
484  //if match has not yet been computed, go on:
485  Int_t nModules = 0;
486  if(fClusterType == 1 || fClusterType == 3) nModules = fGeomEMCAL->GetNumberOfSuperModules();
487  else if(fClusterType == 2) nModules = fGeomPHOS->GetNModules();
488 
489  AliESDEvent *esdev = dynamic_cast<AliESDEvent*>(event);
490  AliAODEvent *aodev = 0;
491  if (!esdev) {
492  aodev = dynamic_cast<AliAODEvent*>(event);
493  if (!aodev) {
494  AliError("Task needs AOD or ESD event, returning");
495  return kFALSE;
496  }
497  }
498 
499  if(!cluster->IsEMCAL() && !cluster->IsPHOS()){AliError("Cluster is neither EMCAL nor PHOS, returning");return kFALSE;}
500 
501  Float_t clusterPosition[3] = {0,0,0};
502  cluster->GetPosition(clusterPosition);
503  Double_t clusterR = TMath::Sqrt( clusterPosition[0]*clusterPosition[0] + clusterPosition[1]*clusterPosition[1] );
504 
505  if(!inSecTrack) return kFALSE;
506  fSecHistControlMatches->Fill(0.,inSecTrack->Pt());
507 
508  AliESDtrack *esdt = dynamic_cast<AliESDtrack*>(inSecTrack);
509  AliAODTrack *aodt = 0;
510  if (!esdt) {
511  aodt = dynamic_cast<AliAODTrack*>(inSecTrack);
512  if (!aodt){
513  AliError("Track is neither ESD nor AOD, continue");
514  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
515  return kFALSE;
516  }
517  }
518 
519  AliExternalTrackParam *trackParam = 0;
520  if (esdt) {
521  const AliExternalTrackParam *in = esdt->GetInnerParam();
522  if (!in){
523  AliDebug(2, "Could not get InnerParam of Track, continue");
524  fSecHistControlMatches->Fill(1.,inSecTrack->Pt());
525  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
526  return kFALSE;
527  }
528  trackParam = new AliExternalTrackParam(*in);
529  } else {
530  // check if tracks should be propagated at all
531  if (fClusterType == 1 || fClusterType == 3){
532  if (TMath::Abs(aodt->GetTrackEtaOnEMCal()) > 0.8)
533  return kFALSE;
534  if( nModules < 13 ){
535  if (( aodt->GetTrackPhiOnEMCal() < 60*TMath::DegToRad() || aodt->GetTrackPhiOnEMCal() > 200*TMath::DegToRad()))
536  return kFALSE;
537  } else if( nModules > 12 ){
538  if (fClusterType == 3 && ( aodt->GetTrackPhiOnEMCal() < 250*TMath::DegToRad() || aodt->GetTrackPhiOnEMCal() > 340*TMath::DegToRad()))
539  return kFALSE;
540  if( fClusterType == 1 && ( aodt->GetTrackPhiOnEMCal() < 60*TMath::DegToRad() || aodt->GetTrackPhiOnEMCal() > 200*TMath::DegToRad()))
541  return kFALSE;
542  }
543  } else {
544  if ( aodt->Phi() < 60*TMath::DegToRad() || aodt->Phi() > 200*TMath::DegToRad())
545  return kFALSE;
546  if (TMath::Abs(aodt->Eta()) > 0.3 )
547  return kFALSE;
548  }
549 
550  Double_t xyz[3] = {0}, pxpypz[3] = {0}, cv[21] = {0};
551  aodt->GetPxPyPz(pxpypz);
552  aodt->GetXYZ(xyz);
553  aodt->GetCovarianceXYZPxPyPz(cv);
554  trackParam = new AliExternalTrackParam(xyz,pxpypz,cv,aodt->Charge());
555  }
556  if(!trackParam){
557  AliError("Could not get TrackParameters, continue");
558  fSecHistControlMatches->Fill(1.,inSecTrack->Pt());
559  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
560  return kFALSE;
561  }
562 
563  Bool_t propagated = kFALSE;
564  AliExternalTrackParam emcParam(*trackParam);
565  Float_t dPhiTemp = 0;
566  Float_t dEtaTemp = 0;
567 
568  if(cluster->IsEMCAL()){
569  if (inSecTrack->Pt() < 0.3 ) {
570  fSecHistControlMatches->Fill(3.,inSecTrack->Pt());
571  }
572  Float_t eta = 0;Float_t phi = 0;Float_t pt = 0;
573  propagated = AliEMCALRecoUtils::ExtrapolateTrackToEMCalSurface(&emcParam, 430, 0.000510999, 20, eta, phi, pt);
574  if(propagated){
575  if( TMath::Abs(eta) > 0.8 ) {
576  delete trackParam;
577  fSecHistControlMatches->Fill(3.,inSecTrack->Pt());
578  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
579  return kFALSE;
580  }
581  // Save some time and memory in case of no DCal present
582  if( nModules < 13 && ( phi < 60*TMath::DegToRad() || phi > 200*TMath::DegToRad())){
583  delete trackParam;
584  fSecHistControlMatches->Fill(3.,inSecTrack->Pt());
585  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
586  return kFALSE;
587  }
588 
589  propagated = AliEMCALRecoUtils::ExtrapolateTrackToCluster(&emcParam, cluster, 0.000510999, 5, dEtaTemp, dPhiTemp);
590  if(!propagated){
591  delete trackParam;
592  fSecHistControlMatches->Fill(4.,inSecTrack->Pt());
593  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
594  return kFALSE;
595  }
596  }else{
597  delete trackParam;
598  fSecHistControlMatches->Fill(2.,inSecTrack->Pt());
599  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
600  return kFALSE;
601  }
602 
603  }else if(cluster->IsPHOS()){
604  propagated = AliTrackerBase::PropagateTrackToBxByBz(&emcParam, clusterR, 0.000510999, 20, kTRUE, 0.8, -1);
605  if (propagated){
606  Double_t trkPos[3] = {0,0,0};
607  emcParam.GetXYZ(trkPos);
608  TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
609  TVector3 clsPosVec(clusterPosition);
610  dPhiTemp = clsPosVec.DeltaPhi(trkPosVec);
611  dEtaTemp = clsPosVec.Eta()-trkPosVec.Eta();
612  }else{
613  delete trackParam;
614  fSecHistControlMatches->Fill(2.,inSecTrack->Pt());
615  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
616  return kFALSE;}
617  }
618 
619  if (propagated){
620  Float_t dR2 = dPhiTemp*dPhiTemp + dEtaTemp*dEtaTemp;
621  //cout << dEtaTemp << " - " << dPhiTemp << " - " << dR2 << endl;
622  if(dR2 > fMatchingResidual){
623  fSecHistControlMatches->Fill(5.,inSecTrack->Pt());
624  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
625  //cout << "NO MATCH! - " << inSecTrack->GetID() << "/" << cluster->GetID() << endl;
626  delete trackParam;
627  return kFALSE;
628  }
629  //cout << "MATCHED!!!!!!!" << endl;
630 
631  if(aodev){
632  //need to search for position in case of AOD
633  Int_t TrackPos = -1;
634  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
635  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
636  if(currTrack->GetID() == inSecTrack->GetID()){
637  TrackPos = iTrack;
638  break;
639  }
640  }
641  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: PropagateV0TrackToClusterAndGetMatchingResidual - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",inSecTrack->GetID()));
642  fSecMapTrackToCluster.insert(make_pair(TrackPos,cluster->GetID()));
643  fSecMapClusterToTrack.insert(make_pair(cluster->GetID(),TrackPos));
644  }else{
645  fSecMapTrackToCluster.insert(make_pair(inSecTrack->GetID(),cluster->GetID()));
646  fSecMapClusterToTrack.insert(make_pair(cluster->GetID(),inSecTrack->GetID()));
647  }
648  fSecVectorDeltaEtaDeltaPhi.push_back(make_pair(dEtaTemp,dPhiTemp));
649  fSecMap_TrID_ClID_ToIndex[make_pair(inSecTrack->GetID(),cluster->GetID())] = fSecNEntries++;
650  if( (Int_t)fSecVectorDeltaEtaDeltaPhi.size() != (fSecNEntries-1)) AliFatal("Fatal error in AliCaloTrackMatcher, vector and map are not in sync!");
651 
652  fSecHistControlMatches->Fill(6.,inSecTrack->Pt());
653  dEta = dEtaTemp;
654  dPhi = dPhiTemp;
655  delete trackParam;
656  return kTRUE;
657  }else AliFatal("Fatal error in AliCaloTrackMatcher, track is labeled as sucessfully propagated although this should be impossible!");
658 
659  fSecMap_TrID_ClID_AlreadyTried[make_pair(inSecTrack->GetID(),cluster->GetID())] = 1.;
660  delete trackParam;
661  return kFALSE;
662 }
663 
664 //________________________________________________________________________
665 //________________________________________________________________________
666 //________________________________________________________________________
667 //________________________________________________________________________
668 //________________________________________________________________________
670  Int_t position = fMap_TrID_ClID_ToIndex[make_pair(trackID,clusterID)];
671  if(position == 0) return kFALSE;
672 
673  pairFloat tempEtaPhi = fVectorDeltaEtaDeltaPhi.at(position-1);
674  dEta = tempEtaPhi.first;
675  dPhi = tempEtaPhi.second;
676  return kTRUE;
677 }
678 //________________________________________________________________________
679 Int_t AliCaloTrackMatcher::GetNMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
680  Int_t matched = 0;
681  multimap<Int_t,Int_t>::iterator it;
682  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it){
683  if(it->first == clusterID){
684  Float_t tempDEta, tempDPhi;
685  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
686  if(!tempTrack) continue;
687  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
688  if(tempTrack->Charge()>0){
689  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) matched++;
690  }else if(tempTrack->Charge()<0){
691  dPhiMin*=-1;
692  dPhiMax*=-1;
693  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) matched++;
694  }
695  }
696  }
697  }
698 
699  return matched;
700 }
701 
702 //________________________________________________________________________
703 Int_t AliCaloTrackMatcher::GetNMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
704  Int_t matched = 0;
705  multimap<Int_t,Int_t>::iterator it;
706  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it){
707  if(it->first == clusterID){
708  Float_t tempDEta, tempDPhi;
709  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
710  if(!tempTrack) continue;
711  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
712  Bool_t match_dEta = kFALSE;
713  Bool_t match_dPhi = kFALSE;
714  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
715  else match_dEta = kFALSE;
716 
717  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
718  else match_dPhi = kFALSE;
719 
720  if (match_dPhi && match_dEta )matched++;
721  }
722  }
723  }
724  return matched;
725 }
726 
727 //________________________________________________________________________
729  Int_t matched = 0;
730  multimap<Int_t,Int_t>::iterator it;
731  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it){
732  if(it->first == clusterID){
733  Float_t tempDEta, tempDPhi;
734  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
735  if(!tempTrack) continue;
736  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
737  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) matched++;
738  }
739  }
740  }
741  return matched;
742 }
743 
744 //________________________________________________________________________
745 Int_t AliCaloTrackMatcher::GetNMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
746 
747  Int_t TrackPos = -1;
748  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
749  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
750  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
751  if(currTrack->GetID() == trackID){
752  TrackPos = iTrack;
753  break;
754  }
755  }
756  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
757  }else TrackPos = trackID; // for ESD just take trackID
758 
759  Int_t matched = 0;
760  multimap<Int_t,Int_t>::iterator it;
761  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
762  if(!tempTrack) return matched;
763  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it){
764  if(it->first == TrackPos){
765  Float_t tempDEta, tempDPhi;
766  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
767  if(tempTrack->Charge()>0){
768  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) matched++;
769  }else if(tempTrack->Charge()<0){
770  dPhiMin*=-1;
771  dPhiMax*=-1;
772  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) matched++;
773  }
774  }
775  }
776  }
777  return matched;
778 }
779 
780 //________________________________________________________________________
781 Int_t AliCaloTrackMatcher::GetNMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
782  Int_t TrackPos = -1;
783  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
784  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
785  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
786  if(currTrack->GetID() == trackID){
787  TrackPos = iTrack;
788  break;
789  }
790  }
791  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
792  }else TrackPos = trackID; // for ESD just take trackID
793 
794  Int_t matched = 0;
795  multimap<Int_t,Int_t>::iterator it;
796  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
797  if(!tempTrack) return matched;
798  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it){
799  if(it->first == TrackPos){
800  Float_t tempDEta, tempDPhi;
801  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
802  Bool_t match_dEta = kFALSE;
803  Bool_t match_dPhi = kFALSE;
804  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
805  else match_dEta = kFALSE;
806 
807  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
808  else match_dPhi = kFALSE;
809 
810  if (match_dPhi && match_dEta )matched++;
811 
812  }
813  }
814  }
815  return matched;
816 }
817 
818 //________________________________________________________________________
820  Int_t TrackPos = -1;
821  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
822  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
823  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
824  if(currTrack->GetID() == trackID){
825  TrackPos = iTrack;
826  break;
827  }
828  }
829  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
830  }else TrackPos = trackID; // for ESD just take trackID
831 
832  Int_t matched = 0;
833  multimap<Int_t,Int_t>::iterator it;
834  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
835  if(!tempTrack) return matched;
836  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it){
837  if(it->first == TrackPos){
838  Float_t tempDEta, tempDPhi;
839  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
840  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) matched++;
841  }
842  }
843  }
844  return matched;
845 }
846 
847 //________________________________________________________________________
848 vector<Int_t> AliCaloTrackMatcher::GetMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
849  vector<Int_t> tempMatchedTracks;
850  multimap<Int_t,Int_t>::iterator it;
851  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it){
852  if(it->first == clusterID){
853  Float_t tempDEta, tempDPhi;
854  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
855  if(!tempTrack) continue;
856  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
857  if(tempTrack->Charge()>0){
858  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) tempMatchedTracks.push_back(it->second);
859  }else if(tempTrack->Charge()<0){
860  dPhiMin*=-1;
861  dPhiMax*=-1;
862  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) tempMatchedTracks.push_back(it->second);
863  }
864  }
865  }
866  }
867  return tempMatchedTracks;
868 }
869 
870 //________________________________________________________________________
871 vector<Int_t> AliCaloTrackMatcher::GetMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
872  vector<Int_t> tempMatchedTracks;
873  multimap<Int_t,Int_t>::iterator it;
874  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it){
875  if(it->first == clusterID){
876  Float_t tempDEta, tempDPhi;
877  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
878  if(!tempTrack) continue;
879  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
880  Bool_t match_dEta = kFALSE;
881  Bool_t match_dPhi = kFALSE;
882  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
883  else match_dEta = kFALSE;
884 
885  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
886  else match_dPhi = kFALSE;
887 
888  if (match_dPhi && match_dEta )tempMatchedTracks.push_back(it->second);
889 
890  }
891  }
892  }
893  return tempMatchedTracks;
894 }
895 
896 //________________________________________________________________________
897 vector<Int_t> AliCaloTrackMatcher::GetMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dR){
898  vector<Int_t> tempMatchedTracks;
899  multimap<Int_t,Int_t>::iterator it;
900  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it){
901  if(it->first == clusterID){
902  Float_t tempDEta, tempDPhi;
903  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
904  if(!tempTrack) continue;
905  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
906  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) tempMatchedTracks.push_back(it->second);
907  }
908  }
909  }
910  return tempMatchedTracks;
911 }
912 
913 //________________________________________________________________________
914 vector<Int_t> AliCaloTrackMatcher::GetMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
915  Int_t TrackPos = -1;
916  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
917  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
918  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
919  if(currTrack->GetID() == trackID){
920  TrackPos = iTrack;
921  break;
922  }
923  }
924  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
925  }else TrackPos = trackID; // for ESD just take trackID
926 
927  vector<Int_t> tempMatchedClusters;
928  multimap<Int_t,Int_t>::iterator it;
929  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
930  if(!tempTrack) return tempMatchedClusters;
931  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it){
932  if(it->first == TrackPos){
933  Float_t tempDEta, tempDPhi;
934  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
935  if(tempTrack->Charge()>0){
936  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) tempMatchedClusters.push_back(it->second);
937  }else if(tempTrack->Charge()<0){
938  dPhiMin*=-1;
939  dPhiMax*=-1;
940  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) tempMatchedClusters.push_back(it->second);
941  }
942  }
943  }
944  }
945 
946  return tempMatchedClusters;
947 }
948 
949 //________________________________________________________________________
950 vector<Int_t> AliCaloTrackMatcher::GetMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
951  Int_t TrackPos = -1;
952  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
953  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
954  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
955  if(currTrack->GetID() == trackID){
956  TrackPos = iTrack;
957  break;
958  }
959  }
960  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
961  }else TrackPos = trackID; // for ESD just take trackID
962 
963  vector<Int_t> tempMatchedClusters;
964  multimap<Int_t,Int_t>::iterator it;
965  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
966  if(!tempTrack) return tempMatchedClusters;
967  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it){
968  if(it->first == TrackPos){
969  Float_t tempDEta, tempDPhi;
970  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
971  Bool_t match_dEta = kFALSE;
972  Bool_t match_dPhi = kFALSE;
973  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
974  else match_dEta = kFALSE;
975 
976  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
977  else match_dPhi = kFALSE;
978 
979  if (match_dPhi && match_dEta )tempMatchedClusters.push_back(it->second);
980  }
981  }
982  }
983  return tempMatchedClusters;
984 }
985 
986 //________________________________________________________________________
987 vector<Int_t> AliCaloTrackMatcher::GetMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, Float_t dR){
988  Int_t TrackPos = -1;
989  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
990  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
991  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
992  if(currTrack->GetID() == trackID){
993  TrackPos = iTrack;
994  break;
995  }
996  }
997  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
998  }else TrackPos = trackID; // for ESD just take trackID
999 
1000  vector<Int_t> tempMatchedClusters;
1001  multimap<Int_t,Int_t>::iterator it;
1002  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1003  if(!tempTrack) return tempMatchedClusters;
1004  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it){
1005  if(it->first == TrackPos){
1006  Float_t tempDEta, tempDPhi;
1007  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1008  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) tempMatchedClusters.push_back(it->second);
1009  }
1010  }
1011  }
1012  return tempMatchedClusters;
1013 }
1014 
1015 //________________________________________________________________________
1016 //________________________________________________________________________
1017 //________________________________________________________________________
1018 //________________________________________________________________________
1019 //________________________________________________________________________
1021  Int_t position = fSecMap_TrID_ClID_ToIndex[make_pair(trackID,clusterID)];
1022  if(position == 0) return kFALSE;
1023 
1024  pairFloat tempEtaPhi = fSecVectorDeltaEtaDeltaPhi.at(position-1);
1025  dEta = tempEtaPhi.first;
1026  dPhi = tempEtaPhi.second;
1027  return kTRUE;
1028 }
1029 //________________________________________________________________________
1031  Int_t position = fSecMap_TrID_ClID_AlreadyTried[make_pair(trackID,clusterID)];
1032  if(position == 0) return kFALSE;
1033  else return kTRUE;
1034 }
1035 //________________________________________________________________________
1036 Int_t AliCaloTrackMatcher::GetNMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
1037  Int_t matched = 0;
1038  multimap<Int_t,Int_t>::iterator it;
1039  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it){
1040  if(it->first == clusterID){
1041  Float_t tempDEta, tempDPhi;
1042  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
1043  if(!tempTrack) continue;
1044  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
1045  if(tempTrack->Charge()>0){
1046  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) matched++;
1047  }else if(tempTrack->Charge()<0){
1048  dPhiMin*=-1;
1049  dPhiMax*=-1;
1050  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) matched++;
1051  }
1052  }
1053  }
1054  }
1055 
1056  return matched;
1057 }
1058 
1059 //________________________________________________________________________
1060 Int_t AliCaloTrackMatcher::GetNMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t clusterID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
1061  Int_t matched = 0;
1062  multimap<Int_t,Int_t>::iterator it;
1063  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it){
1064  if(it->first == clusterID){
1065  Float_t tempDEta, tempDPhi;
1066  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
1067  if(!tempTrack) continue;
1068  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
1069  Bool_t match_dEta = kFALSE;
1070  Bool_t match_dPhi = kFALSE;
1071  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
1072  else match_dEta = kFALSE;
1073 
1074  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
1075  else match_dPhi = kFALSE;
1076 
1077  if (match_dPhi && match_dEta )matched++;
1078  }
1079  }
1080  }
1081 
1082  return matched;
1083 }
1084 
1085 //________________________________________________________________________
1087  Int_t matched = 0;
1088  multimap<Int_t,Int_t>::iterator it;
1089  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it){
1090  if(it->first == clusterID){
1091  Float_t tempDEta, tempDPhi;
1092  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
1093  if(!tempTrack) continue;
1094  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
1095  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) matched++;
1096  }
1097  }
1098  }
1099 
1100  return matched;
1101 }
1102 
1103 //________________________________________________________________________
1104 Int_t AliCaloTrackMatcher::GetNMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
1105  Int_t TrackPos = -1;
1106  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
1107  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
1108  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
1109  if(currTrack->GetID() == trackID){
1110  TrackPos = iTrack;
1111  break;
1112  }
1113  }
1114  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
1115  }else TrackPos = trackID; // for ESD just take trackID
1116 
1117  Int_t matched = 0;
1118  multimap<Int_t,Int_t>::iterator it;
1119  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1120  if(!tempTrack) return matched;
1121  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it){
1122  if(it->first == TrackPos){
1123  Float_t tempDEta, tempDPhi;
1124  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1125  if(tempTrack->Charge()>0){
1126  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) matched++;
1127  }else if(tempTrack->Charge()<0){
1128  dPhiMin*=-1;
1129  dPhiMax*=-1;
1130  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) matched++;
1131  }
1132  }
1133  }
1134  }
1135 
1136  return matched;
1137 }
1138 
1139 //________________________________________________________________________
1140 Int_t AliCaloTrackMatcher::GetNMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t trackID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
1141  Int_t TrackPos = -1;
1142  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
1143  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
1144  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
1145  if(currTrack->GetID() == trackID){
1146  TrackPos = iTrack;
1147  break;
1148  }
1149  }
1150  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
1151  }else TrackPos = trackID; // for ESD just take trackID
1152 
1153  Int_t matched = 0;
1154  multimap<Int_t,Int_t>::iterator it;
1155  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1156  if(!tempTrack) return matched;
1157  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it){
1158  if(it->first == TrackPos){
1159  Float_t tempDEta, tempDPhi;
1160  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1161  Bool_t match_dEta = kFALSE;
1162  Bool_t match_dPhi = kFALSE;
1163  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
1164  else match_dEta = kFALSE;
1165 
1166  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
1167  else match_dPhi = kFALSE;
1168 
1169  if (match_dPhi && match_dEta )matched++;
1170 
1171  }
1172  }
1173  }
1174 
1175  return matched;
1176 }
1177 
1178 //________________________________________________________________________
1180  Int_t TrackPos = -1;
1181  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
1182  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
1183  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
1184  if(currTrack->GetID() == trackID){
1185  TrackPos = iTrack;
1186  break;
1187  }
1188  }
1189  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
1190  }else TrackPos = trackID; // for ESD just take trackID
1191 
1192  Int_t matched = 0;
1193  multimap<Int_t,Int_t>::iterator it;
1194  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1195  if(!tempTrack) return matched;
1196  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it){
1197  if(it->first == TrackPos){
1198  Float_t tempDEta, tempDPhi;
1199  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1200  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) matched++;
1201  }
1202  }
1203  }
1204 
1205  return matched;
1206 }
1207 
1208 //________________________________________________________________________
1209 vector<Int_t> AliCaloTrackMatcher::GetMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
1210  vector<Int_t> tempMatchedTracks;
1211  multimap<Int_t,Int_t>::iterator it;
1212  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it){
1213  if(it->first == clusterID){
1214  Float_t tempDEta, tempDPhi;
1215  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
1216  if(!tempTrack) continue;
1217  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
1218  if(tempTrack->Charge()>0){
1219  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) tempMatchedTracks.push_back(it->second);
1220  }else if(tempTrack->Charge()<0){
1221  dPhiMin*=-1;
1222  dPhiMax*=-1;
1223  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) tempMatchedTracks.push_back(it->second);
1224  }
1225  }
1226  }
1227  }
1228 
1229  return tempMatchedTracks;
1230 }
1231 
1232 //________________________________________________________________________
1233 vector<Int_t> AliCaloTrackMatcher::GetMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t clusterID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
1234  vector<Int_t> tempMatchedTracks;
1235  multimap<Int_t,Int_t>::iterator it;
1236  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it){
1237  if(it->first == clusterID){
1238  Float_t tempDEta, tempDPhi;
1239  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
1240  if(!tempTrack) continue;
1241  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
1242  Bool_t match_dEta = kFALSE;
1243  Bool_t match_dPhi = kFALSE;
1244  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
1245  else match_dEta = kFALSE;
1246 
1247  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
1248  else match_dPhi = kFALSE;
1249 
1250  if (match_dPhi && match_dEta )tempMatchedTracks.push_back(it->second);
1251  }
1252  }
1253  }
1254 
1255  return tempMatchedTracks;
1256 }
1257 
1258 //________________________________________________________________________
1259 vector<Int_t> AliCaloTrackMatcher::GetMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dR){
1260  vector<Int_t> tempMatchedTracks;
1261  multimap<Int_t,Int_t>::iterator it;
1262  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it){
1263  if(it->first == clusterID){
1264  Float_t tempDEta, tempDPhi;
1265  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(it->second));
1266  if(!tempTrack) continue;
1267  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->first,tempDEta,tempDPhi)){
1268  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) tempMatchedTracks.push_back(it->second);
1269  }
1270  }
1271  }
1272 
1273  return tempMatchedTracks;
1274 }
1275 
1276 //________________________________________________________________________
1277 vector<Int_t> AliCaloTrackMatcher::GetMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin){
1278  Int_t TrackPos = -1;
1279  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
1280  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
1281  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
1282  if(currTrack->GetID() == trackID){
1283  TrackPos = iTrack;
1284  break;
1285  }
1286  }
1287  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
1288  }else TrackPos = trackID; // for ESD just take trackID
1289 
1290  vector<Int_t> tempMatchedClusters;
1291  multimap<Int_t,Int_t>::iterator it;
1292  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1293  if(!tempTrack) return tempMatchedClusters;
1294  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it){
1295  if(it->first == TrackPos){
1296  Float_t tempDEta, tempDPhi;
1297  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1298  if(tempTrack->Charge()>0){
1299  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin < tempDPhi) && (tempDPhi < dPhiMax) ) tempMatchedClusters.push_back(it->second);
1300  }else if(tempTrack->Charge()<0){
1301  dPhiMin*=-1;
1302  dPhiMax*=-1;
1303  if( (dEtaMin < tempDEta) && (tempDEta < dEtaMax) && (dPhiMin > tempDPhi) && (tempDPhi > dPhiMax) ) tempMatchedClusters.push_back(it->second);
1304  }
1305  }
1306  }
1307  }
1308 
1309  return tempMatchedClusters;
1310 }
1311 
1312 //________________________________________________________________________
1313 vector<Int_t> AliCaloTrackMatcher::GetMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t trackID, TF1* fFuncPtDepEta, TF1* fFuncPtDepPhi){
1314  Int_t TrackPos = -1;
1315  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
1316  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
1317  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
1318  if(currTrack->GetID() == trackID){
1319  TrackPos = iTrack;
1320  break;
1321  }
1322  }
1323  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
1324  }else TrackPos = trackID; // for ESD just take trackID
1325 
1326  vector<Int_t> tempMatchedClusters;
1327  multimap<Int_t,Int_t>::iterator it;
1328  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1329  if(!tempTrack) return tempMatchedClusters;
1330  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it){
1331  if(it->first == TrackPos){
1332  Float_t tempDEta, tempDPhi;
1333  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1334  Bool_t match_dEta = kFALSE;
1335  Bool_t match_dPhi = kFALSE;
1336  if( TMath::Abs(tempDEta) < fFuncPtDepEta->Eval(tempTrack->Pt())) match_dEta = kTRUE;
1337  else match_dEta = kFALSE;
1338 
1339  if( TMath::Abs(tempDPhi) < fFuncPtDepPhi->Eval(tempTrack->Pt())) match_dPhi = kTRUE;
1340  else match_dPhi = kFALSE;
1341 
1342  if (match_dPhi && match_dEta )tempMatchedClusters.push_back(it->second);
1343  }
1344  }
1345  }
1346 
1347  return tempMatchedClusters;
1348 }
1349 
1350 //________________________________________________________________________
1351 vector<Int_t> AliCaloTrackMatcher::GetMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t trackID, Float_t dR){
1352  Int_t TrackPos = -1;
1353  if(event->IsA()==AliAODEvent::Class()){ // for AOD, we have to look for position of track in the event
1354  for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++){
1355  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(iTrack));
1356  if(currTrack->GetID() == trackID){
1357  TrackPos = iTrack;
1358  break;
1359  }
1360  }
1361  if(TrackPos == -1) AliFatal(Form("AliCaloTrackMatcher: GetNMatchedClusterIDsForTrack - track (ID: '%i') cannot be retrieved from event, should be impossible as it has been used in maim task before!",trackID));
1362  }else TrackPos = trackID; // for ESD just take trackID
1363 
1364  vector<Int_t> tempMatchedClusters;
1365  multimap<Int_t,Int_t>::iterator it;
1366  AliVTrack* tempTrack = dynamic_cast<AliVTrack*>(event->GetTrack(TrackPos));
1367  if(!tempTrack) return tempMatchedClusters;
1368  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it){
1369  if(it->first == TrackPos){
1370  Float_t tempDEta, tempDPhi;
1371  if(GetTrackClusterMatchingResidual(tempTrack->GetID(),it->second,tempDEta,tempDPhi)){
1372  if (TMath::Sqrt(tempDEta*tempDEta + tempDPhi*tempDPhi) < dR ) tempMatchedClusters.push_back(it->second);
1373  }
1374  }
1375  }
1376 
1377  return tempMatchedClusters;
1378 }
1379 
1380 //________________________________________________________________________
1382  Float_t sumTrackEt = 0.;
1383  vector<Int_t> labelsMatched = GetMatchedTrackIDsForCluster(event, clusterID, dR);
1384  if((Int_t) labelsMatched.size()<1) return sumTrackEt;
1385 
1386  TLorentzVector vecTrack;
1387  for (UInt_t i = 0; i < labelsMatched.size(); i++){
1388  AliVTrack* currTrack = dynamic_cast<AliVTrack*>(event->GetTrack(labelsMatched.at(i)));
1389  if(!currTrack) continue;
1390  vecTrack.SetPxPyPzE(currTrack->Px(),currTrack->Py(),currTrack->Pz(),currTrack->E());
1391  sumTrackEt += vecTrack.Et();
1392  }
1393 
1394  return sumTrackEt;
1395 }
1396 
1397 //________________________________________________________________________
1399  TAxis *axisafter = histoRebin->GetYaxis();
1400  Int_t bins = axisafter->GetNbins();
1401  Double_t from = axisafter->GetXmin();
1402  Double_t to = axisafter->GetXmax();
1403  Double_t *newbins = new Double_t[bins+1];
1404  newbins[0] = from;
1405  Double_t factor = TMath::Power(to/from, 1./bins);
1406  for(Int_t i=1; i<=bins; ++i) newbins[i] = factor * newbins[i-1];
1407  axisafter->Set(bins, newbins);
1408  delete [] newbins;
1409  return;
1410 }
1411 
1412 //________________________________________________________________________
1414  if(fSecVectorDeltaEtaDeltaPhi.size()>0){
1415  cout << "******************************" << endl;
1416  cout << "******************************" << endl;
1417  cout << "NEW EVENT !" << endl;
1418  cout << "vector etaphi:" << endl;
1419  cout << fSecVectorDeltaEtaDeltaPhi.size() << endl;
1420  cout << "multimap" << endl;
1421  mapT::iterator iter = fSecMap_TrID_ClID_ToIndex.begin();
1422  for (iter = fSecMap_TrID_ClID_ToIndex.begin(); iter != fSecMap_TrID_ClID_ToIndex.end(); ++iter){
1423  Float_t dEta, dPhi = 0;
1424  if(!GetSecTrackClusterMatchingResidual(iter->first.first,iter->first.second,dEta,dPhi)) continue;
1425  cout << " [" << iter->first.first << "/" << iter->first.second << ", " << iter->second << "] - (" << dEta << "/" << dPhi << ")" << endl;
1426  }
1427  cout << "mapTrackToCluster" << endl;
1428  AliESDEvent *esdev = dynamic_cast<AliESDEvent*>(fInputEvent);
1429  for (Int_t itr=0;itr<esdev->GetNumberOfTracks();itr++){
1430  AliVTrack *inTrack = esdev->GetTrack(itr);
1431  if(!inTrack) continue;
1432  TString tCharge;
1433  if(inTrack->Charge()>0) tCharge = "+";
1434  else if(inTrack->Charge()<0) tCharge = "-";
1435  cout << itr << " (" << tCharge << ") - " << GetNMatchedClusterIDsForSecTrack(fInputEvent,inTrack->GetID(),5,-5,0.2,-0.4) << "\t\t";
1436  }
1437  cout << endl;
1438  multimap<Int_t,Int_t>::iterator it;
1439  for (it=fSecMapTrackToCluster.begin(); it!=fSecMapTrackToCluster.end(); ++it) cout << it->first << " => " << it->second << '\n';
1440  cout << "mapClusterToTrack" << endl;
1441  Int_t tempClus = it->second;
1442  for (it=fSecMapClusterToTrack.begin(); it!=fSecMapClusterToTrack.end(); ++it) cout << it->first << " => " << it->second << '\n';
1443  vector<Int_t> tempTracks = GetMatchedSecTrackIDsForCluster(fInputEvent,tempClus, 5, -5, 0.2, -0.4);
1444  for(UInt_t iJ=0; iJ<tempTracks.size();iJ++){
1445  cout << tempClus << " - " << tempTracks.at(iJ) << endl;
1446  }
1447  }
1448  return;
1449 }
1450 
1451 //________________________________________________________________________
1453  if(fVectorDeltaEtaDeltaPhi.size()>0){
1454  cout << "******************************" << endl;
1455  cout << "******************************" << endl;
1456  cout << "NEW EVENT !" << endl;
1457  cout << "vector etaphi:" << endl;
1458  cout << fVectorDeltaEtaDeltaPhi.size() << endl;
1459  cout << "multimap" << endl;
1460  mapT::iterator iter = fMap_TrID_ClID_ToIndex.begin();
1461  for (iter = fMap_TrID_ClID_ToIndex.begin(); iter != fMap_TrID_ClID_ToIndex.end(); ++iter){
1462  Float_t dEta, dPhi = 0;
1463  if(!GetTrackClusterMatchingResidual(iter->first.first,iter->first.second,dEta,dPhi)) continue;
1464  cout << " [" << iter->first.first << "/" << iter->first.second << ", " << iter->second << "] - (" << dEta << "/" << dPhi << ")" << endl;
1465  }
1466  cout << "mapTrackToCluster" << endl;
1467  AliESDEvent *esdev = dynamic_cast<AliESDEvent*>(fInputEvent);
1468  for (Int_t itr=0;itr<esdev->GetNumberOfTracks();itr++){
1469  AliVTrack *inTrack = esdev->GetTrack(itr);
1470  if(!inTrack) continue;
1471  TString tCharge;
1472  if(inTrack->Charge()>0) tCharge = "+";
1473  else if(inTrack->Charge()<0) tCharge = "-";
1474  cout << itr << " (" << tCharge << ") - " << GetNMatchedClusterIDsForTrack(fInputEvent,inTrack->GetID(),5,-5,0.2,-0.4) << "\t\t";
1475  }
1476  cout << endl;
1477  multimap<Int_t,Int_t>::iterator it;
1478  for (it=fMapTrackToCluster.begin(); it!=fMapTrackToCluster.end(); ++it) cout << it->first << " => " << it->second << '\n';
1479  cout << "mapClusterToTrack" << endl;
1480  Int_t tempClus = it->second;
1481  for (it=fMapClusterToTrack.begin(); it!=fMapClusterToTrack.end(); ++it) cout << it->first << " => " << it->second << '\n';
1482  vector<Int_t> tempTracks = GetMatchedTrackIDsForCluster(fInputEvent,tempClus, 5, -5, 0.2, -0.4);
1483  for(UInt_t iJ=0; iJ<tempTracks.size();iJ++){
1484  cout << tempClus << " - " << tempTracks.at(iJ) << endl;
1485  }
1486  }
1487  return;
1488 }
double Double_t
Definition: External.C:58
Definition: External.C:236
Int_t GetNMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
void ProcessEvent(AliVEvent *event)
Bool_t GetTrackClusterMatchingResidual(Int_t trackID, Int_t clusterID, Float_t &dEta, Float_t &dPhi)
void Initialize(Int_t runNumber)
vector< Int_t > GetMatchedClusterIDsForSecTrack(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
Int_t GetNMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
int Int_t
Definition: External.C:63
unsigned int UInt_t
Definition: External.C:33
float Float_t
Definition: External.C:68
virtual void Terminate(Option_t *)
Bool_t GetSecTrackClusterMatchingResidual(Int_t trackID, Int_t clusterID, Float_t &dEta, Float_t &dPhi)
Float_t SumTrackEtAroundCluster(AliVEvent *event, Int_t clusterID, Float_t dR)
Int_t GetNMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
vector< Int_t > GetMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
virtual void UserExec(Option_t *option)
Definition: External.C:220
Int_t GetNMatchedSecTrackIDsForCluster(AliVEvent *event, Int_t trackID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
Bool_t PropagateV0TrackToClusterAndGetMatchingResidual(AliVTrack *inSecTrack, AliVCluster *cluster, AliVEvent *event, Float_t &dEta, Float_t &dPhi)
pair< Float_t, Float_t > pairFloat
Bool_t IsSecTrackClusterAlreadyTried(Int_t trackID, Int_t clusterID)
const char Option_t
Definition: External.C:48
bool Bool_t
Definition: External.C:53
vector< Int_t > GetMatchedClusterIDsForTrack(AliVEvent *event, Int_t trackID, TF1 *fFuncPtDepEta, TF1 *fFuncPtDepPhi)
vector< Int_t > GetMatchedTrackIDsForCluster(AliVEvent *event, Int_t clusterID, Float_t dEtaMax, Float_t dEtaMin, Float_t dPhiMax, Float_t dPhiMin)
void SetLogBinningYTH2(TH2 *histoRebin)