AliPhysics  6133e27 (6133e27)
 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  fExoticMinCellAmplitude(0),
28  fMaxFcross(0),
29  fCellCrossMaxTimeDiff(0)
30 {
31 }
32 
37 {
38 }
39 
44 {
45  // Initialization
47 
48  GetProperty("createHistos", fCreateHisto);
49 
51  GetProperty("fExoticMinCellAmplitude", fExoticMinCellAmplitude);
52  Float_t fMaxFcross = 0.97;
53  GetProperty("fMaxFcross", fMaxFcross);
55  GetProperty("fCellCrossMaxTimeDiff", fCellCrossMaxTimeDiff);
56 
57  // init reco utils
58  if (!fRecoUtils)
59  fRecoUtils = new AliEMCALRecoUtils;
60  fRecoUtils->SwitchOnRejectExoticCluster();
61  fRecoUtils->SetExoticCellMinAmplitudeCut(fExoticMinCellAmplitude);
62  fRecoUtils->SetExoticCellFractionCut(fMaxFcross);
63  fRecoUtils->SetExoticCellDiffTimeCut(fCellCrossMaxTimeDiff);
64  if (fRecoUtils)
65  fRecoUtils->Print("");
66 
67  return kTRUE;
68 }
69 
74 {
76 
77  // Create my user objects.
78  if (fCreateHisto){
79  fEtaPhiDistBefore = new TH2F("hEtaPhiDistBefore","hEtaPhiDistBefore;#eta;#phi",280,-0.7,0.7,200*3.14,0,2*3.14);
81  fEtaPhiDistAfter = new TH2F("hEtaPhiDistAfter","hEtaPhiDistAfter;#eta;#phi",280,-0.7,0.7,200*3.14,0,2*3.14);
83  fEnergyExoticClusters = new TH1F("fEnergyExoticClusters","fEnergyExoticClusters;E_{ex clus} (GeV)",1500,0,150);
85 
86  // Take ownership of output list
87  fOutput->SetOwner(kTRUE);
88  }
89 }
90 
95 {
97 
98  // loop over clusters
99  AliVCluster *clus = 0;
100  AliClusterContainer * clusCont = 0;
101  TIter nextClusCont(&fClusterCollArray);
102  while ((clusCont = static_cast<AliClusterContainer*>(nextClusCont()))) {
103 
104  if (!clusCont) continue;
105  auto clusItCont = clusCont->all_momentum();
106 
107  for (AliClusterIterableMomentumContainer::iterator clusIterator = clusItCont.begin(); clusIterator != clusItCont.end(); ++clusIterator) {
108  clus = static_cast<AliVCluster *>(clusIterator->second);
109 
110  if (!clus->IsEMCAL()) continue;
111 
112  if (fCreateHisto) {
113  Float_t pos[3] = {0.};
114  clus->GetPosition(pos);
115  TVector3 vec(pos);
116  // Phi needs to be in 0 to 2 Pi
117  fEtaPhiDistBefore->Fill(vec.Eta(), TVector2::Phi_0_2pi(vec.Phi()));
118  }
119 
120  Bool_t exResult = kFALSE;
121 
122  if (fRecoUtils) {
123  if (fRecoUtils->IsRejectExoticCluster()) {
124  Bool_t exRemoval = fRecoUtils->IsRejectExoticCell();
125  fRecoUtils->SwitchOnRejectExoticCell(); //switch on temporarily
126  exResult = fRecoUtils->IsExoticCluster(clus, fCaloCells);
127  if (!exRemoval) fRecoUtils->SwitchOffRejectExoticCell(); //switch back off
128 
129  clus->SetIsExotic(exResult);
130  }
131  }
132 
133  if (fCreateHisto) {
134  if (exResult) {
135  fEnergyExoticClusters->Fill(clus->E());
136  }
137  else {
138  Float_t pos[3] = {0.};
139  clus->GetPosition(pos);
140  TVector3 vec(pos);
141  // Phi needs to be in 0 to 2 Pi
142  fEtaPhiDistAfter->Fill(vec.Eta(), TVector2::Phi_0_2pi(vec.Phi()));
143  }
144  }
145  }
146  }
147 
148  return kTRUE;
149 }
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
ClassImp(AliAnalysisTaskDeltaPt) AliAnalysisTaskDeltaPt
Bool_t fCreateHisto
Flag to make some basic histograms.
const AliClusterIterableMomentumContainer all_momentum() const
static RegisterCorrectionComponent< AliEmcalCorrectionClusterExotics > reg
bool Bool_t
Definition: External.C:53
TH2F * fEtaPhiDistBefore
!eta/phi distribution before
Float_t fCellCrossMaxTimeDiff
Max time difference allowed between leading cell and cross cells (in ns)
Container structure for EMCAL clusters.
Float_t fMaxFcross
Max value of Fcross = 1-Ecross/ecell allowed for clusters to pass exotic cut.
Float_t fExoticMinCellAmplitude
Min energy of leading cell in order for exotic cut to be attempted.
bool GetProperty(std::string propertyName, T &property, bool requiredProperty=true, std::string correctionName="")
Retrieve property.