AliPhysics  97a96ce (97a96ce)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliEmcalJetTask.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-2016, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 
16 #include <vector>
17 
18 #include <TClonesArray.h>
19 #include <TMath.h>
20 #include <TRandom3.h>
21 
22 #include <AliVCluster.h>
23 #include <AliVEvent.h>
24 #include <AliVParticle.h>
25 #include <AliEMCALGeometry.h>
26 
27 #include "AliTLorentzVector.h"
28 #include "AliEmcalJet.h"
29 #include "AliEmcalParticle.h"
30 #include "AliFJWrapper.h"
31 #include "AliEmcalJetUtility.h"
32 #include "AliParticleContainer.h"
33 #include "AliClusterContainer.h"
34 
35 #include "AliEmcalJetTask.h"
36 
37 using std::cout;
38 using std::endl;
39 using std::cerr;
40 
44 
46 
53  fJetsTag(),
54  fJetAlgo(AliJetContainer::antikt_algorithm),
55  fJetType(AliJetContainer::kFullJet),
56  fRecombScheme(AliJetContainer::pt_scheme),
57  fRadius(0.4),
58  fMinJetArea(0.001),
59  fMinJetPt(1.0),
60  fJetPhiMin(-10),
61  fJetPhiMax(+10),
62  fJetEtaMin(-1),
63  fJetEtaMax(+1),
64  fGhostArea(0.005),
65  fTrackEfficiency(1.),
66  fUtilities(0),
67  fLocked(0),
68  fJetsName(),
69  fIsInit(0),
70  fIsPSelSet(0),
71  fIsEmcPart(0),
72  fLegacyMode(kFALSE),
73  fFillGhost(kFALSE),
74  fJets(0),
75  fFastJetWrapper("AliEmcalJetTask","AliEmcalJetTask")
76 {
77 }
78 
85  fJetsTag("Jets"),
86  fJetAlgo(AliJetContainer::antikt_algorithm),
87  fJetType(AliJetContainer::kFullJet),
88  fRecombScheme(AliJetContainer::pt_scheme),
89  fRadius(0.4),
90  fMinJetArea(0.001),
91  fMinJetPt(1.0),
92  fJetPhiMin(-10),
93  fJetPhiMax(+10),
94  fJetEtaMin(-1),
95  fJetEtaMax(+1),
96  fGhostArea(0.005),
97  fTrackEfficiency(1.),
98  fUtilities(0),
99  fLocked(0),
100  fJetsName(),
101  fIsInit(0),
102  fIsPSelSet(0),
103  fIsEmcPart(0),
104  fLegacyMode(kFALSE),
105  fFillGhost(kFALSE),
106  fJets(0),
107  fFastJetWrapper(name,name)
108 {
109 }
110 
115 {
116 }
117 
123 {
124  if (!fUtilities) fUtilities = new TObjArray();
125  if (fUtilities->FindObject(utility)) {
126  Error("AddUtility", "Jet utility %s already connected.", utility->GetName());
127  return utility;
128  }
129  fUtilities->Add(utility);
130  utility->SetJetTask(this);
131 
132  return utility;
133 }
134 
140 {
141  TIter next(fUtilities);
142  AliEmcalJetUtility *utility = 0;
143  while ((utility=static_cast<AliEmcalJetUtility*>(next()))) utility->Init();
144 }
145 
152 {
153  TIter next(fUtilities);
154  AliEmcalJetUtility *utility = 0;
155  while ((utility=static_cast<AliEmcalJetUtility*>(next()))) utility->Prepare(fFastJetWrapper);
156 }
157 
163 {
164  TIter next(fUtilities);
165  AliEmcalJetUtility *utility = 0;
166  while ((utility=static_cast<AliEmcalJetUtility*>(next()))) utility->ProcessJet(jet, ij, fFastJetWrapper);
167 }
168 
174 {
175  TIter next(fUtilities);
176  AliEmcalJetUtility *utility = 0;
177  while ((utility=static_cast<AliEmcalJetUtility*>(next()))) utility->Terminate(fFastJetWrapper);
178 }
179 
185 {
186  // clear the jet array (normally a null operation)
187  fJets->Delete();
188 
189  Int_t n = FindJets();
190 
191  if (n == 0) return kFALSE;
192 
193  FillJetBranch();
194 
195  return kTRUE;
196 }
197 
206 {
207  if (fParticleCollArray.GetEntriesFast() == 0 && fClusterCollArray.GetEntriesFast() == 0){
208  AliError("No tracks or clusters, returning.");
209  return 0;
210  }
211 
213 
214  AliDebug(2,Form("Jet type = %d", fJetType));
215 
216  Int_t iColl = 1;
217  TIter nextPartColl(&fParticleCollArray);
218  AliParticleContainer* tracks = 0;
219  while ((tracks = static_cast<AliParticleContainer*>(nextPartColl()))) {
220  AliDebug(2,Form("Tracks from collection %d: '%s'.", iColl-1, tracks->GetName()));
222  for (AliParticleIterableMomentumContainer::iterator it = itcont.begin(); it != itcont.end(); it++) {
223  // artificial inefficiency
224  if (fTrackEfficiency < 1.) {
225  Double_t rnd = gRandom->Rndm();
226  if (fTrackEfficiency < rnd) {
227  AliDebug(2,Form("Track %d rejected due to artificial tracking inefficiency", it.current_index()));
228  continue;
229  }
230  }
231 
232  AliDebug(2,Form("Track %d accepted (label = %d, pt = %f)", it.current_index(), it->second->GetLabel(), it->first.Pt()));
233  Int_t uid = it.current_index() + fgkConstIndexShift * iColl;
234  fFastJetWrapper.AddInputVector(it->first.Px(), it->first.Py(), it->first.Pz(), it->first.E(), uid);
235  }
236  iColl++;
237  }
238 
239  iColl = 1;
240  TIter nextClusColl(&fClusterCollArray);
241  AliClusterContainer* clusters = 0;
242  while ((clusters = static_cast<AliClusterContainer*>(nextClusColl()))) {
243  AliDebug(2,Form("Clusters from collection %d: '%s'.", iColl-1, clusters->GetName()));
245  for (AliClusterIterableMomentumContainer::iterator it = itcont.begin(); it != itcont.end(); it++) {
246  AliDebug(2,Form("Cluster %d accepted (label = %d, energy = %.3f)", it.current_index(), it->second->GetLabel(), it->first.E()));
247  Int_t uid = -it.current_index() - fgkConstIndexShift * iColl;
248  fFastJetWrapper.AddInputVector(it->first.Px(), it->first.Py(), it->first.Pz(), it->first.E(), uid);
249  }
250  iColl++;
251  }
252 
253  if (fFastJetWrapper.GetInputVectors().size() == 0) return 0;
254 
255  // run jet finder
257 
258  return fFastJetWrapper.GetInclusiveJets().size();
259 }
260 
267 {
269 
270  // loop over fastjet jets
271  std::vector<fastjet::PseudoJet> jets_incl = fFastJetWrapper.GetInclusiveJets();
272  // sort jets according to jet pt
273  static Int_t indexes[9999] = {-1};
274  GetSortedArray(indexes, jets_incl);
275 
276  AliDebug(1,Form("%d jets found", (Int_t)jets_incl.size()));
277  for (UInt_t ijet = 0, jetCount = 0; ijet < jets_incl.size(); ++ijet) {
278  Int_t ij = indexes[ijet];
279  AliDebug(3,Form("Jet pt = %f, area = %f", jets_incl[ij].perp(), fFastJetWrapper.GetJetArea(ij)));
280 
281  if (jets_incl[ij].perp() < fMinJetPt) continue;
282  if (fFastJetWrapper.GetJetArea(ij) < fMinJetArea) continue;
283  if ((jets_incl[ij].eta() < fJetEtaMin) || (jets_incl[ij].eta() > fJetEtaMax) ||
284  (jets_incl[ij].phi() < fJetPhiMin) || (jets_incl[ij].phi() > fJetPhiMax))
285  continue;
286 
287  AliEmcalJet *jet = new ((*fJets)[jetCount])
288  AliEmcalJet(jets_incl[ij].perp(), jets_incl[ij].eta(), jets_incl[ij].phi(), jets_incl[ij].m());
289  jet->SetLabel(ij);
290 
291  fastjet::PseudoJet area(fFastJetWrapper.GetJetAreaVector(ij));
292  jet->SetArea(area.perp());
293  jet->SetAreaEta(area.eta());
294  jet->SetAreaPhi(area.phi());
295  jet->SetAreaE(area.E());
297 
298  // Fill constituent info
299  std::vector<fastjet::PseudoJet> constituents(fFastJetWrapper.GetJetConstituents(ij));
300  FillJetConstituents(jet, constituents, constituents);
301 
302  if (fGeom) {
303  if ((jet->Phi() > fGeom->GetArm1PhiMin() * TMath::DegToRad()) &&
304  (jet->Phi() < fGeom->GetArm1PhiMax() * TMath::DegToRad()) &&
305  (jet->Eta() > fGeom->GetArm1EtaMin()) &&
306  (jet->Eta() < fGeom->GetArm1EtaMax()))
307  jet->SetAxisInEmcal(kTRUE);
308  }
309 
310  ExecuteUtilities(jet, ij);
311 
312  AliDebug(2,Form("Added jet n. %d, pt = %f, area = %f, constituents = %d", jetCount, jet->Pt(), jet->Area(), jet->GetNumberOfConstituents()));
313  jetCount++;
314  }
315 
317 }
318 
325 Bool_t AliEmcalJetTask::GetSortedArray(Int_t indexes[], std::vector<fastjet::PseudoJet> array) const
326 {
327  static Float_t pt[9999] = {0};
328 
329  const Int_t n = (Int_t)array.size();
330 
331  if (n < 1)
332  return kFALSE;
333 
334  for (Int_t i = 0; i < n; i++)
335  pt[i] = array[i].perp();
336 
337  TMath::Sort(n, pt, indexes);
338 
339  return kTRUE;
340 }
341 
348 {
349  if (fTrackEfficiency < 1.) {
350  if (gRandom) delete gRandom;
351  gRandom = new TRandom3(0);
352  }
353 
355 
356  // add jets to event if not yet there
357  if (!(InputEvent()->FindListObject(fJetsName))) {
358  fJets = new TClonesArray("AliEmcalJet");
359  fJets->SetName(fJetsName);
360  ::Info("AliEmcalJetTask::ExecOnce", "Jet collection with name '%s' has been added to the event.", fJetsName.Data());
361  InputEvent()->AddObject(fJets);
362  }
363  else {
364  AliError(Form("%s: Object with name %s already in event! Returning", GetName(), fJetsName.Data()));
365  return;
366  }
367 
368  // setup fj wrapper
369  fFastJetWrapper.SetAreaType(fastjet::active_area_explicit_ghosts);
375 
376  // setting legacy mode
377  if (fLegacyMode) {
379  }
380 
381  InitUtilities();
382 
384 }
385 
395 void AliEmcalJetTask::FillJetConstituents(AliEmcalJet *jet, std::vector<fastjet::PseudoJet>& constituents,
396  std::vector<fastjet::PseudoJet>& constituents_unsub, Int_t flag, TClonesArray *particles_sub)
397 {
398  Int_t nt = 0;
399  Int_t nc = 0;
400  Double_t neutralE = 0.;
401  Double_t maxCh = 0.;
402  Double_t maxNe = 0.;
403  Int_t gall = 0;
404  Int_t gemc = 0;
405  Int_t cemc = 0;
406  Int_t ncharged = 0;
407  Int_t nneutral = 0;
408  Double_t mcpt = 0.;
409  Double_t emcpt = 0.;
410 
411  Int_t uid = -1;
412 
413  jet->SetNumberOfTracks(constituents.size());
414  jet->SetNumberOfClusters(constituents.size());
415 
416  for (UInt_t ic = 0; ic < constituents.size(); ++ic) {
417 
418  if (flag == 0) {
419  uid = constituents[ic].user_index();
420  }
421  else {
422  if (constituents[ic].perp()<1.e-10) {
423  uid=-1;
424  }
425  else {
426  uid = GetIndexSub(constituents[ic].phi(), constituents_unsub);
427  }
428  if (uid==0) {
429  AliError("correspondence between un/subtracted constituent not found");
430  continue;
431  }
432  }
433 
434  AliDebug(3,Form("Processing constituent %d", uid));
435  if (uid == -1) { //ghost particle
436  ++gall;
437  if (fGeom) {
438  Double_t gphi = constituents[ic].phi();
439  if (gphi < 0) gphi += TMath::TwoPi();
440  gphi *= TMath::RadToDeg();
441  Double_t geta = constituents[ic].eta();
442  if ((gphi > fGeom->GetArm1PhiMin()) && (gphi < fGeom->GetArm1PhiMax()) &&
443  (geta > fGeom->GetArm1EtaMin()) && (geta < fGeom->GetArm1EtaMax()))
444  ++gemc;
445  }
446 
447  if (fFillGhost) jet->AddGhost(constituents[ic].px(),
448  constituents[ic].py(),
449  constituents[ic].pz(),
450  constituents[ic].e());
451  }
452  else if (uid >= fgkConstIndexShift) { // track constituent
453  Int_t iColl = uid / fgkConstIndexShift;
454  Int_t tid = uid - iColl * fgkConstIndexShift;
455  iColl--;
456  AliDebug(3,Form("Constituent %d is a track from collection %d and with ID %d", uid, iColl, tid));
457  AliParticleContainer* partCont = GetParticleContainer(iColl);
458  if (!partCont) {
459  AliError(Form("Could not find particle container %d",iColl));
460  continue;
461  }
462  AliVParticle *t = partCont->GetParticle(tid);
463  if (!t) {
464  AliError(Form("Could not find track %d",tid));
465  continue;
466  }
467 
468  Double_t cEta = t->Eta();
469  Double_t cPhi = t->Phi();
470  Double_t cPt = t->Pt();
471  Double_t cP = t->P();
472  if (t->Charge() == 0) {
473  neutralE += cP;
474  ++nneutral;
475  if (cPt > maxNe) maxNe = cPt;
476  } else {
477  ++ncharged;
478  if (cPt > maxCh) maxCh = cPt;
479  }
480 
481  // check if MC particle
482  if (TMath::Abs(t->GetLabel()) > fMinMCLabel) mcpt += cPt;
483 
484  if (fGeom) {
485  if (cPhi < 0) cPhi += TMath::TwoPi();
486  cPhi *= TMath::RadToDeg();
487  if ((cPhi > fGeom->GetArm1PhiMin()) && (cPhi < fGeom->GetArm1PhiMax()) &&
488  (cEta > fGeom->GetArm1EtaMin()) && (cEta < fGeom->GetArm1EtaMax())) {
489  emcpt += cPt;
490  ++cemc;
491  }
492  }
493 
494  if (flag == 0 || particles_sub == 0) {
495  jet->AddTrackAt(tid, nt);
496  }
497  else {
498  Int_t part_sub_id = particles_sub->GetEntriesFast();
499  AliEmcalParticle* part_sub = new ((*particles_sub)[part_sub_id]) AliEmcalParticle(dynamic_cast<AliVTrack*>(t)); // SA: probably need to be fixed!!
500  part_sub->SetPtEtaPhiM(constituents[ic].perp(),constituents[ic].eta(),constituents[ic].phi(),constituents[ic].m());
501  jet->AddTrackAt(part_sub_id, nt);
502  }
503 
504  ++nt;
505  }
506  else if (uid <= -fgkConstIndexShift) { // cluster constituent
507  Int_t iColl = -uid / fgkConstIndexShift;
508  Int_t cid = -uid - iColl * fgkConstIndexShift;
509  iColl--;
510  AliClusterContainer* clusCont = GetClusterContainer(iColl);
511  AliVCluster *c = clusCont->GetCluster(cid);
512 
513  if (!c) continue;
514 
516  clusCont->GetMomentum(nP, cid);
517 
518  Double_t cEta = nP.Eta();
519  Double_t cPhi = nP.Phi_0_2pi();
520  Double_t cPt = nP.Pt();
521  Double_t cP = nP.P();
522 
523  neutralE += cP;
524  if (cPt > maxNe) maxNe = cPt;
525 
526  // MC particle
527  if (TMath::Abs(c->GetLabel()) > fMinMCLabel) mcpt += c->GetMCEnergyFraction() > 1e-6 ? cPt * c->GetMCEnergyFraction() : cPt;
528 
529  if (fGeom) {
530  if (cPhi<0) cPhi += TMath::TwoPi();
531  cPhi *= TMath::RadToDeg();
532  if ((cPhi > fGeom->GetArm1PhiMin()) && (cPhi < fGeom->GetArm1PhiMax()) &&
533  (cEta > fGeom->GetArm1EtaMin()) && (cEta < fGeom->GetArm1EtaMax())) {
534  emcpt += cPt;
535  ++cemc;
536  }
537  }
538 
539  if (flag == 0 || particles_sub == 0) {
540  jet->AddClusterAt(cid, nc);
541  }
542  else {
543  Int_t part_sub_id = particles_sub->GetEntriesFast();
544  AliEmcalParticle* part_sub = new ((*particles_sub)[part_sub_id]) AliEmcalParticle(c);
545  part_sub->SetPtEtaPhiM(constituents[ic].perp(),constituents[ic].eta(),constituents[ic].phi(),constituents[ic].m());
546  jet->AddTrackAt(part_sub_id, nt);
547  }
548 
549  ++nc;
550  ++nneutral;
551  }
552  else {
553  AliError(Form("%s: No logical way to end up here (uid = %d).", GetName(), uid));
554  continue;
555  }
556  }
557 
558  jet->SetNumberOfTracks(nt);
559  jet->SetNumberOfClusters(nc);
560  jet->SetNEF(neutralE / jet->E());
561  jet->SetMaxChargedPt(maxCh);
562  jet->SetMaxNeutralPt(maxNe);
563  if (gall > 0) jet->SetAreaEmc(jet->Area() * gemc / gall);
564  else jet->SetAreaEmc(-1);
565  jet->SetNEmc(cemc);
566  jet->SetNumberOfCharged(ncharged);
567  jet->SetNumberOfNeutrals(nneutral);
568  jet->SetMCPt(mcpt);
569  jet->SetPtEmc(emcpt);
570  jet->SortConstituents();
571 }
572 
579 Int_t AliEmcalJetTask::GetIndexSub(Double_t phi_sub, std::vector<fastjet::PseudoJet>& constituents_unsub)
580 {
581  Double_t dphi=0;
582  Double_t phimin=100;
583  Int_t index=0;
584  for (UInt_t ii = 0; ii < constituents_unsub.size(); ii++) {
585  dphi=0;
586  Double_t phi_unsub = constituents_unsub[ii].phi();
587  dphi=phi_unsub-phi_sub;
588  if (dphi < -1*TMath::Pi()) dphi += (2*TMath::Pi());
589  else if (dphi > TMath::Pi()) dphi -= (2*TMath::Pi());
590  if(TMath::Abs(dphi)<phimin){ phimin=TMath::Abs(dphi);
591  index=ii;} }
592  if (constituents_unsub[index].user_index()!=-1) return constituents_unsub[index].user_index();
593 
594  return 0;
595 }
596 
604 {
605  if (fLocked) {
606  AliFatal("Jet finder task is locked! Changing properties is not allowed.");
607  return kTRUE;
608  }
609  else {
610  return kFALSE;
611  }
612 }
613 
621 {
622  if(!fIsPSelSet) {
623  fIsPSelSet = kTRUE;
624  fOfflineTriggerMask = offlineTriggerMask;
625  }
626  else {
627  AliWarning("Manually setting the event selection for jet finders is not allowed! Using trigger=old_trigger | your_trigger");
628  fOfflineTriggerMask = fOfflineTriggerMask | offlineTriggerMask;
629  }
630 }
631 
638 {
639  if (IsLocked()) return;
640 
641  TIter nextPartColl(&fParticleCollArray);
642  AliParticleContainer* tracks = 0;
643  while ((tracks = static_cast<AliParticleContainer*>(nextPartColl()))) {
644  tracks->SetParticleEtaLimits(emi, ema);
645  }
646 }
647 
653 {
654  if (IsLocked()) return;
655 
656  TIter nextClusColl(&fClusterCollArray);
657  AliClusterContainer* clusters = 0;
658  while ((clusters = static_cast<AliClusterContainer*>(nextClusColl()))) {
659  clusters->SetClusPtCut(min);
660  }
661 }
662 
668 {
669  if (IsLocked()) return;
670 
671  TIter nextClusColl(&fClusterCollArray);
672  AliClusterContainer* clusters = 0;
673  while ((clusters = static_cast<AliClusterContainer*>(nextClusColl()))) {
674  clusters->SetClusECut(min);
675  }
676 }
677 
683 {
684  if (IsLocked()) return;
685 
686  TIter nextPartColl(&fParticleCollArray);
687  AliParticleContainer* tracks = 0;
688  while ((tracks = static_cast<AliParticleContainer*>(nextPartColl()))) {
689  tracks->SetParticlePtCut(min);
690  }
691 }
692 
699 {
700  if (IsLocked()) return;
701 
702  TIter nextPartColl(&fParticleCollArray);
703  AliParticleContainer* tracks = 0;
704  while ((tracks = static_cast<AliParticleContainer*>(nextPartColl()))) {
705  tracks->SetParticlePhiLimits(pmi, pma);
706  }
707 }
708 
715 fastjet::JetAlgorithm AliEmcalJetTask::ConvertToFJAlgo(EJetAlgo_t algo)
716 {
717  switch(algo) {
719  return fastjet::kt_algorithm;
721  return fastjet::antikt_algorithm;
723  return fastjet::cambridge_algorithm;
725  return fastjet::genkt_algorithm;
727  return fastjet::cambridge_for_passive_algorithm;
729  return fastjet::genkt_for_passive_algorithm;
731  return fastjet::plugin_algorithm;
733  return fastjet::undefined_jet_algorithm;
734 
735  default:
736  ::Error("AliEmcalJetTask::ConvertToFJAlgo", "Jet algorithm %d not recognized!!!", algo);
737  return fastjet::undefined_jet_algorithm;
738  }
739 }
740 
747 fastjet::RecombinationScheme AliEmcalJetTask::ConvertToFJRecoScheme(ERecoScheme_t reco)
748 {
749  switch(reco) {
751  return fastjet::E_scheme;
752 
754  return fastjet::pt_scheme;
755 
757  return fastjet::pt2_scheme;
758 
760  return fastjet::Et_scheme;
761 
763  return fastjet::Et2_scheme;
764 
766  return fastjet::BIpt_scheme;
767 
769  return fastjet::BIpt2_scheme;
770 
772  return fastjet::external_scheme;
773 
774  default:
775  ::Error("AliEmcalJetTask::ConvertToFJRecoScheme", "Recombination scheme %d not recognized!!!", reco);
776  return fastjet::external_scheme;
777  }
778 }
779 
785 
786  //This method has to be called after the run number is known because it needs the EMCal geometry object.
787 
788  UInt_t jetAcceptanceType = AliEmcalJet::kUser; // all jets satify the "no acceptance cut" condition
789 
790  // Check if TPC
791  if( eta < 0.9 && eta > -0.9 ) {
792  jetAcceptanceType |= AliEmcalJet::kTPC;
793  // Check if TPCfid
794  if (eta < 0.9 - r && eta > -0.9 + r)
795  jetAcceptanceType |= AliEmcalJet::kTPCfid;
796  }
797 
798  // Check if EMCAL
799  if( IsJetInEmcal(eta, phi, 0) ) {
800  jetAcceptanceType |= AliEmcalJet::kEMCAL;
801  // Check if EMCALfid
802  if( IsJetInEmcal(eta, phi, r) )
803  jetAcceptanceType |= AliEmcalJet::kEMCALfid;
804  }
805 
806  // Check if DCAL
807  if( IsJetInDcal(eta, phi, 0) ) {
808  jetAcceptanceType |= AliEmcalJet::kDCAL;
809  // Check if DCALfid
810  if( IsJetInDcal(eta, phi, r) )
811  jetAcceptanceType |= AliEmcalJet::kDCALfid;
812  }
813 
814  return jetAcceptanceType;
815 }
816 
821 {
822  if (!fGeom) return kFALSE;
823 
824  if ( eta < fGeom->GetArm1EtaMax() - r && eta > fGeom->GetArm1EtaMin() + r ) {
825  if(fRunNumber >= 177295 && fRunNumber <= 197470) {//small SM masked in 2012 and 2013
826  if ( phi < 3.135 - r && phi > 1.405 + r )
827  return kTRUE;
828  }
829  else {
830  if ( phi < fGeom->GetEMCALPhiMax() * TMath::DegToRad() - r && phi > fGeom->GetArm1PhiMin() * TMath::DegToRad() + r)
831  return kTRUE;
832  }
833  }
834 
835  return kFALSE;
836 }
837 
842 {
843  if (!fGeom) return kFALSE;
844  if (eta < fGeom->GetArm1EtaMax() - r && eta > fGeom->GetArm1EtaMin() + r ) {
845  if ( phi < fGeom->GetDCALPhiMax() * TMath::DegToRad() - r && phi > fGeom->GetDCALPhiMin() * TMath::DegToRad() + r)
846  return kTRUE;
847  }
848  return kFALSE;
849 }
void SetMaxNeutralPt(Double32_t t)
Definition: AliEmcalJet.h:182
fastjet::PseudoJet GetJetAreaVector(UInt_t idx) const
Double_t Area() const
Definition: AliEmcalJet.h:117
TObjArray fClusterCollArray
cluster collection array
void AddTrackAt(Int_t track, Int_t idx)
Definition: AliEmcalJet.h:195
void SetParticlePtCut(Double_t cut)
TClonesArray * fJets
=true ghost particles will be filled in AliEmcalJet obj
virtual Int_t Run()
double Double_t
Definition: External.C:58
void SetEtaRange(Double_t emi, Double_t ema)
Bool_t IsLocked() const
DCal acceptance.
Definition: AliEmcalJet.h:62
Double_t Eta() const
Definition: AliEmcalJet.h:108
const AliClusterIterableMomentumContainer accepted_momentum() const
Base task in the EMCAL framework.
bidirectional stl iterator over the EMCAL iterable container
EJetType_t fJetType
Bool_t IsJetInEmcal(Double_t eta, Double_t phi, Double_t r)
Double_t Phi() const
Definition: AliEmcalJet.h:104
void FillJetConstituents(AliEmcalJet *jet, std::vector< fastjet::PseudoJet > &constituents, std::vector< fastjet::PseudoJet > &constituents_sub, Int_t flag=0, TClonesArray *particles_sub=0)
void AddClusterAt(Int_t clus, Int_t idx)
Definition: AliEmcalJet.h:194
Bool_t fFillGhost
=true to enable FJ 2.x behavior
AliEmcalJetUtility * AddUtility(AliEmcalJetUtility *utility)
Declaration of class AliTLorentzVector.
virtual void ProcessJet(AliEmcalJet *jet, Int_t ij, AliFJWrapper &fjw)=0
void SetArea(Double_t a)
Definition: AliEmcalJet.h:176
void ExecuteUtilities(AliEmcalJet *jet, Int_t ij)
Full acceptance, i.e. no acceptance cut applied – left to user.
Definition: AliEmcalJet.h:64
TCanvas * c
Definition: TestFitELoss.C:172
void SetMaxRap(Double_t maxrap)
Definition: AliFJWrapper.h:100
virtual void Terminate(AliFJWrapper &fjw)=0
void SetAreaE(Double_t a)
Definition: AliEmcalJet.h:179
ERecoScheme_t fRecombScheme
void SetRecombScheme(const fastjet::RecombinationScheme &scheme)
Definition: AliFJWrapper.h:96
Double_t fTrackEfficiency
void SetMCPt(Double_t p)
Definition: AliEmcalJet.h:189
Double_t E() const
Definition: AliEmcalJet.h:106
virtual ~AliEmcalJetTask()
TRandom * gRandom
static FJRecoScheme ConvertToFJRecoScheme(ERecoScheme_t reco)
void SetLabel(Int_t l)
Definition: AliEmcalJet.h:175
General jet finder task implementing a wrapper for FastJet.
static TString GenerateJetName(EJetType_t jetType, EJetAlgo_t jetAlgo, ERecoScheme_t recoScheme, Double_t radius, AliParticleContainer *partCont, AliClusterContainer *clusCont, TString tag)
UShort_t GetNumberOfConstituents() const
Definition: AliEmcalJet.h:127
Container for particles within the EMCAL framework.
Bool_t fIsPSelSet
=true if already initialized
void SetR(Double_t r)
Definition: AliFJWrapper.h:101
TObjArray fParticleCollArray
particle/track collection array
AliParticleContainer * GetParticleContainer(Int_t i=0) const
TPC fiducial acceptance (each eta edge narrowed by jet R)
Definition: AliEmcalJet.h:59
const std::vector< fastjet::PseudoJet > & GetInclusiveJets() const
Definition: AliFJWrapper.h:33
void SetJetAcceptanceType(UInt_t type)
Definition: AliEmcalJet.h:252
UInt_t FindJetAcceptanceType(Double_t eta, Double_t phi, Double_t r)
int Int_t
Definition: External.C:63
Bool_t IsJetInDcal(Double_t eta, Double_t phi, Double_t r)
TPC acceptance.
Definition: AliEmcalJet.h:58
unsigned int UInt_t
Definition: External.C:33
float Float_t
Definition: External.C:68
Double_t GetJetArea(UInt_t idx) const
AliEMCALGeometry * fGeom
!emcal geometry
Double_t Phi_0_2pi() const
EMCal acceptance.
Definition: AliEmcalJet.h:60
void SetLegacyMode(Bool_t mode)
Definition: AliFJWrapper.h:111
virtual AliVParticle * GetParticle(Int_t i=-1) const
TObjArray * fUtilities
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
void SetAlgorithm(const fastjet::JetAlgorithm &algor)
Definition: AliFJWrapper.h:95
AliClusterContainer * GetClusterContainer(Int_t i=0) const
Double_t Phi_0_2pi() const
Definition: AliEmcalJet.h:116
void SetMaxChargedPt(Double32_t t)
Definition: AliEmcalJet.h:183
virtual void Clear(const Option_t *="")
void SetNEF(Double_t nef)
Definition: AliEmcalJet.h:184
Bool_t GetSortedArray(Int_t indexes[], std::vector< fastjet::PseudoJet > array) const
void AddGhost(const Double_t dPx, const Double_t dPy, const Double_t dPz, const Double_t dE)
void SetNEmc(Int_t n)
Definition: AliEmcalJet.h:190
AliVCluster * GetCluster(Int_t i) const
virtual void Prepare(AliFJWrapper &fjw)=0
Int_t fMinMCLabel
minimum MC label value for the tracks/clusters being considered MC particles
void SetJetTask(AliEmcalJetTask *jetTask)
void SetMinJetClusE(Double_t min)
Double_t Pt() const
Definition: AliEmcalJet.h:96
void SelectCollisionCandidates(UInt_t offlineTriggerMask=AliVEvent::kMB)
void SetGhostArea(Double_t gharea)
Definition: AliFJWrapper.h:99
EJetAlgo_t fJetAlgo
void SetNumberOfCharged(Int_t n)
Definition: AliEmcalJet.h:187
std::vector< fastjet::PseudoJet > GetJetConstituents(UInt_t idx) const
void SetAreaEta(Double_t a)
Definition: AliEmcalJet.h:177
static FJJetAlgo ConvertToFJAlgo(EJetAlgo_t algo)
void SetParticleEtaLimits(Double_t min, Double_t max)
const std::vector< fastjet::PseudoJet > & GetInputVectors() const
Definition: AliFJWrapper.h:31
void SetMinJetTrackPt(Double_t min)
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
const AliParticleIterableMomentumContainer accepted_momentum() const
Bool_t GetMomentum(TLorentzVector &mom, const AliVCluster *vc, Double_t mass) const
Represent a jet reconstructed using the EMCal jet framework.
Definition: AliEmcalJet.h:44
void SetNumberOfTracks(Int_t n)
Definition: AliEmcalJet.h:186
Bool_t fLegacyMode
=true if emcal particles are given as input (for clusters)
void SetClusPtCut(Double_t cut)
Int_t fRunNumber
!run number (triggering RunChanged()
virtual void AddInputVector(Double_t px, Double_t py, Double_t pz, Double_t E, Int_t index=-99999)
void SortConstituents()
void SetMinJetClusPt(Double_t min)
DCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:63
Int_t GetIndexSub(Double_t phi_sub, std::vector< fastjet::PseudoJet > &constituents_unsub)
bool Bool_t
Definition: External.C:53
void SetPhiRange(Double_t pmi, Double_t pma)
void SetClusECut(Double_t cut)
void SetAreaPhi(Double_t a)
Definition: AliEmcalJet.h:178
void SetAreaType(const fastjet::AreaType &atype)
Definition: AliFJWrapper.h:97
void SetParticlePhiLimits(Double_t min, Double_t max)
void SetPtEmc(Double_t pt)
Definition: AliEmcalJet.h:191
Container structure for EMCAL clusters.
EMCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:61
Container for jet within the EMCAL jet framework.
AliFJWrapper fFastJetWrapper
jet collection
void SetNumberOfClusters(Int_t n)
Definition: AliEmcalJet.h:185
virtual void Init()=0
const Double_t phimin
void SetAxisInEmcal(Bool_t b)
Definition: AliEmcalJet.h:181
static const Int_t fgkConstIndexShift
fastjet wrapper
void SetAreaEmc(Double_t a)
Definition: AliEmcalJet.h:180
void SetNumberOfNeutrals(Int_t n)
Definition: AliEmcalJet.h:188