AliPhysics  095eea3 (095eea3)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliEmcalCorrectionClusterExotics.cxx
Go to the documentation of this file.
1 // AliEmcalCorrectionClusterExotics
2 //
3 
5 
6 #include <TH2F.h>
7 #include <TList.h>
8 
9 #include "AliClusterContainer.h"
10 #include "AliEMCALRecoUtils.h"
11 
15 
16 // Actually registers the class with the base class
18 
23  AliEmcalCorrectionComponent("AliEmcalCorrectionClusterExotics"),
24  fEtaPhiDistBefore(0),
25  fEtaPhiDistAfter(0),
26  fEnergyExoticClusters(0)
27 {
28 }
29 
34 {
35 }
36 
41 {
42  // Initialization
44 
45  GetProperty("createHistos", fCreateHisto);
46 
47  // init reco utils
48  if (!fRecoUtils)
49  fRecoUtils = new AliEMCALRecoUtils;
50  fRecoUtils->SwitchOnRejectExoticCluster();
51  if (fRecoUtils)
52  fRecoUtils->Print("");
53 
54  return kTRUE;
55 }
56 
61 {
63 
64  // Create my user objects.
65  if (fCreateHisto){
66  fEtaPhiDistBefore = new TH2F("hEtaPhiDistBefore","hEtaPhiDistBefore;#eta;#phi",280,-0.7,0.7,200*3.14,0,2*3.14);
68  fEtaPhiDistAfter = new TH2F("hEtaPhiDistAfter","hEtaPhiDistAfter;#eta;#phi",280,-0.7,0.7,200*3.14,0,2*3.14);
70  fEnergyExoticClusters = new TH1F("fEnergyExoticClusters","fEnergyExoticClusters;E_{ex clus} (GeV)",1500,0,150);
72 
73  // Take ownership of output list
74  fOutput->SetOwner(kTRUE);
75  }
76 }
77 
82 {
84 
85  // loop over clusters
86  AliVCluster *clus = 0;
87  AliClusterContainer * clusCont = 0;
88  TIter nextClusCont(&fClusterCollArray);
89  while ((clusCont = static_cast<AliClusterContainer*>(nextClusCont()))) {
90 
91  if (!clusCont) continue;
92  auto clusItCont = clusCont->all_momentum();
93 
94  for (AliClusterIterableMomentumContainer::iterator clusIterator = clusItCont.begin(); clusIterator != clusItCont.end(); ++clusIterator) {
95  clus = static_cast<AliVCluster *>(clusIterator->second);
96 
97  if (!clus->IsEMCAL()) continue;
98 
99  if (fCreateHisto) {
100  Float_t pos[3] = {0.};
101  clus->GetPosition(pos);
102  TVector3 vec(pos);
103  // Phi needs to be in 0 to 2 Pi
104  fEtaPhiDistBefore->Fill(vec.Eta(), TVector2::Phi_0_2pi(vec.Phi()));
105  }
106 
107  Bool_t exResult = kFALSE;
108 
109  if (fRecoUtils) {
110  if (fRecoUtils->IsRejectExoticCluster()) {
111  Bool_t exRemoval = fRecoUtils->IsRejectExoticCell();
112  fRecoUtils->SwitchOnRejectExoticCell(); //switch on temporarily
113  exResult = fRecoUtils->IsExoticCluster(clus, fCaloCells);
114  if (!exRemoval) fRecoUtils->SwitchOffRejectExoticCell(); //switch back off
115 
116  clus->SetIsExotic(exResult);
117  }
118  }
119 
120  if (fCreateHisto) {
121  if (exResult) {
122  fEnergyExoticClusters->Fill(clus->E());
123  }
124  else {
125  Float_t pos[3] = {0.};
126  clus->GetPosition(pos);
127  TVector3 vec(pos);
128  // Phi needs to be in 0 to 2 Pi
129  fEtaPhiDistAfter->Fill(vec.Eta(), TVector2::Phi_0_2pi(vec.Phi()));
130  }
131  }
132  }
133  }
134 
135  return kTRUE;
136 }
Exotic cluster removal in the EMCal correction framework.
Definition: External.C:236
bidirectional stl iterator over the EMCAL iterable container
AliVCaloCells * fCaloCells
! Pointer to CaloCells
TH2F * fEtaPhiDistAfter
!eta/phi distribution after
AliEMCALRecoUtils * fRecoUtils
Pointer to RecoUtils.
TObjArray fClusterCollArray
Cluster collection array.
TH1F * fEnergyExoticClusters
!energy of exotic clusters
float Float_t
Definition: External.C:68
Base class for correction components in the EMCal correction framework.
TList * fOutput
! List of output histograms
Bool_t fCreateHisto
Flag to make some basic histograms.
const AliClusterIterableMomentumContainer all_momentum() const
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
static RegisterCorrectionComponent< AliEmcalCorrectionClusterExotics > reg
bool Bool_t
Definition: External.C:53
TH2F * fEtaPhiDistBefore
!eta/phi distribution before
Container structure for EMCAL clusters.
bool GetProperty(std::string propertyName, T &property, bool requiredProperty=true, std::string correctionName="")
Retrieve property.