]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/UserTasks/AliAnalysisTaskJetMatching.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / UserTasks / AliAnalysisTaskJetMatching.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 // 
17 // General task to match two sets of jets
18 //
19 // This task takes two TClonesArray's as input: 
20 // [1] fSourceJets - e.g. pythia jets
21 // [2] fTargetJets - e.g. a samle containing pythia jets embedded in a pbpb event
22 // the task will try to match jets from the source array to the target array, and
23 // save the found TARGET jets in a new array 'fMatchedJets'
24 // the purpose of this task is twofold
25 // [1] matching for embedding studies
26 // [2] matching to create a detector response function
27 //
28 // matching is done in three steps
29 // [1] geometric matching, where jets are matched by R = sqrt(dphi*dphi-deta*deta) or directly via dphi and deta
30 // [2] optional injection / bijection check 
31 //     in this check, it is made sure that fSourceJets -> fMatchedJets is either an injective non-surjective 
32 //     or bijective map, depending on the matching resolution which is chosen for step [1]
33 //     so that each source jet has a unique target and vice-versa.
34 //     if R (or dphi, deta) are proportional to, or larger than, the jet radius, matching tends to be 
35 //     bijective (each source has a target), if R is chosen to be very small, source jets might be lost 
36 //     as no target can be found.
37 //     the mapping is done in such a way that each target is matched to its closest source and each source
38 //     is mapped to its closest target
39 // [3] constituent matching
40 //     - how many constituents of the source jet are present in the matched jet? 
41 //       a cut on the constituent fraction can be performed (not recommended)
42 //     - how much of the original pt is recovered in the matched jet? 
43 //       a cut on the fraction of recovered / original pt can be performed (recommended)
44 //
45 // detector response
46 //     to obtain a detector respose function, supply
47 // [1] fSourceJets - particle level jets
48 // [2] fTargetJets - detector level jets
49 // 
50 // Author: Redmer Alexander Bertens, Utrecht Univeristy, Utrecht, Netherlands
51 // rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl 
52
53 // root includes
54 #include <TClonesArray.h>
55 #include <TChain.h>
56 #include <TMath.h>
57 #include <TF1.h>
58 #include <TF2.h>
59 #include <TH1F.h>
60 #include <TH2F.h>
61 #include <TProfile.h>
62 #include <TFile.h>
63 #include <TTree.h>
64 #include <TKey.h>
65 #include <TSystem.h>
66 // aliroot includes
67 #include <AliAnalysisTask.h>
68 #include <AliAnalysisManager.h>
69 #include <AliLog.h>
70 #include <AliVEvent.h>
71 #include <AliVParticle.h>
72 // emcal jet framework includes
73 #include <AliEmcalJet.h>
74 #include <AliAnalysisTaskJetMatching.h>
75 #include <AliLocalRhoParameter.h>
76
77 class AliAnalysisTaskJetMatching;
78 using namespace std;
79
80 ClassImp(AliAnalysisTaskJetMatching)
81
82 AliAnalysisTaskJetMatching::AliAnalysisTaskJetMatching() : AliAnalysisTaskEmcalJet("AliAnalysisTaskJetMatching", kTRUE), 
83     fDebug(0), fSourceJets(0), fSourceJetsName(0), fTargetJets(0), fTargetJetsName(0), fMatchedJets(0), fMatchedJetsName(GetName()), fSourceRho(0), fSourceRhoName(0), fTargetRho(0), fTargetRhoName(0), fUseScaledRho(0), fSourceRadius(0.3), fTargetRadius(0.3), fMatchingScheme(kGeoEtaPhi), fUseEmcalBaseJetCuts(kFALSE), fSourceBKG(kNoSourceBKG), fTargetBKG(kNoTargetBKG), fOutputList(0), fHistUnsortedCorrelation(0), fHistMatchedCorrelation(0), fHistSourceJetPt(0), fHistMatchedSourceJetPt(0), fHistTargetJetPt(0), fHistMatchedJetPt(0), fHistSourceMatchedJetPt(0), fHistDetectorResponseProb(0), fHistNoConstSourceJet(0), fHistNoConstTargetJet(0), fHistNoConstMatchJet(0), fProfFracPtMatched(0), fProfFracPtJets(0), fProfFracNoMatched(0), fProfFracNoJets(0), fHistDiJet(0), fHistDiJetLeadingJet(0), fHistDiJetDPhi(0), fHistAnalysisSummary(0), fProfQAMatched(0), fProfQA(0), fNoMatchedJets(200), fMatchEta(.3), fMatchPhi(.3), fMatchR(.08), fDoDetectorResponse(kFALSE), fMatchConstituents(kTRUE), fMinFracRecoveredConstituents(.5), fMinFracRecoveredConstituentPt(0.5), fGetBijection(kTRUE), fh1Trials(0x0), fh1AvgTrials(0x0), fh1Xsec(0x0), fAvgTrials(0) {
84     // default constructor
85     ClearMatchedJetsCache();
86 }
87 //_____________________________________________________________________________
88 AliAnalysisTaskJetMatching::AliAnalysisTaskJetMatching(const char* name) : AliAnalysisTaskEmcalJet(name, kTRUE),
89     fDebug(0), fSourceJets(0), fSourceJetsName(0), fTargetJets(0), fTargetJetsName(0), fMatchedJets(0), fMatchedJetsName(GetName()), fSourceRho(0), fSourceRhoName(0), fTargetRho(0), fTargetRhoName(0), fUseScaledRho(0), fSourceRadius(0.3), fTargetRadius(0.3), fMatchingScheme(kGeoEtaPhi), fUseEmcalBaseJetCuts(kFALSE), fSourceBKG(kNoSourceBKG), fTargetBKG(kNoTargetBKG), fOutputList(0), fHistUnsortedCorrelation(0), fHistMatchedCorrelation(0), fHistSourceJetPt(0), fHistMatchedSourceJetPt(0), fHistTargetJetPt(0), fHistMatchedJetPt(0), fHistSourceMatchedJetPt(0), fHistDetectorResponseProb(0), fHistNoConstSourceJet(0), fHistNoConstTargetJet(0), fHistNoConstMatchJet(0), fProfFracPtMatched(0), fProfFracPtJets(0), fProfFracNoMatched(0), fProfFracNoJets(0), fHistDiJet(0), fHistDiJetLeadingJet(0), fHistDiJetDPhi(0), fHistAnalysisSummary(0), fProfQAMatched(0), fProfQA(0), fNoMatchedJets(200), fMatchEta(.3), fMatchPhi(.3), fMatchR(.08), fDoDetectorResponse(kFALSE), fMatchConstituents(kTRUE), fMinFracRecoveredConstituents(0.5), fMinFracRecoveredConstituentPt(0.5), fGetBijection(kTRUE), fh1Trials(0x0), fh1AvgTrials(0x0), fh1Xsec(0x0), fAvgTrials(0) {
90     // constructor
91     ClearMatchedJetsCache();
92     DefineInput(0, TChain::Class());
93     DefineOutput(1, TList::Class());
94 }
95 //_____________________________________________________________________________
96 AliAnalysisTaskJetMatching::~AliAnalysisTaskJetMatching()
97 {
98     // destructor
99     if(fOutputList)             delete fOutputList;
100 }
101 //_____________________________________________________________________________
102 void AliAnalysisTaskJetMatching::ExecOnce() 
103 {
104     // initialize the anaysis
105     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
106     // get the stand alone jets from the input event
107     fSourceJets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fSourceJetsName.Data()));
108     if(!fSourceJets) AliFatal(Form("%s: Container with name %s not found. Aborting", GetName(), fSourceJetsName.Data()));
109     fTargetJets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTargetJetsName.Data()));
110     if(!fTargetJets) AliFatal(Form("%s: Container with name %s not found. Aborting", GetName(), fSourceJetsName.Data()));
111     // append the list of matched jets to the event
112     fMatchedJets->Delete();
113     if (!(InputEvent()->FindListObject(fMatchedJetsName))) InputEvent()->AddObject(fMatchedJets);
114     else AliFatal(Form("%s: Object with name %s already in event! Aborting", GetName(), fMatchedJetsName.Data()));
115     FillAnalysisSummaryHistogram();
116     switch (fSourceBKG) {
117         case kSourceLocalRho : {
118             fSourceRho = dynamic_cast<AliLocalRhoParameter*>(InputEvent()->FindListObject(fSourceRhoName));
119             if(!fSourceRho) AliFatal(Form("%s: Object with name %s requested but not found! Aborting", GetName(), fSourceRhoName.Data()));
120         } break;
121         default : break;
122     }
123     switch (fTargetBKG) {
124         case kTargetLocalRho : {
125             fTargetRho = dynamic_cast<AliLocalRhoParameter*>(InputEvent()->FindListObject(fTargetRhoName));
126             if(!fTargetRho) AliFatal(Form("%s: Object with name %s requested but not found! Aborting", GetName(), fTargetRhoName.Data()));
127         } break;
128         default : break;
129     }
130     AliAnalysisTaskEmcalJet::ExecOnce(); // init base class
131     if(fDoDetectorResponse) fMatchConstituents = kFALSE;
132 }
133 //_____________________________________________________________________________
134 void AliAnalysisTaskJetMatching::UserCreateOutputObjects()
135 {
136     // create output objects
137     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
138      // add the matched jets array to the event
139     fMatchedJets = new TClonesArray("AliEmcalJet");
140     fMatchedJets->SetName(fMatchedJetsName);
141     fOutputList = new TList();
142     fOutputList->SetOwner(kTRUE);
143     // add analysis histograms
144     fHistUnsortedCorrelation = BookTH1F("fHistUnsortedCorrelation", "#Delta #varphi unsorted", 50, 0, TMath::Pi());
145     fHistMatchedCorrelation = BookTH1F("fHistMatchedCorrelation", "#Delta #varphi matched", 50, 0, TMath::Pi());
146     fHistSourceJetPt = (fDoDetectorResponse) ? BookTH1F("fHistParticleLevelJetPt", "p_{t}^{gen} [GeV/c]", 150, -50, 250) :  BookTH1F("fHistSourceJetPt", "p_{t} [GeV/c]", 150, -50, 250);      
147     fHistMatchedSourceJetPt = (fDoDetectorResponse) ? BookTH1F("fHistMatchedParticleLevelJetPt", "p_{t}^{gen} [GeV/c]", 150, -50, 250) : BookTH1F("fHistMatchedSourceJetPt", "p_{t} [GeV/c]", 150, -50, 250);
148     fHistTargetJetPt = BookTH1F("fHistTargetJetPt", "p_{t} [GeV/c]", 150, -50, 250);
149     fHistMatchedJetPt = (fDoDetectorResponse) ? BookTH1F("fHistDetectorLevelJet", "p_{t}^{rec}", 150, -50, 250) : BookTH1F("fHistMatchedJetPt", "p_{t} [GeV/c]", 150, -50, 250);
150     fHistSourceMatchedJetPt = (fDoDetectorResponse) ? BookTH2F("fHistDetectorResponse", "particle level jet p_{t}^{gen} [GeV/c]", "detector level jet p_{t}^{rec} [GeV/c]", 300, -50, 250, 300, -50, 250) : BookTH2F("fHistSourceMatchedJetPt", "source jet p_{t} [GeV/c]", "matched jet p_{t} [GeV/c]", 300, -50, 250, 300, -50, 250);
151     if(fDoDetectorResponse) {
152         fHistDetectorResponseProb = BookTH2F("fHistDetectorResponseProb", "(p_{t}^{det} - p_{t}^{part})/p_{t}^{part}", "p_{t}^{part}", 100, -1.5, 1., 20, 0, 200);
153         fh1Xsec = new TProfile("fh1Xsec","xsec from pyxsec.root",1,0,1);
154         fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
155         fOutputList->Add(fh1Xsec);
156         fh1Trials = new TH1F("fh1Trials","trials root file",1,0,1);
157         fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
158         fOutputList->Add(fh1Trials);
159         fh1AvgTrials = new TH1F("fh1AvgTrials","avg trials root file",1,0,1);
160         fh1AvgTrials->GetXaxis()->SetBinLabel(1,"#sum{avg ntrials}");
161         fOutputList->Add(fh1AvgTrials);
162     }
163     fHistNoConstSourceJet = BookTH1F("fHistNoConstSourceJet", "number of constituents source jets", 50, 0, 50);
164     fHistNoConstTargetJet = BookTH1F("fHistNoConstTargetJet", "number of constituents target jets", 50, 0, 50);
165     fHistNoConstMatchJet = BookTH1F("fHistNoConstMatchJet", "number of constituents matched jets", 50, 0, 50);
166     if(!fDoDetectorResponse) { // these observables cannot be measured in current detector response mode
167         fProfFracPtMatched = new TProfile("fProfFracPtMatched", "recovered target p_{T} / source p_{T}", 15, -50, 250);
168         fOutputList->Add(fProfFracPtMatched);
169         fProfFracNoMatched = new TProfile("fProfFracNoMatched", "recovered target constituents / source constituents", 15, -50, 250);
170         fOutputList->Add(fProfFracNoMatched);
171     }
172     // the analysis summary histo which stores all the analysis flags is always written to file
173     fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 51, -0.5, 15.5);
174     fProfQAMatched = new TProfile("fProfQAMatched", "fProfQAMatched", 3, 0, 3);
175     fProfQAMatched->GetXaxis()->SetBinLabel(1, "<#delta p_{t} >");
176     fProfQAMatched->GetXaxis()->SetBinLabel(2, "<#delta #eta>");
177     fProfQAMatched->GetXaxis()->SetBinLabel(3, "<#delta #varphi>");
178     fOutputList->Add(fProfQAMatched);
179     fProfQA = new TProfile("fProfQA", "fProfQA", 3, 0, 3);
180     fProfQA->GetXaxis()->SetBinLabel(1, "<#delta p_{t} >");
181     fProfQA->GetXaxis()->SetBinLabel(2, "<#delta #eta>");
182     fProfQA->GetXaxis()->SetBinLabel(3, "<#delta #varphi>");
183     fOutputList->Add(fProfQA);
184     fProfFracPtJets = new TProfile("fProfFracPtJets", "source p_{T} / target p_{T}", 15, -50, 250);
185     fOutputList->Add(fProfFracPtJets);
186     fProfFracNoJets = new TProfile("fProfFracNoJets", "source constituents / target constituents", 15, -50, 250);
187     fOutputList->Add(fProfFracNoJets);
188     switch (fMatchingScheme) {
189         case kDiJet : {
190             fHistDiJet = BookTH2F("fHistDiJet", "matched di-jet #varphi", "matched di-jet #eta", 100, 0., TMath::TwoPi(), 100, -5., 5.);
191             fHistDiJetLeadingJet = BookTH2F("fHistDiJetLeadingJet", "leading jet #varphi", "leadingd jet #eta", 100, 0., TMath::TwoPi(), 100, -5., 5.);
192             fHistDiJetDPhi = BookTH1F("fHistDiJetDPhi", "leading jet #varphi - (matched jet #varphi - #pi)", 100, -1.*TMath::Pi(), TMath::Pi());
193         } break;
194         default : break;
195     }
196     fOutputList->Sort();
197     PostData(1, fOutputList);
198 }
199 //_____________________________________________________________________________
200 TH1F* AliAnalysisTaskJetMatching::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Bool_t append)
201 {
202     // book a TH1F and connect it to the output container
203     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
204     if(!fOutputList) return 0x0;
205     TString title(name);
206     title += Form(";%s;[counts]", x);
207     TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
208     histogram->Sumw2();
209     if(append) fOutputList->Add(histogram);
210     return histogram;   
211 }
212 //_____________________________________________________________________________
213 TH2F* AliAnalysisTaskJetMatching::BookTH2F(const char* name, const char* x, const char*y, Int_t binsx, Double_t minx, Double_t maxx, Int_t binsy, Double_t miny, Double_t maxy, Bool_t append)
214 {
215     // book a TH2F and connect it to the output container
216     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
217     if(!fOutputList) return 0x0;
218     TString title(name);
219     title += Form(";%s;%s", x, y);
220     TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
221     histogram->Sumw2();
222     if(append) fOutputList->Add(histogram);
223     return histogram;   
224 }
225 //_____________________________________________________________________________
226 Bool_t AliAnalysisTaskJetMatching::Notify()
227 {
228     // for each file get the number of trials and pythia cross section
229     // see  $ALICE_ROOT/PWGJE/AliAnalysisTaskJetProperties.cxx
230     //      $ALICE_ROOT/PWG/Tools/AliAnalysisHelperJetTasks.cxx
231     // this function is implenented here temporarily to avoid introducing a dependency
232     // later on this could just be a call to a static helper function
233     if(!fDoDetectorResponse) return kTRUE;
234     Float_t xsection(0), ftrials(1);
235     fAvgTrials = ftrials;
236     TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
237     if(tree) {
238         TFile *curfile = tree->GetCurrentFile();
239         if (!curfile) return kFALSE;
240         TString file(curfile->GetName()); 
241         if(file.Contains("root_archive.zip#")) file.Replace(file.Index("#",1,file.Index("root_archive",12,0,TString::kExact),TString::kExact)+1,file.Index(".root",5,TString::kExact)-file.Index("root_archive",12,0,TString::kExact),"");
242         else file.ReplaceAll(gSystem->BaseName(file.Data()),"");
243         TFile *fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec.root"));
244         if(!fxsec) {
245             fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
246             if(fxsec) {
247                 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0); 
248                 if(key) {
249                       TList *list = dynamic_cast<TList*>(key->ReadObj());
250                       if(list) {
251                           xsection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
252                           ftrials  = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
253                       }
254                 }
255                 fxsec->Close();
256             }
257         } else {
258             TTree *xtree = (TTree*)fxsec->Get("Xsection");
259             if(xtree) {
260                 UInt_t _ftrials  = 0;
261                 Double_t _xsection  = 0;
262                 xtree->SetBranchAddress("xsection",&_xsection);
263                 xtree->SetBranchAddress("ntrials",&_ftrials);
264                 xtree->GetEntry(0);
265                 ftrials = _ftrials;
266                 xsection = _xsection;
267             }
268             fxsec->Close();
269         }
270         fh1Xsec->Fill("<#sigma>",xsection);
271         Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
272         if(ftrials >= nEntries && nEntries > 0.) fAvgTrials = ftrials/nEntries;
273         fh1Trials->Fill("#sum{ntrials}",ftrials); 
274     }  
275     return kTRUE;
276 }
277 //_____________________________________________________________________________
278 Bool_t AliAnalysisTaskJetMatching::Run()
279 {
280     // execute once for each event
281     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
282     if(!(InputEvent() && fSourceJets && fTargetJets && IsEventSelected())) return kFALSE;
283     if(fh1AvgTrials) fh1AvgTrials->Fill("#sum{avg ntrials}", fAvgTrials);
284     // step one: do geometric matching 
285     switch (fMatchingScheme) {
286         // FIXME having separate dedicated functions is not necessary and historical
287         // cluttered code - should be cleaned up and merged into one
288         case kGeoEtaPhi : {
289             DoGeometricMatchingEtaPhi();
290         } break;
291         case kGeoR : {
292             DoGeometricMatchingR();
293         } break;
294         case kDiJet : {
295             DoDiJetMatching();
296         } break;
297         default : break;
298     }
299     // break if no jet was matched
300     if(!fMatchedJetContainer[1][0]) return kTRUE;
301     // optional step two: get a bijection (avoid duplicate matches)
302     if(fGetBijection)           GetBijection();
303     // optional step three: match constituents within matched jets
304     if(fMatchConstituents)      DoConstituentMatching();
305     // stream data to output
306     PostMatchedJets();
307     FillMatchedJetHistograms();
308     if(fDebug > 0) PrintInfo();
309     PostData(1, fOutputList);
310     return kTRUE;
311 }
312 //_____________________________________________________________________________
313 void AliAnalysisTaskJetMatching::DoGeometricMatchingEtaPhi()
314 {
315     // do geometric matching based on eta phi distance between jet centers
316     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
317     fNoMatchedJets = 0; // reset the matched jet counter
318     Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
319     for(Int_t i(0); i < iSource; i++) {
320         AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
321         if(!PassesCuts(sourceJet, 0)) continue;
322         for(Int_t j(0); j < iTarget; j++) {
323             AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
324             if(!PassesCuts(targetJet, 1)) continue;
325             if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
326             if((TMath::Abs(sourceJet->Eta() - targetJet->Eta()) < fMatchEta )) {
327                 Double_t sourcePhi(sourceJet->Phi()), targetPhi(targetJet->Phi());
328                 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
329                 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
330                 if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) {       // accept the jets as matching 
331                 Bool_t isBestMatch(kTRUE);
332                     if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
333                         if(fDebug > 0) printf(" > Entering first bijection test \n");
334                         for(Int_t k(i); k < iSource; k++) {
335                             AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
336                             if(PassesCuts(candidateSourceJet, 0)) continue;
337                             if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
338                             if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
339                                 isBestMatch = kFALSE;
340                                 break;
341                             }
342                         }
343                         if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
344                     }
345                     if(isBestMatch) {
346                         fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
347                         fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
348                         fNoMatchedJets++;
349                     }
350                     if(fNoMatchedJets > 199) {   // maximum to keep the cache at reasonable size
351                         AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
352                         return;
353                     }
354                 }
355             }
356         }
357     }
358 }
359 //_____________________________________________________________________________
360 void AliAnalysisTaskJetMatching::DoGeometricMatchingR()
361 {
362     // do geometric matching based on shortest path between jet centers
363     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
364     fNoMatchedJets = 0; // reset the matched jet counter
365     Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
366     for(Int_t i(0); i < iSource; i++) {
367         AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
368         if(!PassesCuts(sourceJet, 0)) continue;
369         for(Int_t j(0); j < iTarget; j++) {
370             AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
371             if(!PassesCuts(targetJet, 1)) continue;
372             else if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
373             if(GetR(sourceJet, targetJet) <= fMatchR) {
374                 Bool_t isBestMatch(kTRUE);
375                 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
376                     if(fDebug > 0) printf(" > Entering first bijection test \n");
377                     for(Int_t k(i); k < iSource; k++) {
378                         AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
379                         if(!PassesCuts(candidateSourceJet, 0)) continue;
380                         if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
381                         if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
382                             isBestMatch = kFALSE;
383                             break;
384                         }
385                     }
386                     if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
387                 }
388                 if(isBestMatch) {
389                     fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
390                     fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
391                     fNoMatchedJets++;
392                 }
393                 if(fNoMatchedJets > 99) {
394                     AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
395                     return;
396                 }
397             }
398         }
399     }
400 }
401 //_____________________________________________________________________________
402 void AliAnalysisTaskJetMatching::DoDiJetMatching()
403 {
404     // match dijets. this is in a sense a 'special' mode of the task as both jet supplied jet arrays 
405     // (target and source) are the same jet collection, matching will be performed on distribution in
406     // azimuth
407     // no ouptut array is produced, dedicated dijet histo's are filled in this function instead
408     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
409     // get the leading jet in a given acceptance (TPC is default)
410     Int_t leadingJetIndex(-1), subLeadingJetIndex(-1);
411     // retrieve the leading jet, leadingJetIndex points to the leading jet
412     AliEmcalJet* leadingJet(GetLeadingJet(fSourceJets, leadingJetIndex));
413     if(!leadingJet) return;
414     // fill phi and eta of leading jet (should always be in selected acceptance)
415     fHistDiJetLeadingJet->Fill(leadingJet->Phi(), leadingJet->Eta());
416     Double_t sourcePhi(leadingJet->Phi()), targetPhi(-1);
417     // get the sub-leading jet - faster when leading jet is also provided
418     AliEmcalJet* subLeadingJet(GetSubLeadingJet(fSourceJets, leadingJetIndex, subLeadingJetIndex));
419     if(!subLeadingJet) return;
420     else { // check if the sub-leading jet is actually the away-side jet
421         targetPhi = subLeadingJet->Phi() + TMath::Pi();
422         // rotate jets to common phase
423         if(TMath::Abs(sourcePhi) > TMath::Abs(sourcePhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
424         if(TMath::Abs(sourcePhi) > TMath::Abs(sourcePhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
425         if(TMath::Abs(targetPhi) > TMath::Abs(targetPhi+TMath::TwoPi())) targetPhi+=TMath::TwoPi();
426         if(TMath::Abs(targetPhi) > TMath::Abs(targetPhi-TMath::TwoPi())) targetPhi-=TMath::TwoPi();
427         if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) {
428             fHistDiJet->Fill(subLeadingJet->Phi(), subLeadingJet->Eta());
429             fHistDiJetDPhi->Fill(sourcePhi-targetPhi);
430         }
431     }
432 }
433 //_____________________________________________________________________________
434 void AliAnalysisTaskJetMatching::DoConstituentMatching()
435 {
436     // match constituents within matched jets 
437     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
438     if(!fTracks) {
439         AliFatal(Form("%s Fatal error! To do deep matching, supply jet constituents ! \n", GetName()));
440         return; // coverity ...
441     }
442     for(Int_t i(0); i < fNoMatchedJets; i++) {
443         AliEmcalJet* sourceJet = fMatchedJetContainer[i][0];
444         AliEmcalJet* targetJet = fMatchedJetContainer[i][1];
445         if(sourceJet && targetJet) {    // duplicate check: slot migth be NULL
446             Double_t targetPt(0);
447             Int_t iSJ(sourceJet->GetNumberOfTracks());
448             Int_t iTJ(targetJet->GetNumberOfTracks());
449             Int_t overlap(0), alreadyFound(0);
450             for(Int_t j(0); j < iSJ; j++) {
451                 alreadyFound = 0;
452                 Int_t idSource((Int_t)sourceJet->TrackAt(j));
453                 for(Int_t k(0); k < iTJ; k++) { // compare all target tracks to the source track
454                     if(idSource == targetJet->TrackAt(k) && alreadyFound == 0) {
455                         overlap++;
456                         alreadyFound++; // avoid possible duplicate matching
457                         AliVParticle* vp(static_cast<AliVParticle*>(targetJet->TrackAt(k, fTracks)));
458                         if(vp) targetPt += vp->Pt();
459                         continue;
460                     }
461                 }
462             }
463             if((float)overlap/(float)iSJ < fMinFracRecoveredConstituents || targetPt/sourceJet->Pt() < fMinFracRecoveredConstituentPt) {
464                 if(fDebug > 0) printf("  \n > Purging jet, recovered constituents ratio  %i / %i = %.2f <  or pt ratio %.2f / %.2f = %.2f < %.2f", 
465                         overlap, iSJ, (float)overlap/(float)iSJ, targetPt, sourceJet->Pt(), targetPt/sourceJet->Pt(), fMinFracRecoveredConstituentPt);
466                 fMatchedJetContainer[i][0] = 0x0;
467                 fMatchedJetContainer[i][1] = 0x0;
468                 continue;
469             }
470             if(sourceJet->Pt() > 0) {
471                 Double_t sourceRho(0), targetRho(0);
472                 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(sourceJet->Phi(), fSourceRadius)*sourceJet->Area();
473                 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(targetJet->Phi(), fTargetRadius)*targetJet->Area();
474                 fProfFracPtMatched->Fill(sourceJet->Pt()-sourceRho, (targetPt-targetRho) / (sourceJet->Pt()-sourceRho));
475                 fProfFracPtJets->Fill(sourceJet->Pt()-sourceRho, (targetJet->Pt()-targetRho) / (sourceJet->Pt()-sourceRho));
476                 fProfFracNoMatched->Fill(sourceJet->Pt()-sourceRho, (double)overlap / (double)sourceJet->GetNumberOfTracks());
477                 fProfFracNoJets->Fill(sourceJet->Pt()-sourceRho, (double)targetJet->GetNumberOfTracks() / (double)sourceJet->GetNumberOfTracks());
478             }
479             if(fDebug > 0) {
480                 printf("\n > Jet A: %i const\t", iSJ);
481                 printf(" > Jet B %i const\t", iTJ);
482                 printf(" -> OVERLAP: %i tracks <- \n", overlap);
483             }
484         }
485     }
486 }
487 //_____________________________________________________________________________
488 void AliAnalysisTaskJetMatching::GetBijection()
489 {
490     // bijection of source and matched jets, based on closest distance between jets
491     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
492     for(Int_t i(0); i < fNoMatchedJets; i++) {
493         for(Int_t j(i+1); j < fNoMatchedJets; j++) {
494             if(fMatchedJetContainer[i][0] == fMatchedJetContainer[j][0]) {
495                 // found source with two targets, now see which target is closer to the source
496                 if(!(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1] && fMatchedJetContainer[j][0] && fMatchedJetContainer[j][1] )) continue;
497                 Double_t rA(GetR(fMatchedJetContainer[i][0], fMatchedJetContainer[i][1]));      // distance between connection A = i
498                 Double_t rB(GetR(fMatchedJetContainer[j][0], fMatchedJetContainer[j][1]));      // distance between connection B = j
499                 if (rB > rA) {  // jet two is far away, clear it from both target and source list
500                     fMatchedJetContainer[j][0] = 0x0;
501                     fMatchedJetContainer[j][1] = 0x0;
502                 } else {                // jet one is far away, clear it from both target and source list
503                     fMatchedJetContainer[i][0] = fMatchedJetContainer[j][0];    // switch to avoid breaking loop
504                     fMatchedJetContainer[i][1] = fMatchedJetContainer[j][1];    
505                     fMatchedJetContainer[j][0] = 0x0;                           // clear duplicate jet from cache
506                     fMatchedJetContainer[j][1] = 0x0;
507                 }
508                 if(fDebug > 0) printf(" found duplicate jet, chose %.2f over %.2f \n" , (rB > rA) ? rA : rB, (rB > rA) ? rB : rA);
509             }
510         }
511     }
512 }
513 //_____________________________________________________________________________
514 void AliAnalysisTaskJetMatching::PostMatchedJets()
515 {
516     // post matched jets
517     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
518     fMatchedJets->Delete();     // should be a NULL operation, but added just in case
519     for(Int_t i(0), p(0); i < fNoMatchedJets; i++) {
520         if(fMatchedJetContainer[i][1]) {        // duplicate jets will have NULL value here and are skipped
521             new((*fMatchedJets)[p]) AliEmcalJet(*fMatchedJetContainer[i][1]);
522             p++; 
523         }
524     }
525 }
526 //_____________________________________________________________________________
527 void AliAnalysisTaskJetMatching::FillAnalysisSummaryHistogram() const
528 {
529     // fill the analysis summary histrogram, saves all relevant analysis settigns
530     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
531     fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fUseScaledRho");
532     fHistAnalysisSummary->SetBinContent(1, (int)fUseScaledRho);
533     fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fMatchingScheme");
534     fHistAnalysisSummary->SetBinContent(2, (int)fMatchingScheme);
535     fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fSourceBKG");
536     fHistAnalysisSummary->SetBinContent(3, (int)fSourceBKG);
537     fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fTargetBKG");
538     fHistAnalysisSummary->SetBinContent(4, (int)fTargetBKG);
539     fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fMatchEta");
540     fHistAnalysisSummary->SetBinContent(5, fMatchEta);
541     fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fMatchPhi");
542     fHistAnalysisSummary->SetBinContent(6, fMatchPhi);
543     fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fMatchR");
544     fHistAnalysisSummary->SetBinContent(7, fMatchR);
545     fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "iter");
546     fHistAnalysisSummary->SetBinContent(8, 1.);
547 }
548 //_____________________________________________________________________________
549 void AliAnalysisTaskJetMatching::FillMatchedJetHistograms() 
550 {
551     // fill matched jet histograms
552     if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
553     for(Int_t i(0); i < fSourceJets->GetEntriesFast(); i++) {
554         AliEmcalJet* source = static_cast<AliEmcalJet*>(fSourceJets->At(i));
555         if(!source) continue;
556         else if(fUseEmcalBaseJetCuts && !AcceptJet(source, 0)) continue;
557         Double_t sourceRho(0), targetRho(0);
558         if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(source->Phi(), fSourceRadius)*source->Area();
559         fHistSourceJetPt->Fill(source->Pt()-sourceRho);
560         fHistNoConstSourceJet->Fill(source->GetNumberOfConstituents());
561         for(Int_t j(0); j < fTargetJets->GetEntriesFast(); j++) {
562             AliEmcalJet* target = static_cast<AliEmcalJet*>(fTargetJets->At(j));
563             if(target) {
564                 if(fUseEmcalBaseJetCuts && !AcceptJet(target, 1)) continue;
565                 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(target->Phi(), fTargetRadius)*target->Area();
566                 fProfQA->Fill(0.5, TMath::Abs((source->Pt()-sourceRho)-(target->Pt()-targetRho)));  
567                 fProfQA->Fill(1.5, TMath::Abs(source->Eta()-target->Eta()));
568                 fProfQA->Fill(2.5, TMath::Abs(source->Phi()-target->Phi()));
569                 fHistUnsortedCorrelation->Fill(PhaseShift(source->Phi()-target->Phi(), 2));
570                 if(j==0) {
571                     fHistTargetJetPt->Fill(target->Pt()-targetRho);
572                     fHistNoConstTargetJet->Fill(target->GetNumberOfConstituents());
573                 }
574             }
575         }
576     }
577     for(Int_t i(0); i < fNoMatchedJets; i++) {
578         if(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1]) {
579             Double_t sourceRho(0), targetRho(0);
580             if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(fMatchedJetContainer[i][0]->Phi(), fSourceRadius)*fMatchedJetContainer[i][0]->Area();
581             if(fTargetRho) targetRho = fTargetRho->GetLocalVal(fMatchedJetContainer[i][1]->Phi(), fTargetRadius)*fMatchedJetContainer[i][1]->Area();
582             fHistMatchedCorrelation->Fill(PhaseShift(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi(), 2));
583             fHistMatchedSourceJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho);
584             fHistMatchedJetPt->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
585             fHistNoConstMatchJet->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
586             fProfQAMatched->Fill(0.5, TMath::Abs((fMatchedJetContainer[i][0]->Pt()-sourceRho)-(fMatchedJetContainer[i][1]->Pt()-targetRho)));
587             fProfQAMatched->Fill(1.5, TMath::Abs(fMatchedJetContainer[i][0]->Eta()-fMatchedJetContainer[i][1]->Eta()));
588             fProfQAMatched->Fill(2.5, TMath::Abs(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi()));
589             
590             fHistSourceMatchedJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, fMatchedJetContainer[i][1]->Pt()-targetRho);
591             if(fDoDetectorResponse) {
592                 fProfFracPtJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (fMatchedJetContainer[i][1]->Pt()-targetRho) / (fMatchedJetContainer[i][0]->Pt()-sourceRho));
593                 fProfFracNoJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (double)fMatchedJetContainer[i][1]->GetNumberOfTracks() / (double)fMatchedJetContainer[i][0]->GetNumberOfTracks());
594                 fHistDetectorResponseProb->Fill((fMatchedJetContainer[i][1]->Pt()-fMatchedJetContainer[i][0]->Pt())/fMatchedJetContainer[i][0]->Pt(), fMatchedJetContainer[i][0]->Pt());
595             }
596         }
597     }
598 }
599 //_____________________________________________________________________________
600 AliEmcalJet* AliAnalysisTaskJetMatching::GetLeadingJet(TClonesArray* source, Int_t &leadingJetIndex, Double_t etaMin, Double_t etaMax) 
601 {
602     // return the leading jet within given acceptance
603     Int_t iJets(source->GetEntriesFast());
604     Double_t pt(0);
605     AliEmcalJet* leadingJet(0x0);
606     for(Int_t i(0); i < iJets; i++) {
607         AliEmcalJet* jet = static_cast<AliEmcalJet*>(source->At(i));
608         if(jet->Eta() < etaMin || jet->Eta() > etaMax) continue;
609         if(jet->Pt() > pt) {
610            leadingJet = jet;
611            pt = leadingJet->Pt();
612            leadingJetIndex = i;
613         }
614     }
615     return leadingJet;
616 }
617 //_____________________________________________________________________________
618 AliEmcalJet* AliAnalysisTaskJetMatching::GetSubLeadingJet(TClonesArray* source, Int_t leadingJetIndex, Int_t &subLeadingJetIndex, Double_t etaMin, Double_t etaMax) 
619 {
620     // return the leading jet within given acceptance
621     // same as GetLeadingJet() but skips the leading jet (so returned jet is 
622     // sub-leading by design)
623     Int_t iJets(source->GetEntriesFast());
624     // if the leading jet isn't given, retrieve it
625     if(leadingJetIndex < 0) GetLeadingJet(source, leadingJetIndex, etaMin, etaMax);
626     Double_t pt(0);
627     AliEmcalJet* leadingJet(0x0);
628     for(Int_t i(0); i < iJets; i++) {
629         if(i == leadingJetIndex) continue;
630         AliEmcalJet* jet = static_cast<AliEmcalJet*>(source->At(i));
631         if(jet->Eta() < etaMin || jet->Eta() > etaMax) continue;
632         if(jet->Pt() > pt) {
633            leadingJet = jet;
634            pt = leadingJet->Pt();
635            subLeadingJetIndex = i;
636         }
637     }
638     return leadingJet;
639 }
640 //_____________________________________________________________________________
641 void AliAnalysisTaskJetMatching::PrintInfo() const
642 {
643     // print some info 
644     printf("\n > No. of source jets from %s \n \t %i \n ", fSourceJetsName.Data(), fSourceJets->GetEntriesFast());
645     printf(" > No. of target jets from %s \n \t %i \n ", fTargetJetsName.Data(), fTargetJets->GetEntriesFast());
646     printf(" > No. of matched jets from %s \n \t %i \n ", fMatchedJetsName.Data(), fMatchedJets->GetEntriesFast());
647     if(fDebug > 3) InputEvent()->GetList()->ls();
648 }
649 //_____________________________________________________________________________
650 void AliAnalysisTaskJetMatching::Terminate(Option_t *)
651 {
652     // terminate
653 }
654 //_____________________________________________________________________________