1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 // General task to match two sets of jets
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
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)
46 // to obtain a detector respose function, supply
47 // [1] fSourceJets - particle level jets
48 // [2] fTargetJets - detector level jets
50 // Author: Redmer Alexander Bertens, Utrecht Univeristy, Utrecht, Netherlands
51 // rbertens@cern.ch, rbertens@nikhef.nl, r.a.bertens@uu.nl
54 #include <TClonesArray.h>
63 #include <AliAnalysisTask.h>
64 #include <AliAnalysisManager.h>
66 #include <AliVEvent.h>
67 #include <AliVParticle.h>
68 // emcal jet framework includes
69 #include <AliEmcalJet.h>
70 #include <AliAnalysisTaskJetMatching.h>
71 #include <AliLocalRhoParameter.h>
73 class AliAnalysisTaskJetMatching;
76 ClassImp(AliAnalysisTaskJetMatching)
78 AliAnalysisTaskJetMatching::AliAnalysisTaskJetMatching() : AliAnalysisTaskEmcalJet("AliAnalysisTaskJetMatching", kTRUE),
79 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), fHistNoConstSourceJet(0), fHistNoConstTargetJet(0), fHistNoConstMatchJet(0), fProfFracPtMatched(0), fProfFracPtJets(0), fProfFracNoMatched(0), fProfFracNoJets(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) {
80 // default constructor
81 ClearMatchedJetsCache();
83 //_____________________________________________________________________________
84 AliAnalysisTaskJetMatching::AliAnalysisTaskJetMatching(const char* name) : AliAnalysisTaskEmcalJet(name, kTRUE),
85 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), fHistNoConstSourceJet(0), fHistNoConstTargetJet(0), fHistNoConstMatchJet(0), fProfFracPtMatched(0), fProfFracPtJets(0), fProfFracNoMatched(0), fProfFracNoJets(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) {
87 ClearMatchedJetsCache();
88 DefineInput(0, TChain::Class());
89 DefineOutput(1, TList::Class());
91 //_____________________________________________________________________________
92 AliAnalysisTaskJetMatching::~AliAnalysisTaskJetMatching()
95 if(fOutputList) delete fOutputList;
97 //_____________________________________________________________________________
98 void AliAnalysisTaskJetMatching::ExecOnce()
100 // initialize the anaysis
101 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
102 // get the stand alone jets from the input event
103 fSourceJets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fSourceJetsName.Data()));
104 if(!fSourceJets) AliFatal(Form("%s: Container with name %s not found. Aborting", GetName(), fSourceJetsName.Data()));
105 fTargetJets = dynamic_cast<TClonesArray*>(InputEvent()->FindListObject(fTargetJetsName.Data()));
106 if(!fTargetJets) AliFatal(Form("%s: Container with name %s not found. Aborting", GetName(), fSourceJetsName.Data()));
107 // append the list of matched jets to the event
108 fMatchedJets->Delete();
109 if (!(InputEvent()->FindListObject(fMatchedJetsName))) InputEvent()->AddObject(fMatchedJets);
110 else AliFatal(Form("%s: Object with name %s already in event! Aborting", GetName(), fMatchedJetsName.Data()));
111 FillAnalysisSummaryHistogram();
112 switch (fSourceBKG) {
113 case kSourceLocalRho : {
114 fSourceRho = dynamic_cast<AliLocalRhoParameter*>(InputEvent()->FindListObject(fSourceRhoName));
115 if(!fSourceRho) AliFatal(Form("%s: Object with name %s requested but not found! Aborting", GetName(), fSourceRhoName.Data()));
119 switch (fTargetBKG) {
120 case kTargetLocalRho : {
121 fTargetRho = dynamic_cast<AliLocalRhoParameter*>(InputEvent()->FindListObject(fTargetRhoName));
122 if(!fTargetRho) AliFatal(Form("%s: Object with name %s requested but not found! Aborting", GetName(), fTargetRhoName.Data()));
126 AliAnalysisTaskEmcalJet::ExecOnce(); // init base class
127 if(fDoDetectorResponse) fMatchConstituents = kFALSE;
129 //_____________________________________________________________________________
130 void AliAnalysisTaskJetMatching::UserCreateOutputObjects()
132 // create output objects
133 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
134 // add the matched jets array to the event
135 fMatchedJets = new TClonesArray("AliEmcalJet");
136 fMatchedJets->SetName(fMatchedJetsName);
137 fOutputList = new TList();
138 fOutputList->SetOwner(kTRUE);
139 // add analysis histograms
140 fHistUnsortedCorrelation = BookTH1F("fHistUnsortedCorrelation", "#Delta #varphi unsorted", 50, 0, TMath::Pi());
141 fHistMatchedCorrelation = BookTH1F("fHistMatchedCorrelation", "#Delta #varphi matched", 50, 0, TMath::Pi());
142 fHistSourceJetPt = (fDoDetectorResponse) ? BookTH1F("fHistParticleLevelJetPt", "p_{t}^{gen} [GeV/c]", 150, -50, 250) : BookTH1F("fHistSourceJetPt", "p_{t} [GeV/c]", 150, -50, 250);
143 fHistMatchedSourceJetPt = (fDoDetectorResponse) ? BookTH1F("fHistMatchedParticleLevelJetPt", "p_{t}^{gen} [GeV/c]", 150, -50, 250) : BookTH1F("fHistMatchedSourceJetPt", "p_{t} [GeV/c]", 150, -50, 250);
144 fHistTargetJetPt = BookTH1F("fHistTargetJetPt", "p_{t} [GeV/c]", 150, -50, 250);
145 fHistMatchedJetPt = (fDoDetectorResponse) ? BookTH1F("fHistDetectorLevelJet", "p_{t}^{rec}", 150, -50, 250) : BookTH1F("fHistMatchedJetPt", "p_{t} [GeV/c]", 150, -50, 250);
146 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);
147 fHistNoConstSourceJet = BookTH1F("fHistNoConstSourceJet", "number of constituents source jets", 50, 0, 50);
148 fHistNoConstTargetJet = BookTH1F("fHistNoConstTargetJet", "number of constituents target jets", 50, 0, 50);
149 fHistNoConstMatchJet = BookTH1F("fHistNoConstMatchJet", "number of constituents matched jets", 50, 0, 50);
150 if(!fDoDetectorResponse) { // these observables cannot be measured in current detector response mode
151 fProfFracPtMatched = new TProfile("fProfFracPtMatched", "recovered target p_{T} / source p_{T}", 15, -50, 250);
152 fOutputList->Add(fProfFracPtMatched);
153 fProfFracNoMatched = new TProfile("fProfFracNoMatched", "recovered target constituents / source constituents", 15, -50, 250);
154 fOutputList->Add(fProfFracNoMatched);
156 // the analysis summary histo which stores all the analysis flags is always written to file
157 fHistAnalysisSummary = BookTH1F("fHistAnalysisSummary", "flag", 51, -0.5, 15.5);
158 fProfQAMatched = new TProfile("fProfQAMatched", "fProfQAMatched", 3, 0, 3);
159 fProfQAMatched->GetXaxis()->SetBinLabel(1, "<#delta p_{t} >");
160 fProfQAMatched->GetXaxis()->SetBinLabel(2, "<#delta #eta>");
161 fProfQAMatched->GetXaxis()->SetBinLabel(3, "<#delta #varphi>");
162 fOutputList->Add(fProfQAMatched);
163 fProfQA = new TProfile("fProfQA", "fProfQA", 3, 0, 3);
164 fProfQA->GetXaxis()->SetBinLabel(1, "<#delta p_{t} >");
165 fProfQA->GetXaxis()->SetBinLabel(2, "<#delta #eta>");
166 fProfQA->GetXaxis()->SetBinLabel(3, "<#delta #varphi>");
167 fOutputList->Add(fProfQA);
168 fProfFracPtJets = new TProfile("fProfFracPtJets", "source p_{T} / target p_{T}", 15, -50, 250);
169 fOutputList->Add(fProfFracPtJets);
170 fProfFracNoJets = new TProfile("fProfFracNoJets", "source constituents / target constituents", 15, -50, 250);
171 fOutputList->Add(fProfFracNoJets);
173 PostData(1, fOutputList);
175 //_____________________________________________________________________________
176 TH1F* AliAnalysisTaskJetMatching::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Bool_t append)
178 // book a TH1F and connect it to the output container
179 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
180 if(!fOutputList) return 0x0;
182 title += Form(";%s;[counts]", x);
183 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
185 if(append) fOutputList->Add(histogram);
188 //_____________________________________________________________________________
189 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)
191 // book a TH2F and connect it to the output container
192 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
193 if(!fOutputList) return 0x0;
195 title += Form(";%s;%s", x, y);
196 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
198 if(append) fOutputList->Add(histogram);
201 //_____________________________________________________________________________
202 Bool_t AliAnalysisTaskJetMatching::Run()
204 // execute once for each event
205 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
206 if(!(InputEvent() && fSourceJets && fTargetJets && IsEventSelected())) return kFALSE;
207 // step one: do geometric matching
208 switch (fMatchingScheme) {
210 DoGeometricMatchingEtaPhi();
213 DoGeometricMatchingR();
217 // break if no jet was matched
218 if(!fMatchedJetContainer[1][0]) return kTRUE;
219 // optional step two: get a bijection (avoid duplicate matches)
220 if(fGetBijection) GetBijection();
221 // optional step three: match constituents within matched jets
222 if(fMatchConstituents) DoConstituentMatching();
223 // stream data to output
225 FillMatchedJetHistograms();
226 if(fDebug > 0) PrintInfo();
227 PostData(1, fOutputList);
230 //_____________________________________________________________________________
231 void AliAnalysisTaskJetMatching::DoGeometricMatchingEtaPhi()
233 // do geometric matching based on eta phi distance between jet centers
234 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
235 fNoMatchedJets = 0; // reset the matched jet counter
236 Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
237 for(Int_t i(0); i < iSource; i++) {
238 AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
239 if(!PassesCuts(sourceJet)) continue;
240 if(fUseEmcalBaseJetCuts && !AcceptJet(sourceJet, 0)) continue;
241 for(Int_t j(0); j < iTarget; j++) {
242 AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
243 if(!PassesCuts(targetJet)) continue;
244 if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
245 if((TMath::Abs(sourceJet->Eta() - targetJet->Eta()) < fMatchEta )) {
246 Double_t sourcePhi(sourceJet->Phi()), targetPhi(targetJet->Phi());
247 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
248 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
249 if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) { // accept the jets as matching
250 Bool_t isBestMatch(kTRUE);
251 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
252 if(fDebug > 0) printf(" > Entering first bbijection test \n");
253 for(Int_t k(i); k < iSource; k++) {
254 AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
255 if(!PassesCuts(candidateSourceJet)) continue;
256 if(fUseEmcalBaseJetCuts && !AcceptJet(candidateSourceJet, 0)) continue;
257 if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
258 if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
259 isBestMatch = kFALSE;
263 if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
266 fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
267 fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
270 if(fNoMatchedJets > 199) { // maximum to keep the cache at reasonable size
271 AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
279 //_____________________________________________________________________________
280 void AliAnalysisTaskJetMatching::DoGeometricMatchingR()
282 // do geometric matching based on shortest path between jet centers
283 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
284 fNoMatchedJets = 0; // reset the matched jet counter
285 Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
286 for(Int_t i(0); i < iSource; i++) {
287 AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
288 if(!PassesCuts(sourceJet)) continue;
289 else if (fUseEmcalBaseJetCuts && !AcceptJet(sourceJet, 0)) continue;
290 for(Int_t j(0); j < iTarget; j++) {
291 AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
292 if(!PassesCuts(targetJet)) continue;
293 else if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
294 if(GetR(sourceJet, targetJet) <= fMatchR) {
295 Bool_t isBestMatch(kTRUE);
296 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
297 if(fDebug > 0) printf(" > Entering first bijection test \n");
298 for(Int_t k(i); k < iSource; k++) {
299 AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
300 if(!PassesCuts(candidateSourceJet)) continue;
301 if(fUseEmcalBaseJetCuts && !AcceptJet(candidateSourceJet, 0)) continue;
302 if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
303 if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
304 isBestMatch = kFALSE;
308 if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
311 fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
312 fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
315 if(fNoMatchedJets > 99) {
316 AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
323 //_____________________________________________________________________________
324 void AliAnalysisTaskJetMatching::DoConstituentMatching()
326 // match constituents within matched jets
327 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
329 AliFatal(Form("%s Fatal error! To do deep matching, supply jet constituents ! \n", GetName()));
330 return; // coverity ...
332 for(Int_t i(0); i < fNoMatchedJets; i++) {
333 AliEmcalJet* sourceJet = fMatchedJetContainer[i][0];
334 AliEmcalJet* targetJet = fMatchedJetContainer[i][1];
335 if(sourceJet && targetJet) { // duplicate check: slot migth be NULL
336 Double_t targetPt(0);
337 Int_t iSJ(sourceJet->GetNumberOfTracks());
338 Int_t iTJ(targetJet->GetNumberOfTracks());
339 Int_t overlap(0), alreadyFound(0);
340 for(Int_t j(0); j < iSJ; j++) {
342 Int_t idSource((Int_t)sourceJet->TrackAt(j));
343 for(Int_t k(0); k < iTJ; k++) { // compare all target tracks to the source track
344 if(idSource == targetJet->TrackAt(k) && alreadyFound == 0) {
346 alreadyFound++; // avoid possible duplicate matching
347 AliVParticle* vp(static_cast<AliVParticle*>(targetJet->TrackAt(k, fTracks)));
348 if(vp) targetPt += vp->Pt();
353 if((float)overlap/(float)iSJ < fMinFracRecoveredConstituents || targetPt/sourceJet->Pt() < fMinFracRecoveredConstituentPt) {
354 if(fDebug > 0) printf(" \n > Purging jet, recovered constituents ratio %i / %i = %.2f < or pt ratio %.2f / %.2f = %.2f < %.2f",
355 overlap, iSJ, (float)overlap/(float)iSJ, targetPt, sourceJet->Pt(), targetPt/sourceJet->Pt(), fMinFracRecoveredConstituentPt);
356 fMatchedJetContainer[i][0] = 0x0;
357 fMatchedJetContainer[i][1] = 0x0;
360 if(sourceJet->Pt() > 0) {
361 Double_t sourceRho(0), targetRho(0);
362 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(sourceJet->Phi(), fSourceRadius)*sourceJet->Area();
363 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(targetJet->Phi(), fTargetRadius)*targetJet->Area();
364 fProfFracPtMatched->Fill(sourceJet->Pt()-sourceRho, (targetPt-targetRho) / (sourceJet->Pt()-sourceRho));
365 fProfFracPtJets->Fill(sourceJet->Pt()-sourceRho, (targetJet->Pt()-targetRho) / (sourceJet->Pt()-sourceRho));
366 fProfFracNoMatched->Fill(sourceJet->Pt()-sourceRho, (double)overlap / (double)sourceJet->GetNumberOfTracks());
367 fProfFracNoJets->Fill(sourceJet->Pt()-sourceRho, (double)targetJet->GetNumberOfTracks() / (double)sourceJet->GetNumberOfTracks());
370 printf("\n > Jet A: %i const\t", iSJ);
371 printf(" > Jet B %i const\t", iTJ);
372 printf(" -> OVERLAP: %i tracks <- \n", overlap);
377 //_____________________________________________________________________________
378 void AliAnalysisTaskJetMatching::GetBijection()
380 // bijection of source and matched jets, based on closest distance between jets
381 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
382 for(Int_t i(0); i < fNoMatchedJets; i++) {
383 for(Int_t j(i+1); j < fNoMatchedJets; j++) {
384 if(fMatchedJetContainer[i][0] == fMatchedJetContainer[j][0]) {
385 // found source with two targets, now see which target is closer to the source
386 if(!(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1] && fMatchedJetContainer[j][0] && fMatchedJetContainer[j][1] )) continue;
387 Double_t rA(GetR(fMatchedJetContainer[i][0], fMatchedJetContainer[i][1])); // distance between connection A = i
388 Double_t rB(GetR(fMatchedJetContainer[j][0], fMatchedJetContainer[j][1])); // distance between connection B = j
389 if (rB > rA) { // jet two is far away, clear it from both target and source list
390 fMatchedJetContainer[j][0] = 0x0;
391 fMatchedJetContainer[j][1] = 0x0;
392 } else { // jet one is far away, clear it from both target and source list
393 fMatchedJetContainer[i][0] = fMatchedJetContainer[j][0]; // switch to avoid breaking loop
394 fMatchedJetContainer[i][1] = fMatchedJetContainer[j][1];
395 fMatchedJetContainer[j][0] = 0x0; // clear duplicate jet from cache
396 fMatchedJetContainer[j][1] = 0x0;
398 if(fDebug > 0) printf(" found duplicate jet, chose %.2f over %.2f \n" , (rB > rA) ? rA : rB, (rB > rA) ? rB : rA);
403 //_____________________________________________________________________________
404 void AliAnalysisTaskJetMatching::PostMatchedJets()
407 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
408 fMatchedJets->Delete(); // should be a NULL operation, but added just in case
409 for(Int_t i(0), p(0); i < fNoMatchedJets; i++) {
410 if(fMatchedJetContainer[i][1]) { // duplicate jets will have NULL value here and are skipped
411 new((*fMatchedJets)[p]) AliEmcalJet(*fMatchedJetContainer[i][1]);
416 //_____________________________________________________________________________
417 void AliAnalysisTaskJetMatching::FillAnalysisSummaryHistogram() const
419 // fill the analysis summary histrogram, saves all relevant analysis settigns
420 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
421 fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fUseScaledRho");
422 fHistAnalysisSummary->SetBinContent(1, (int)fUseScaledRho);
423 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fMatchingScheme");
424 fHistAnalysisSummary->SetBinContent(2, (int)fMatchingScheme);
425 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fSourceBKG");
426 fHistAnalysisSummary->SetBinContent(3, (int)fSourceBKG);
427 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fTargetBKG");
428 fHistAnalysisSummary->SetBinContent(4, (int)fTargetBKG);
429 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fMatchEta");
430 fHistAnalysisSummary->SetBinContent(5, fMatchEta);
431 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fMatchPhi");
432 fHistAnalysisSummary->SetBinContent(6, fMatchPhi);
433 fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fMatchR");
434 fHistAnalysisSummary->SetBinContent(7, fMatchR);
435 fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "iter");
436 fHistAnalysisSummary->SetBinContent(8, 1.);
438 //_____________________________________________________________________________
439 void AliAnalysisTaskJetMatching::FillMatchedJetHistograms()
441 // fill matched jet histograms
442 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
443 for(Int_t i(0); i < fSourceJets->GetEntriesFast(); i++) {
444 AliEmcalJet* source = static_cast<AliEmcalJet*>(fSourceJets->At(i));
445 if(!source) continue;
446 else if(fUseEmcalBaseJetCuts && !AcceptJet(source, 0)) continue;
447 Double_t sourceRho(0), targetRho(0);
448 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(source->Phi(), fSourceRadius)*source->Area();
449 fHistSourceJetPt->Fill(source->Pt()-sourceRho);
450 fHistNoConstSourceJet->Fill(source->GetNumberOfConstituents());
451 for(Int_t j(0); j < fTargetJets->GetEntriesFast(); j++) {
452 AliEmcalJet* target = static_cast<AliEmcalJet*>(fTargetJets->At(j));
454 if(fUseEmcalBaseJetCuts && !AcceptJet(target, 1)) continue;
455 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(target->Phi(), fTargetRadius)*target->Area();
456 fProfQA->Fill(0.5, TMath::Abs((source->Pt()-sourceRho)-(target->Pt()-targetRho)));
457 fProfQA->Fill(1.5, TMath::Abs(source->Eta()-target->Eta()));
458 fProfQA->Fill(2.5, TMath::Abs(source->Phi()-target->Phi()));
459 fHistUnsortedCorrelation->Fill(PhaseShift(source->Phi()-target->Phi(), 2));
461 fHistTargetJetPt->Fill(target->Pt()-targetRho);
462 fHistNoConstTargetJet->Fill(target->GetNumberOfConstituents());
467 for(Int_t i(0); i < fNoMatchedJets; i++) {
468 if(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1]) {
469 Double_t sourceRho(0), targetRho(0);
470 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(fMatchedJetContainer[i][0]->Phi(), fSourceRadius)*fMatchedJetContainer[i][0]->Area();
471 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(fMatchedJetContainer[i][1]->Phi(), fTargetRadius)*fMatchedJetContainer[i][1]->Area();
472 fHistMatchedCorrelation->Fill(PhaseShift(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi(), 2));
473 fHistMatchedSourceJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho);
474 fHistMatchedJetPt->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
475 fHistNoConstMatchJet->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
476 fProfQAMatched->Fill(0.5, TMath::Abs((fMatchedJetContainer[i][0]->Pt()-sourceRho)-(fMatchedJetContainer[i][1]->Pt()-targetRho)));
477 fProfQAMatched->Fill(1.5, TMath::Abs(fMatchedJetContainer[i][0]->Eta()-fMatchedJetContainer[i][1]->Eta()));
478 fProfQAMatched->Fill(2.5, TMath::Abs(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi()));
479 fHistSourceMatchedJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, fMatchedJetContainer[i][1]->Pt()-targetRho);
480 if(fDoDetectorResponse) {
481 fProfFracPtJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (fMatchedJetContainer[i][1]->Pt()-targetRho) / (fMatchedJetContainer[i][0]->Pt()-sourceRho));
482 fProfFracNoJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (double)fMatchedJetContainer[i][1]->GetNumberOfTracks() / (double)fMatchedJetContainer[i][0]->GetNumberOfTracks());
487 //_____________________________________________________________________________
488 void AliAnalysisTaskJetMatching::PrintInfo() const
491 printf("\n > No. of source jets from %s \n \t %i \n ", fSourceJetsName.Data(), fSourceJets->GetEntriesFast());
492 printf(" > No. of target jets from %s \n \t %i \n ", fTargetJetsName.Data(), fTargetJets->GetEntriesFast());
493 printf(" > No. of matched jets from %s \n \t %i \n ", fMatchedJetsName.Data(), fMatchedJets->GetEntriesFast());
494 if(fDebug > 3) InputEvent()->GetList()->ls();
496 //_____________________________________________________________________________
497 void AliAnalysisTaskJetMatching::Terminate(Option_t *)
501 //_____________________________________________________________________________