AliPhysics  8b695ca (8b695ca)
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 
49  GetProperty("fExoticMinCellAmplitude", fExoticMinCellAmplitude);
50  Float_t fMaxFcross = 0.97;
51  GetProperty("fMaxFcross", fMaxFcross);
53  GetProperty("fCellCrossMaxTimeDiff", fCellCrossMaxTimeDiff);
54 
55  // init reco utils
56  if (!fRecoUtils)
59  fRecoUtils->SetExoticCellMinAmplitudeCut(fExoticMinCellAmplitude);
61  fRecoUtils->SetExoticCellDiffTimeCut(fCellCrossMaxTimeDiff);
62  if (fRecoUtils)
63  fRecoUtils->Print("");
64 
65  return kTRUE;
66 }
67 
72 {
74 
75  // Create my user objects.
76  if (fCreateHisto){
77  fEtaPhiDistBefore = new TH2F("hEtaPhiDistBefore","hEtaPhiDistBefore;#eta;#phi",280,-0.7,0.7,200*3.14,0,2*3.14);
79  fEtaPhiDistAfter = new TH2F("hEtaPhiDistAfter","hEtaPhiDistAfter;#eta;#phi",280,-0.7,0.7,200*3.14,0,2*3.14);
81  fEnergyExoticClusters = new TH1F("fEnergyExoticClusters","fEnergyExoticClusters;E_{ex clus} (GeV)",1500,0,150);
83 
84  // Take ownership of output list
85  fOutput->SetOwner(kTRUE);
86  }
87 }
88 
93 {
95 
96  // loop over clusters
97  AliVCluster *clus = 0;
98  AliClusterContainer * clusCont = 0;
99  TIter nextClusCont(&fClusterCollArray);
100  while ((clusCont = static_cast<AliClusterContainer*>(nextClusCont()))) {
101 
102  if (!clusCont) continue;
103  auto clusItCont = clusCont->all_momentum();
104 
105  for (AliClusterIterableMomentumContainer::iterator clusIterator = clusItCont.begin(); clusIterator != clusItCont.end(); ++clusIterator) {
106  clus = static_cast<AliVCluster *>(clusIterator->second);
107 
108  if (!clus->IsEMCAL()) continue;
109 
110  if (fCreateHisto) {
111  Float_t pos[3] = {0.};
112  clus->GetPosition(pos);
113  TVector3 vec(pos);
114  // Phi needs to be in 0 to 2 Pi
115  fEtaPhiDistBefore->Fill(vec.Eta(), TVector2::Phi_0_2pi(vec.Phi()));
116  }
117 
118  Bool_t exResult = kFALSE;
119 
120  if (fRecoUtils) {
122  Bool_t exRemoval = fRecoUtils->IsRejectExoticCell();
123  fRecoUtils->SwitchOnRejectExoticCell(); //switch on temporarily
124  exResult = fRecoUtils->IsExoticCluster(clus, fCaloCells);
125  if (!exRemoval) fRecoUtils->SwitchOffRejectExoticCell(); //switch back off
126 
127  clus->SetIsExotic(exResult);
128  }
129  }
130 
131  if (fCreateHisto) {
132  if (exResult) {
133  fEnergyExoticClusters->Fill(clus->E());
134  }
135  else {
136  Float_t pos[3] = {0.};
137  clus->GetPosition(pos);
138  TVector3 vec(pos);
139  // Phi needs to be in 0 to 2 Pi
140  fEtaPhiDistAfter->Fill(vec.Eta(), TVector2::Phi_0_2pi(vec.Phi()));
141  }
142  }
143  }
144  }
145 
146  return kTRUE;
147 }
Exotic cluster removal in the EMCal correction framework.
Definition: External.C:236
Bool_t IsRejectExoticCluster() const
bidirectional stl iterator over the EMCAL iterable container
Bool_t IsExoticCluster(const AliVCluster *cluster, AliVCaloCells *cells, Int_t bc=0)
AliVCaloCells * fCaloCells
! Pointer to CaloCells
Some utilities for cluster and cell treatment.
TH2F * fEtaPhiDistAfter
!eta/phi distribution after
AliEMCALRecoUtils * fRecoUtils
Pointer to RecoUtils.
TObjArray fClusterCollArray
Cluster collection array.
void Print(const Option_t *) const
TH1F * fEnergyExoticClusters
!energy of exotic clusters
float Float_t
Definition: External.C:68
Base class for correction components in the EMCal correction framework.
void SetExoticCellDiffTimeCut(Float_t dt)
TList * fOutput
! List of output histograms
Bool_t fCreateHisto
Flag to make some basic histograms.
const AliClusterIterableMomentumContainer all_momentum() const
void SwitchOnRejectExoticCluster()
Bool_t IsRejectExoticCell() const
void SetExoticCellMinAmplitudeCut(Float_t ma)
static RegisterCorrectionComponent< AliEmcalCorrectionClusterExotics > reg
void SetExoticCellFractionCut(Float_t f)
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.