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>
67 #include <AliAnalysisTask.h>
68 #include <AliAnalysisManager.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>
77 class AliAnalysisTaskJetMatching;
80 ClassImp(AliAnalysisTaskJetMatching)
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), 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();
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), 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) {
91 ClearMatchedJetsCache();
92 DefineInput(0, TChain::Class());
93 DefineOutput(1, TList::Class());
95 //_____________________________________________________________________________
96 AliAnalysisTaskJetMatching::~AliAnalysisTaskJetMatching()
99 if(fOutputList) delete fOutputList;
101 //_____________________________________________________________________________
102 void AliAnalysisTaskJetMatching::ExecOnce()
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()));
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()));
130 AliAnalysisTaskEmcalJet::ExecOnce(); // init base class
131 if(fDoDetectorResponse) fMatchConstituents = kFALSE;
133 //_____________________________________________________________________________
134 void AliAnalysisTaskJetMatching::UserCreateOutputObjects()
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);
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);
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) {
190 fHistDiJet = BookTH2F("fHistDiJet", "matched jet #varphi", "matched jet #eta", 100, 0., TMath::TwoPi(), 100, -2.5, 2.5);
191 fHistDiJetLeadingJet = BookTH2F("fHistDiJetLeadingJet", "matched jet #varphi", "matched jet #eta", 100, 0., TMath::TwoPi(), 100, -2.5, 2.5);
196 PostData(1, fOutputList);
198 //_____________________________________________________________________________
199 TH1F* AliAnalysisTaskJetMatching::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Bool_t append)
201 // book a TH1F and connect it to the output container
202 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
203 if(!fOutputList) return 0x0;
205 title += Form(";%s;[counts]", x);
206 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
208 if(append) fOutputList->Add(histogram);
211 //_____________________________________________________________________________
212 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 // book a TH2F and connect it to the output container
215 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
216 if(!fOutputList) return 0x0;
218 title += Form(";%s;%s", x, y);
219 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
221 if(append) fOutputList->Add(histogram);
224 //_____________________________________________________________________________
225 Bool_t AliAnalysisTaskJetMatching::Notify()
227 // for each file get the number of trials and pythia cross section
228 // see $ALICE_ROOT/PWGJE/AliAnalysisTaskJetProperties.cxx
229 // $ALICE_ROOT/PWG/Tools/AliAnalysisHelperJetTasks.cxx
230 // this function is implenented here temporarily to avoid introducing a dependency
231 // later on this could just be a call to a static helper function
232 if(!fDoDetectorResponse) return kTRUE;
233 Float_t xsection(0), ftrials(1);
234 fAvgTrials = ftrials;
235 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
237 TFile *curfile = tree->GetCurrentFile();
238 if (!curfile) return kFALSE;
239 TString file(curfile->GetName());
240 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),"");
241 else file.ReplaceAll(gSystem->BaseName(file.Data()),"");
242 TFile *fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec.root"));
244 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
246 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
248 TList *list = dynamic_cast<TList*>(key->ReadObj());
250 xsection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
251 ftrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
257 TTree *xtree = (TTree*)fxsec->Get("Xsection");
260 Double_t _xsection = 0;
261 xtree->SetBranchAddress("xsection",&_xsection);
262 xtree->SetBranchAddress("ntrials",&_ftrials);
265 xsection = _xsection;
269 fh1Xsec->Fill("<#sigma>",xsection);
270 Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
271 if(ftrials >= nEntries && nEntries > 0.) fAvgTrials = ftrials/nEntries;
272 fh1Trials->Fill("#sum{ntrials}",ftrials);
276 //_____________________________________________________________________________
277 Bool_t AliAnalysisTaskJetMatching::Run()
279 // execute once for each event
280 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
281 if(!(InputEvent() && fSourceJets && fTargetJets && IsEventSelected())) return kFALSE;
282 if(fh1AvgTrials) fh1AvgTrials->Fill("#sum{avg ntrials}", fAvgTrials);
283 // step one: do geometric matching
284 switch (fMatchingScheme) {
285 // FIXME having separate dedicated functions is not necessary and historical
286 // cluttered code - should be cleaned up and merged into one
288 DoGeometricMatchingEtaPhi();
291 DoGeometricMatchingR();
298 // break if no jet was matched
299 if(!fMatchedJetContainer[1][0]) return kTRUE;
300 // optional step two: get a bijection (avoid duplicate matches)
301 if(fGetBijection) GetBijection();
302 // optional step three: match constituents within matched jets
303 if(fMatchConstituents) DoConstituentMatching();
304 // stream data to output
306 FillMatchedJetHistograms();
307 if(fDebug > 0) PrintInfo();
308 PostData(1, fOutputList);
311 //_____________________________________________________________________________
312 void AliAnalysisTaskJetMatching::DoGeometricMatchingEtaPhi()
314 // do geometric matching based on eta phi distance between jet centers
315 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
316 fNoMatchedJets = 0; // reset the matched jet counter
317 Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
318 for(Int_t i(0); i < iSource; i++) {
319 AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
320 if(!PassesCuts(sourceJet, 0)) continue;
321 for(Int_t j(0); j < iTarget; j++) {
322 AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
323 if(!PassesCuts(targetJet, 1)) continue;
324 if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
325 if((TMath::Abs(sourceJet->Eta() - targetJet->Eta()) < fMatchEta )) {
326 Double_t sourcePhi(sourceJet->Phi()), targetPhi(targetJet->Phi());
327 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
328 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
329 if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) { // accept the jets as matching
330 Bool_t isBestMatch(kTRUE);
331 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
332 if(fDebug > 0) printf(" > Entering first bbijection test \n");
333 for(Int_t k(i); k < iSource; k++) {
334 AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
335 if(PassesCuts(candidateSourceJet, 0)) continue;
336 if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
337 if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
338 isBestMatch = kFALSE;
342 if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
345 fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
346 fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
349 if(fNoMatchedJets > 199) { // maximum to keep the cache at reasonable size
350 AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
358 //_____________________________________________________________________________
359 void AliAnalysisTaskJetMatching::DoGeometricMatchingR()
361 // do geometric matching based on shortest path between jet centers
362 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
363 fNoMatchedJets = 0; // reset the matched jet counter
364 Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
365 for(Int_t i(0); i < iSource; i++) {
366 AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
367 if(!PassesCuts(sourceJet, 0)) continue;
368 for(Int_t j(0); j < iTarget; j++) {
369 AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
370 if(!PassesCuts(targetJet, 1)) continue;
371 else if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
372 if(GetR(sourceJet, targetJet) <= fMatchR) {
373 Bool_t isBestMatch(kTRUE);
374 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
375 if(fDebug > 0) printf(" > Entering first bijection test \n");
376 for(Int_t k(i); k < iSource; k++) {
377 AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
378 if(!PassesCuts(candidateSourceJet, 0)) continue;
379 if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
380 if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
381 isBestMatch = kFALSE;
385 if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
388 fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
389 fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
392 if(fNoMatchedJets > 99) {
393 AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
400 //_____________________________________________________________________________
401 void AliAnalysisTaskJetMatching::DoDiJetMatching()
403 // match dijets. this is in a sense a 'special' mode of the task as both jet supplied jet arrays
404 // (target and source) are the same jet collection, matching will be performed on distribution in
406 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
407 // get the leading jet in a given acceptance
408 Int_t leadingJetIndex(-1), subLeadingJetIndex(-1);
409 AliEmcalJet* leadingJet(GetLeadingJet(fSourceJets, leadingJetIndex));
410 if(!leadingJet) return;
411 fHistDiJetLeadingJet->Fill(leadingJet->Eta(), leadingJet->Phi());
412 Double_t sourcePhi(leadingJet->Phi()), targetPhi(-1);
413 // grab the sub-leading jet
414 AliEmcalJet* subLeadingJet(GetSubLeadingJet(fSourceJets, leadingJetIndex, subLeadingJetIndex));
415 if(!subLeadingJet) return;
416 else { // check if the sub-leading jet is actually the away-side jet
417 targetPhi = subLeadingJet->Phi() + TMath::Pi();
418 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
419 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
420 if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) fHistDiJet->Fill(subLeadingJet->Eta(), subLeadingJet->Phi());
423 //_____________________________________________________________________________
424 void AliAnalysisTaskJetMatching::DoConstituentMatching()
426 // match constituents within matched jets
427 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
429 AliFatal(Form("%s Fatal error! To do deep matching, supply jet constituents ! \n", GetName()));
430 return; // coverity ...
432 for(Int_t i(0); i < fNoMatchedJets; i++) {
433 AliEmcalJet* sourceJet = fMatchedJetContainer[i][0];
434 AliEmcalJet* targetJet = fMatchedJetContainer[i][1];
435 if(sourceJet && targetJet) { // duplicate check: slot migth be NULL
436 Double_t targetPt(0);
437 Int_t iSJ(sourceJet->GetNumberOfTracks());
438 Int_t iTJ(targetJet->GetNumberOfTracks());
439 Int_t overlap(0), alreadyFound(0);
440 for(Int_t j(0); j < iSJ; j++) {
442 Int_t idSource((Int_t)sourceJet->TrackAt(j));
443 for(Int_t k(0); k < iTJ; k++) { // compare all target tracks to the source track
444 if(idSource == targetJet->TrackAt(k) && alreadyFound == 0) {
446 alreadyFound++; // avoid possible duplicate matching
447 AliVParticle* vp(static_cast<AliVParticle*>(targetJet->TrackAt(k, fTracks)));
448 if(vp) targetPt += vp->Pt();
453 if((float)overlap/(float)iSJ < fMinFracRecoveredConstituents || targetPt/sourceJet->Pt() < fMinFracRecoveredConstituentPt) {
454 if(fDebug > 0) printf(" \n > Purging jet, recovered constituents ratio %i / %i = %.2f < or pt ratio %.2f / %.2f = %.2f < %.2f",
455 overlap, iSJ, (float)overlap/(float)iSJ, targetPt, sourceJet->Pt(), targetPt/sourceJet->Pt(), fMinFracRecoveredConstituentPt);
456 fMatchedJetContainer[i][0] = 0x0;
457 fMatchedJetContainer[i][1] = 0x0;
460 if(sourceJet->Pt() > 0) {
461 Double_t sourceRho(0), targetRho(0);
462 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(sourceJet->Phi(), fSourceRadius)*sourceJet->Area();
463 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(targetJet->Phi(), fTargetRadius)*targetJet->Area();
464 fProfFracPtMatched->Fill(sourceJet->Pt()-sourceRho, (targetPt-targetRho) / (sourceJet->Pt()-sourceRho));
465 fProfFracPtJets->Fill(sourceJet->Pt()-sourceRho, (targetJet->Pt()-targetRho) / (sourceJet->Pt()-sourceRho));
466 fProfFracNoMatched->Fill(sourceJet->Pt()-sourceRho, (double)overlap / (double)sourceJet->GetNumberOfTracks());
467 fProfFracNoJets->Fill(sourceJet->Pt()-sourceRho, (double)targetJet->GetNumberOfTracks() / (double)sourceJet->GetNumberOfTracks());
470 printf("\n > Jet A: %i const\t", iSJ);
471 printf(" > Jet B %i const\t", iTJ);
472 printf(" -> OVERLAP: %i tracks <- \n", overlap);
477 //_____________________________________________________________________________
478 void AliAnalysisTaskJetMatching::GetBijection()
480 // bijection of source and matched jets, based on closest distance between jets
481 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
482 for(Int_t i(0); i < fNoMatchedJets; i++) {
483 for(Int_t j(i+1); j < fNoMatchedJets; j++) {
484 if(fMatchedJetContainer[i][0] == fMatchedJetContainer[j][0]) {
485 // found source with two targets, now see which target is closer to the source
486 if(!(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1] && fMatchedJetContainer[j][0] && fMatchedJetContainer[j][1] )) continue;
487 Double_t rA(GetR(fMatchedJetContainer[i][0], fMatchedJetContainer[i][1])); // distance between connection A = i
488 Double_t rB(GetR(fMatchedJetContainer[j][0], fMatchedJetContainer[j][1])); // distance between connection B = j
489 if (rB > rA) { // jet two is far away, clear it from both target and source list
490 fMatchedJetContainer[j][0] = 0x0;
491 fMatchedJetContainer[j][1] = 0x0;
492 } else { // jet one is far away, clear it from both target and source list
493 fMatchedJetContainer[i][0] = fMatchedJetContainer[j][0]; // switch to avoid breaking loop
494 fMatchedJetContainer[i][1] = fMatchedJetContainer[j][1];
495 fMatchedJetContainer[j][0] = 0x0; // clear duplicate jet from cache
496 fMatchedJetContainer[j][1] = 0x0;
498 if(fDebug > 0) printf(" found duplicate jet, chose %.2f over %.2f \n" , (rB > rA) ? rA : rB, (rB > rA) ? rB : rA);
503 //_____________________________________________________________________________
504 void AliAnalysisTaskJetMatching::PostMatchedJets()
507 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
508 fMatchedJets->Delete(); // should be a NULL operation, but added just in case
509 for(Int_t i(0), p(0); i < fNoMatchedJets; i++) {
510 if(fMatchedJetContainer[i][1]) { // duplicate jets will have NULL value here and are skipped
511 new((*fMatchedJets)[p]) AliEmcalJet(*fMatchedJetContainer[i][1]);
516 //_____________________________________________________________________________
517 void AliAnalysisTaskJetMatching::FillAnalysisSummaryHistogram() const
519 // fill the analysis summary histrogram, saves all relevant analysis settigns
520 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
521 fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fUseScaledRho");
522 fHistAnalysisSummary->SetBinContent(1, (int)fUseScaledRho);
523 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fMatchingScheme");
524 fHistAnalysisSummary->SetBinContent(2, (int)fMatchingScheme);
525 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fSourceBKG");
526 fHistAnalysisSummary->SetBinContent(3, (int)fSourceBKG);
527 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fTargetBKG");
528 fHistAnalysisSummary->SetBinContent(4, (int)fTargetBKG);
529 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fMatchEta");
530 fHistAnalysisSummary->SetBinContent(5, fMatchEta);
531 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fMatchPhi");
532 fHistAnalysisSummary->SetBinContent(6, fMatchPhi);
533 fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fMatchR");
534 fHistAnalysisSummary->SetBinContent(7, fMatchR);
535 fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "iter");
536 fHistAnalysisSummary->SetBinContent(8, 1.);
538 //_____________________________________________________________________________
539 void AliAnalysisTaskJetMatching::FillMatchedJetHistograms()
541 // fill matched jet histograms
542 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
543 for(Int_t i(0); i < fSourceJets->GetEntriesFast(); i++) {
544 AliEmcalJet* source = static_cast<AliEmcalJet*>(fSourceJets->At(i));
545 if(!source) continue;
546 else if(fUseEmcalBaseJetCuts && !AcceptJet(source, 0)) continue;
547 Double_t sourceRho(0), targetRho(0);
548 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(source->Phi(), fSourceRadius)*source->Area();
549 fHistSourceJetPt->Fill(source->Pt()-sourceRho);
550 fHistNoConstSourceJet->Fill(source->GetNumberOfConstituents());
551 for(Int_t j(0); j < fTargetJets->GetEntriesFast(); j++) {
552 AliEmcalJet* target = static_cast<AliEmcalJet*>(fTargetJets->At(j));
554 if(fUseEmcalBaseJetCuts && !AcceptJet(target, 1)) continue;
555 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(target->Phi(), fTargetRadius)*target->Area();
556 fProfQA->Fill(0.5, TMath::Abs((source->Pt()-sourceRho)-(target->Pt()-targetRho)));
557 fProfQA->Fill(1.5, TMath::Abs(source->Eta()-target->Eta()));
558 fProfQA->Fill(2.5, TMath::Abs(source->Phi()-target->Phi()));
559 fHistUnsortedCorrelation->Fill(PhaseShift(source->Phi()-target->Phi(), 2));
561 fHistTargetJetPt->Fill(target->Pt()-targetRho);
562 fHistNoConstTargetJet->Fill(target->GetNumberOfConstituents());
567 for(Int_t i(0); i < fNoMatchedJets; i++) {
568 if(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1]) {
569 Double_t sourceRho(0), targetRho(0);
570 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(fMatchedJetContainer[i][0]->Phi(), fSourceRadius)*fMatchedJetContainer[i][0]->Area();
571 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(fMatchedJetContainer[i][1]->Phi(), fTargetRadius)*fMatchedJetContainer[i][1]->Area();
572 fHistMatchedCorrelation->Fill(PhaseShift(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi(), 2));
573 fHistMatchedSourceJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho);
574 fHistMatchedJetPt->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
575 fHistNoConstMatchJet->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
576 fProfQAMatched->Fill(0.5, TMath::Abs((fMatchedJetContainer[i][0]->Pt()-sourceRho)-(fMatchedJetContainer[i][1]->Pt()-targetRho)));
577 fProfQAMatched->Fill(1.5, TMath::Abs(fMatchedJetContainer[i][0]->Eta()-fMatchedJetContainer[i][1]->Eta()));
578 fProfQAMatched->Fill(2.5, TMath::Abs(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi()));
580 fHistSourceMatchedJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, fMatchedJetContainer[i][1]->Pt()-targetRho);
581 if(fDoDetectorResponse) {
582 fProfFracPtJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (fMatchedJetContainer[i][1]->Pt()-targetRho) / (fMatchedJetContainer[i][0]->Pt()-sourceRho));
583 fProfFracNoJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (double)fMatchedJetContainer[i][1]->GetNumberOfTracks() / (double)fMatchedJetContainer[i][0]->GetNumberOfTracks());
584 fHistDetectorResponseProb->Fill((fMatchedJetContainer[i][1]->Pt()-fMatchedJetContainer[i][0]->Pt())/fMatchedJetContainer[i][0]->Pt(), fMatchedJetContainer[i][0]->Pt());
589 //_____________________________________________________________________________
590 AliEmcalJet* AliAnalysisTaskJetMatching::GetLeadingJet(TClonesArray* source, Int_t &leadingJetIndex, Double_t etaMin, Double_t etaMax)
592 // return the leading jet within given acceptance
593 Int_t iJets(source->GetEntriesFast());
595 AliEmcalJet* leadingJet(0x0);
596 for(Int_t i(0); i < iJets; i++) {
597 AliEmcalJet* jet = static_cast<AliEmcalJet*>(source->At(i));
598 if(jet->Eta() < etaMin || jet->Eta() > etaMax) continue;
601 pt = leadingJet->Pt();
607 //_____________________________________________________________________________
608 AliEmcalJet* AliAnalysisTaskJetMatching::GetSubLeadingJet(TClonesArray* source, Int_t leadingJetIndex, Int_t &subLeadingJetIndex, Double_t etaMin, Double_t etaMax)
610 // return the leading jet within given acceptance
611 // same as GetLeadingJet() but skips the leading jet (so returned jet is
612 // sub-leading by design)
613 Int_t iJets(source->GetEntriesFast());
615 AliEmcalJet* leadingJet(0x0);
616 for(Int_t i(0); i < iJets; i++) {
617 if(i == leadingJetIndex) continue;
618 AliEmcalJet* jet = static_cast<AliEmcalJet*>(source->At(i));
619 if(jet->Eta() < etaMin || jet->Eta() > etaMax) continue;
622 pt = leadingJet->Pt();
623 subLeadingJetIndex = i;
628 //_____________________________________________________________________________
629 void AliAnalysisTaskJetMatching::PrintInfo() const
632 printf("\n > No. of source jets from %s \n \t %i \n ", fSourceJetsName.Data(), fSourceJets->GetEntriesFast());
633 printf(" > No. of target jets from %s \n \t %i \n ", fTargetJetsName.Data(), fTargetJets->GetEntriesFast());
634 printf(" > No. of matched jets from %s \n \t %i \n ", fMatchedJetsName.Data(), fMatchedJets->GetEntriesFast());
635 if(fDebug > 3) InputEvent()->GetList()->ls();
637 //_____________________________________________________________________________
638 void AliAnalysisTaskJetMatching::Terminate(Option_t *)
642 //_____________________________________________________________________________