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), fHistDiJetDPhi(0), fHistDiJetDPt(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), fHistDiJetDPhi(0), fHistDiJetDPt(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 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 fHistDiJetDPt = BookTH1F("fHistDiJetDPt", "leading jet p_{T} - sub leading jet p_{T} (GeV/c)", 100, -25, 25);
198 PostData(1, fOutputList);
200 //_____________________________________________________________________________
201 TH1F* AliAnalysisTaskJetMatching::BookTH1F(const char* name, const char* x, Int_t bins, Double_t min, Double_t max, Bool_t append)
203 // book a TH1F and connect it to the output container
204 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
205 if(!fOutputList) return 0x0;
207 title += Form(";%s;[counts]", x);
208 TH1F* histogram = new TH1F(name, title.Data(), bins, min, max);
210 if(append) fOutputList->Add(histogram);
213 //_____________________________________________________________________________
214 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)
216 // book a TH2F and connect it to the output container
217 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
218 if(!fOutputList) return 0x0;
220 title += Form(";%s;%s", x, y);
221 TH2F* histogram = new TH2F(name, title.Data(), binsx, minx, maxx, binsy, miny, maxy);
223 if(append) fOutputList->Add(histogram);
226 //_____________________________________________________________________________
227 Bool_t AliAnalysisTaskJetMatching::Notify()
229 // for each file get the number of trials and pythia cross section
230 // see $ALICE_ROOT/PWGJE/AliAnalysisTaskJetProperties.cxx
231 // $ALICE_ROOT/PWG/Tools/AliAnalysisHelperJetTasks.cxx
232 // this function is implenented here temporarily to avoid introducing a dependency
233 // later on this could just be a call to a static helper function
234 if(!fDoDetectorResponse) return kTRUE;
235 Float_t xsection(0), ftrials(1);
236 fAvgTrials = ftrials;
237 TTree *tree = AliAnalysisManager::GetAnalysisManager()->GetTree();
239 TFile *curfile = tree->GetCurrentFile();
240 if (!curfile) return kFALSE;
241 TString file(curfile->GetName());
242 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),"");
243 else file.ReplaceAll(gSystem->BaseName(file.Data()),"");
244 TFile *fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec.root"));
246 fxsec = TFile::Open(Form("%s%s",file.Data(),"pyxsec_hists.root"));
248 TKey* key = (TKey*)fxsec->GetListOfKeys()->At(0);
250 TList *list = dynamic_cast<TList*>(key->ReadObj());
252 xsection = ((TProfile*)list->FindObject("h1Xsec"))->GetBinContent(1);
253 ftrials = ((TH1F*)list->FindObject("h1Trials"))->GetBinContent(1);
259 TTree *xtree = (TTree*)fxsec->Get("Xsection");
262 Double_t _xsection = 0;
263 xtree->SetBranchAddress("xsection",&_xsection);
264 xtree->SetBranchAddress("ntrials",&_ftrials);
267 xsection = _xsection;
271 fh1Xsec->Fill("<#sigma>",xsection);
272 Float_t nEntries = (Float_t)tree->GetTree()->GetEntries();
273 if(ftrials >= nEntries && nEntries > 0.) fAvgTrials = ftrials/nEntries;
274 fh1Trials->Fill("#sum{ntrials}",ftrials);
278 //_____________________________________________________________________________
279 Bool_t AliAnalysisTaskJetMatching::Run()
281 // execute once for each event
282 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
283 if(!(InputEvent() && fSourceJets && fTargetJets && IsEventSelected())) return kFALSE;
284 if(fh1AvgTrials) fh1AvgTrials->Fill("#sum{avg ntrials}", fAvgTrials);
285 // step one: do geometric matching
286 switch (fMatchingScheme) {
287 // FIXME having separate dedicated functions is not necessary and historical
288 // cluttered code - should be cleaned up and merged into one
290 DoGeometricMatchingEtaPhi();
291 // break if no jet was matched
292 if(!fMatchedJetContainer[1][0]) return kTRUE;
295 DoGeometricMatchingR();
296 // break if no jet was matched
297 if(!fMatchedJetContainer[1][0]) return kTRUE;
304 // optional step two: get a bijection (avoid duplicate matches)
305 if(fGetBijection) GetBijection();
306 // optional step three: match constituents within matched jets
307 if(fMatchConstituents) DoConstituentMatching();
308 // stream data to output
310 FillMatchedJetHistograms();
311 if(fDebug > 0) PrintInfo();
312 PostData(1, fOutputList);
315 //_____________________________________________________________________________
316 void AliAnalysisTaskJetMatching::DoGeometricMatchingEtaPhi()
318 // do geometric matching based on eta phi distance between jet centers
319 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
320 fNoMatchedJets = 0; // reset the matched jet counter
321 Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
322 for(Int_t i(0); i < iSource; i++) {
323 AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
324 if(!PassesCuts(sourceJet, 0)) continue;
325 for(Int_t j(0); j < iTarget; j++) {
326 AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
327 if(!PassesCuts(targetJet, 1)) continue;
328 if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
329 if((TMath::Abs(sourceJet->Eta() - targetJet->Eta()) < fMatchEta )) {
330 Double_t sourcePhi(sourceJet->Phi()), targetPhi(targetJet->Phi());
331 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
332 if(TMath::Abs(sourcePhi-targetPhi) > TMath::Abs(sourcePhi-targetPhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
333 if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) { // accept the jets as matching
334 Bool_t isBestMatch(kTRUE);
335 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
336 if(fDebug > 0) printf(" > Entering first bijection test \n");
337 for(Int_t k(i); k < iSource; k++) {
338 AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
339 if(PassesCuts(candidateSourceJet, 0)) continue;
340 if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
341 if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
342 isBestMatch = kFALSE;
346 if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
349 fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
350 fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
353 if(fNoMatchedJets > 199) { // maximum to keep the cache at reasonable size
354 AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
362 //_____________________________________________________________________________
363 void AliAnalysisTaskJetMatching::DoGeometricMatchingR()
365 // do geometric matching based on shortest path between jet centers
366 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
367 fNoMatchedJets = 0; // reset the matched jet counter
368 Int_t iSource(fSourceJets->GetEntriesFast()), iTarget(fTargetJets->GetEntriesFast());
369 for(Int_t i(0); i < iSource; i++) {
370 AliEmcalJet* sourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(i)));
371 if(!PassesCuts(sourceJet, 0)) continue;
372 for(Int_t j(0); j < iTarget; j++) {
373 AliEmcalJet* targetJet(static_cast<AliEmcalJet*>(fTargetJets->At(j)));
374 if(!PassesCuts(targetJet, 1)) continue;
375 else if (fUseEmcalBaseJetCuts && !AcceptJet(targetJet, 1)) continue;
376 if(GetR(sourceJet, targetJet) <= fMatchR) {
377 Bool_t isBestMatch(kTRUE);
378 if(fGetBijection) { // match has been found, for bijection only store it if there's no better match
379 if(fDebug > 0) printf(" > Entering first bijection test \n");
380 for(Int_t k(i); k < iSource; k++) {
381 AliEmcalJet* candidateSourceJet(static_cast<AliEmcalJet*>(fSourceJets->At(k)));
382 if(!PassesCuts(candidateSourceJet, 0)) continue;
383 if(fDebug > 0) printf("source distance %.2f \t candidate distance %.2f \n", GetR(sourceJet, targetJet),GetR(candidateSourceJet, targetJet));
384 if(GetR(sourceJet, targetJet) > GetR(candidateSourceJet, targetJet)) {
385 isBestMatch = kFALSE;
389 if(fDebug > 0) (isBestMatch) ? printf(" kept source \n ") : printf(" we can do better (rejected source) \n");
392 fMatchedJetContainer[fNoMatchedJets][0] = sourceJet;
393 fMatchedJetContainer[fNoMatchedJets][1] = targetJet;
396 if(fNoMatchedJets > 99) {
397 AliError(Form("%s: Found too many matched jets (> 100). Adjust matching criteria !", GetName()));
404 //_____________________________________________________________________________
405 void AliAnalysisTaskJetMatching::DoDiJetMatching()
407 // match dijets. this is in a sense a 'special' mode of the task as both jet supplied jet arrays
408 // (target and source) are the same jet collection, matching will be performed on distribution in
410 // no ouptut array is produced, dedicated dijet histo's are filled in this function instead
411 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
412 // get the leading jet in a given acceptance (TPC is default)
413 Int_t leadingJetIndex(-1), subLeadingJetIndex(-1);
414 // retrieve the leading jet, leadingJetIndex points to the leading jet
415 AliEmcalJet* leadingJet(GetLeadingJet(fSourceJets, leadingJetIndex));
416 if(!leadingJet) return;
417 // fill phi and eta of leading jet (should always be in selected acceptance)
418 fHistDiJetLeadingJet->Fill(leadingJet->Phi(), leadingJet->Eta());
419 Double_t sourcePhi(leadingJet->Phi()), targetPhi(-1);
420 // get the sub-leading jet - faster when leading jet is also provided
421 AliEmcalJet* subLeadingJet(GetSubLeadingJet(fSourceJets, leadingJetIndex, subLeadingJetIndex));
422 if(!subLeadingJet) return;
423 else { // check if the sub-leading jet is actually the away-side jet
424 targetPhi = subLeadingJet->Phi() + TMath::Pi();
425 // rotate jets to common phase
426 if(TMath::Abs(sourcePhi) > TMath::Abs(sourcePhi+TMath::TwoPi())) sourcePhi+=TMath::TwoPi();
427 if(TMath::Abs(sourcePhi) > TMath::Abs(sourcePhi-TMath::TwoPi())) sourcePhi-=TMath::TwoPi();
428 if(TMath::Abs(targetPhi) > TMath::Abs(targetPhi+TMath::TwoPi())) targetPhi+=TMath::TwoPi();
429 if(TMath::Abs(targetPhi) > TMath::Abs(targetPhi-TMath::TwoPi())) targetPhi-=TMath::TwoPi();
430 if(TMath::Abs(sourcePhi-targetPhi) < fMatchPhi) {
431 fHistDiJet->Fill(subLeadingJet->Phi(), subLeadingJet->Eta());
432 fHistDiJetDPhi->Fill(sourcePhi-targetPhi);
433 fHistDiJetDPt->Fill(leadingJet->Pt() - subLeadingJet->Pt());
437 //_____________________________________________________________________________
438 void AliAnalysisTaskJetMatching::DoConstituentMatching()
440 // match constituents within matched jets
441 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
443 AliFatal(Form("%s Fatal error! To do deep matching, supply jet constituents ! \n", GetName()));
444 return; // coverity ...
446 for(Int_t i(0); i < fNoMatchedJets; i++) {
447 AliEmcalJet* sourceJet = fMatchedJetContainer[i][0];
448 AliEmcalJet* targetJet = fMatchedJetContainer[i][1];
449 if(sourceJet && targetJet) { // duplicate check: slot migth be NULL
450 Double_t targetPt(0);
451 Int_t iSJ(sourceJet->GetNumberOfTracks());
452 Int_t iTJ(targetJet->GetNumberOfTracks());
453 Int_t overlap(0), alreadyFound(0);
454 for(Int_t j(0); j < iSJ; j++) {
456 Int_t idSource((Int_t)sourceJet->TrackAt(j));
457 for(Int_t k(0); k < iTJ; k++) { // compare all target tracks to the source track
458 if(idSource == targetJet->TrackAt(k) && alreadyFound == 0) {
460 alreadyFound++; // avoid possible duplicate matching
461 AliVParticle* vp(static_cast<AliVParticle*>(targetJet->TrackAt(k, fTracks)));
462 if(vp) targetPt += vp->Pt();
467 if((float)overlap/(float)iSJ < fMinFracRecoveredConstituents || targetPt/sourceJet->Pt() < fMinFracRecoveredConstituentPt) {
468 if(fDebug > 0) printf(" \n > Purging jet, recovered constituents ratio %i / %i = %.2f < or pt ratio %.2f / %.2f = %.2f < %.2f",
469 overlap, iSJ, (float)overlap/(float)iSJ, targetPt, sourceJet->Pt(), targetPt/sourceJet->Pt(), fMinFracRecoveredConstituentPt);
470 fMatchedJetContainer[i][0] = 0x0;
471 fMatchedJetContainer[i][1] = 0x0;
474 if(sourceJet->Pt() > 0) {
475 Double_t sourceRho(0), targetRho(0);
476 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(sourceJet->Phi(), fSourceRadius)*sourceJet->Area();
477 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(targetJet->Phi(), fTargetRadius)*targetJet->Area();
478 fProfFracPtMatched->Fill(sourceJet->Pt()-sourceRho, (targetPt-targetRho) / (sourceJet->Pt()-sourceRho));
479 fProfFracPtJets->Fill(sourceJet->Pt()-sourceRho, (targetJet->Pt()-targetRho) / (sourceJet->Pt()-sourceRho));
480 fProfFracNoMatched->Fill(sourceJet->Pt()-sourceRho, (double)overlap / (double)sourceJet->GetNumberOfTracks());
481 fProfFracNoJets->Fill(sourceJet->Pt()-sourceRho, (double)targetJet->GetNumberOfTracks() / (double)sourceJet->GetNumberOfTracks());
484 printf("\n > Jet A: %i const\t", iSJ);
485 printf(" > Jet B %i const\t", iTJ);
486 printf(" -> OVERLAP: %i tracks <- \n", overlap);
491 //_____________________________________________________________________________
492 void AliAnalysisTaskJetMatching::GetBijection()
494 // bijection of source and matched jets, based on closest distance between jets
495 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
496 for(Int_t i(0); i < fNoMatchedJets; i++) {
497 for(Int_t j(i+1); j < fNoMatchedJets; j++) {
498 if(fMatchedJetContainer[i][0] == fMatchedJetContainer[j][0]) {
499 // found source with two targets, now see which target is closer to the source
500 if(!(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1] && fMatchedJetContainer[j][0] && fMatchedJetContainer[j][1] )) continue;
501 Double_t rA(GetR(fMatchedJetContainer[i][0], fMatchedJetContainer[i][1])); // distance between connection A = i
502 Double_t rB(GetR(fMatchedJetContainer[j][0], fMatchedJetContainer[j][1])); // distance between connection B = j
503 if (rB > rA) { // jet two is far away, clear it from both target and source list
504 fMatchedJetContainer[j][0] = 0x0;
505 fMatchedJetContainer[j][1] = 0x0;
506 } else { // jet one is far away, clear it from both target and source list
507 fMatchedJetContainer[i][0] = fMatchedJetContainer[j][0]; // switch to avoid breaking loop
508 fMatchedJetContainer[i][1] = fMatchedJetContainer[j][1];
509 fMatchedJetContainer[j][0] = 0x0; // clear duplicate jet from cache
510 fMatchedJetContainer[j][1] = 0x0;
512 if(fDebug > 0) printf(" found duplicate jet, chose %.2f over %.2f \n" , (rB > rA) ? rA : rB, (rB > rA) ? rB : rA);
517 //_____________________________________________________________________________
518 void AliAnalysisTaskJetMatching::PostMatchedJets()
521 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
522 fMatchedJets->Delete(); // should be a NULL operation, but added just in case
523 for(Int_t i(0), p(0); i < fNoMatchedJets; i++) {
524 if(fMatchedJetContainer[i][1]) { // duplicate jets will have NULL value here and are skipped
525 new((*fMatchedJets)[p]) AliEmcalJet(*fMatchedJetContainer[i][1]);
530 //_____________________________________________________________________________
531 void AliAnalysisTaskJetMatching::FillAnalysisSummaryHistogram() const
533 // fill the analysis summary histrogram, saves all relevant analysis settigns
534 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
535 fHistAnalysisSummary->GetXaxis()->SetBinLabel(1, "fUseScaledRho");
536 fHistAnalysisSummary->SetBinContent(1, (int)fUseScaledRho);
537 fHistAnalysisSummary->GetXaxis()->SetBinLabel(2, "fMatchingScheme");
538 fHistAnalysisSummary->SetBinContent(2, (int)fMatchingScheme);
539 fHistAnalysisSummary->GetXaxis()->SetBinLabel(3, "fSourceBKG");
540 fHistAnalysisSummary->SetBinContent(3, (int)fSourceBKG);
541 fHistAnalysisSummary->GetXaxis()->SetBinLabel(4, "fTargetBKG");
542 fHistAnalysisSummary->SetBinContent(4, (int)fTargetBKG);
543 fHistAnalysisSummary->GetXaxis()->SetBinLabel(5, "fMatchEta");
544 fHistAnalysisSummary->SetBinContent(5, fMatchEta);
545 fHistAnalysisSummary->GetXaxis()->SetBinLabel(6, "fMatchPhi");
546 fHistAnalysisSummary->SetBinContent(6, fMatchPhi);
547 fHistAnalysisSummary->GetXaxis()->SetBinLabel(7, "fMatchR");
548 fHistAnalysisSummary->SetBinContent(7, fMatchR);
549 fHistAnalysisSummary->GetXaxis()->SetBinLabel(8, "iter");
550 fHistAnalysisSummary->SetBinContent(8, 1.);
552 //_____________________________________________________________________________
553 void AliAnalysisTaskJetMatching::FillMatchedJetHistograms()
555 // fill matched jet histograms
556 if(fDebug > 0) printf("__FILE__ = %s \n __LINE __ %i , __FUNC__ %s \n ", __FILE__, __LINE__, __func__);
557 for(Int_t i(0); i < fSourceJets->GetEntriesFast(); i++) {
558 AliEmcalJet* source = static_cast<AliEmcalJet*>(fSourceJets->At(i));
559 if(!source) continue;
560 else if(fUseEmcalBaseJetCuts && !AcceptJet(source, 0)) continue;
561 Double_t sourceRho(0), targetRho(0);
562 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(source->Phi(), fSourceRadius)*source->Area();
563 fHistSourceJetPt->Fill(source->Pt()-sourceRho);
564 fHistNoConstSourceJet->Fill(source->GetNumberOfConstituents());
565 for(Int_t j(0); j < fTargetJets->GetEntriesFast(); j++) {
566 AliEmcalJet* target = static_cast<AliEmcalJet*>(fTargetJets->At(j));
568 if(fUseEmcalBaseJetCuts && !AcceptJet(target, 1)) continue;
569 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(target->Phi(), fTargetRadius)*target->Area();
570 fProfQA->Fill(0.5, TMath::Abs((source->Pt()-sourceRho)-(target->Pt()-targetRho)));
571 fProfQA->Fill(1.5, TMath::Abs(source->Eta()-target->Eta()));
572 fProfQA->Fill(2.5, TMath::Abs(source->Phi()-target->Phi()));
573 fHistUnsortedCorrelation->Fill(PhaseShift(source->Phi()-target->Phi(), 2));
575 fHistTargetJetPt->Fill(target->Pt()-targetRho);
576 fHistNoConstTargetJet->Fill(target->GetNumberOfConstituents());
581 for(Int_t i(0); i < fNoMatchedJets; i++) {
582 if(fMatchedJetContainer[i][0] && fMatchedJetContainer[i][1]) {
583 Double_t sourceRho(0), targetRho(0);
584 if(fSourceRho) sourceRho = fSourceRho->GetLocalVal(fMatchedJetContainer[i][0]->Phi(), fSourceRadius)*fMatchedJetContainer[i][0]->Area();
585 if(fTargetRho) targetRho = fTargetRho->GetLocalVal(fMatchedJetContainer[i][1]->Phi(), fTargetRadius)*fMatchedJetContainer[i][1]->Area();
586 fHistMatchedCorrelation->Fill(PhaseShift(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi(), 2));
587 fHistMatchedSourceJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho);
588 fHistMatchedJetPt->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
589 fHistNoConstMatchJet->Fill(fMatchedJetContainer[i][1]->Pt()-targetRho);
590 fProfQAMatched->Fill(0.5, TMath::Abs((fMatchedJetContainer[i][0]->Pt()-sourceRho)-(fMatchedJetContainer[i][1]->Pt()-targetRho)));
591 fProfQAMatched->Fill(1.5, TMath::Abs(fMatchedJetContainer[i][0]->Eta()-fMatchedJetContainer[i][1]->Eta()));
592 fProfQAMatched->Fill(2.5, TMath::Abs(fMatchedJetContainer[i][0]->Phi()-fMatchedJetContainer[i][1]->Phi()));
594 fHistSourceMatchedJetPt->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, fMatchedJetContainer[i][1]->Pt()-targetRho);
595 if(fDoDetectorResponse) {
596 fProfFracPtJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (fMatchedJetContainer[i][1]->Pt()-targetRho) / (fMatchedJetContainer[i][0]->Pt()-sourceRho));
597 fProfFracNoJets->Fill(fMatchedJetContainer[i][0]->Pt()-sourceRho, (double)fMatchedJetContainer[i][1]->GetNumberOfTracks() / (double)fMatchedJetContainer[i][0]->GetNumberOfTracks());
598 fHistDetectorResponseProb->Fill((fMatchedJetContainer[i][1]->Pt()-fMatchedJetContainer[i][0]->Pt())/fMatchedJetContainer[i][0]->Pt(), fMatchedJetContainer[i][0]->Pt());
603 //_____________________________________________________________________________
604 AliEmcalJet* AliAnalysisTaskJetMatching::GetLeadingJet(TClonesArray* source, Int_t &leadingJetIndex, Double_t etaMin, Double_t etaMax)
606 // return the leading jet within giiven acceptance
607 Int_t iJets(source->GetEntriesFast());
609 AliEmcalJet* leadingJet(0x0);
610 for(Int_t i(0); i < iJets; i++) {
611 AliEmcalJet* jet = static_cast<AliEmcalJet*>(source->At(i));
612 if(jet->Eta() < etaMin || jet->Eta() > etaMax) continue;
615 pt = leadingJet->Pt();
621 //_____________________________________________________________________________
622 AliEmcalJet* AliAnalysisTaskJetMatching::GetSubLeadingJet(TClonesArray* source, Int_t leadingJetIndex, Int_t &subLeadingJetIndex)
624 // return the sub-leading jet within given acceptance
625 // same as GetLeadingJet() but skips the leading jet (so returned jet is
626 // sub-leading by design)
627 // there is no eta requirement on the location of the sub-leading jet
628 Int_t iJets(source->GetEntriesFast());
629 // if the leading jet isn't given, retrieve it
630 if(leadingJetIndex < 0) GetLeadingJet(source, leadingJetIndex);
631 AliEmcalJet *leadingJet(0x0), *prevLeadingJet(static_cast<AliEmcalJet*>(source->At(leadingJetIndex)));
632 if(!prevLeadingJet) return 0x0;
633 Double_t pt(0), leadingPt(prevLeadingJet->Pt());
634 for(Int_t i(0); i < iJets; i++) {
635 if(i == leadingJetIndex) continue;
636 AliEmcalJet* jet = static_cast<AliEmcalJet*>(source->At(i));
637 // check if jet is actually sub-leading
638 if(jet->Pt() > pt && jet->Pt() <= leadingPt) {
640 pt = leadingJet->Pt();
641 subLeadingJetIndex = i;
646 //_____________________________________________________________________________
647 void AliAnalysisTaskJetMatching::PrintInfo() const
650 printf("\n > No. of source jets from %s \n \t %i \n ", fSourceJetsName.Data(), fSourceJets->GetEntriesFast());
651 printf(" > No. of target jets from %s \n \t %i \n ", fTargetJetsName.Data(), fTargetJets->GetEntriesFast());
652 printf(" > No. of matched jets from %s \n \t %i \n ", fMatchedJetsName.Data(), fMatchedJets->GetEntriesFast());
653 if(fDebug > 3) InputEvent()->GetList()->ls();
655 //_____________________________________________________________________________
656 void AliAnalysisTaskJetMatching::Terminate(Option_t *)
660 //_____________________________________________________________________________