AliPhysics  bdbde52 (bdbde52)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliEmcalCorrectionClusterNonLinearity.cxx
Go to the documentation of this file.
1 // AliEmcalCorrectionClusterNonLinearity
2 //
3 
5 
6 #include <TList.h>
7 
8 #include "AliClusterContainer.h"
9 
13 
14 // Actually registers the class with the base class
16 
17 const std::map <std::string, AliEMCALRecoUtils::NonlinearityFunctions> AliEmcalCorrectionClusterNonLinearity::fgkNonlinearityFunctionMap = {
18  { "kPi0MC", AliEMCALRecoUtils::kPi0MC },
19  { "kPi0GammaGamma", AliEMCALRecoUtils::kPi0GammaGamma },
20  { "kPi0GammaConversion", AliEMCALRecoUtils::kPi0GammaConversion },
21  { "kNoCorrection", AliEMCALRecoUtils::kNoCorrection },
22  { "kBeamTest", AliEMCALRecoUtils::kBeamTest },
23  { "kBeamTestCorrected", AliEMCALRecoUtils::kBeamTestCorrected },
24  { "kPi0MCv2", AliEMCALRecoUtils::kPi0MCv2 },
25  { "kPi0MCv3", AliEMCALRecoUtils::kPi0MCv3 },
26  { "kBeamTestCorrectedv2", AliEMCALRecoUtils::kBeamTestCorrectedv2 },
27  { "kSDMv5", AliEMCALRecoUtils::kSDMv5 },
28  { "kPi0MCv5", AliEMCALRecoUtils::kPi0MCv5 },
29  { "kSDMv6", AliEMCALRecoUtils::kSDMv6 },
30  { "kPi0MCv6", AliEMCALRecoUtils::kPi0MCv6 },
31  { "kBeamTestCorrectedv3", AliEMCALRecoUtils::kBeamTestCorrectedv3 }
32 };
33 
38  AliEmcalCorrectionComponent("AliEmcalCorrectionClusterNonLinearity"),
39  fEnergyDistBefore(0),
40  fEnergyTimeHistBefore(0),
41  fEnergyDistAfter(0),
42  fEnergyTimeHistAfter(0)
43 
44 {
45 }
46 
51 {
52 }
53 
58 {
59  // Initialization
61 
62  GetProperty("createHistos", fCreateHisto);
63 
64  std::string nonLinFunctStr = "";
65  GetProperty("nonLinFunct", nonLinFunctStr);
66  UInt_t nonLinFunct = fgkNonlinearityFunctionMap.at(nonLinFunctStr);
67 
68  // init reco utils
69  if (!fRecoUtils)
70  fRecoUtils = new AliEMCALRecoUtils;
71  fRecoUtils->SetNonLinearityFunction(nonLinFunct);
72 
73  if (fRecoUtils) {
74  fRecoUtils->InitNonLinearityParam();
75  fRecoUtils->Print("");
76  }
77 
78  return kTRUE;
79 }
80 
85 {
87 
88  // Create my user objects.
89  if (fCreateHisto){
90  fEnergyDistBefore = new TH1F("hEnergyDistBefore","hEnergyDistBefore;E_{clus} (GeV)",1500,0,150);
92  fEnergyTimeHistBefore = new TH2F("hEnergyTimeDistBefore","hEnergyTimeDistBefore;E_{clus} (GeV);time",1500,0,150,500,0,1e-6);
94  fEnergyDistAfter = new TH1F("hEnergyDistAfter","hEnergyDistAfter;E_{clus} (GeV)",1500,0,150);
96  fEnergyTimeHistAfter = new TH2F("hEnergyTimeDistAfter","hEnergyTimeDistAfter;E_{clus} (GeV);time",1500,0,150,500,0,1e-6);
98 
99  // Take ownership of output list
100  fOutput->SetOwner(kTRUE);
101  }
102 }
103 
108 {
110 
111  // loop over clusters
112  AliVCluster *clus = 0;
113  AliClusterContainer * clusCont = 0;
114  TIter nextClusCont(&fClusterCollArray);
115  while ((clusCont = static_cast<AliClusterContainer*>(nextClusCont()))) {
116 
117  if (!clusCont) return kFALSE;
118  auto clusItCont = clusCont->all_momentum();
119 
120  for (AliClusterIterableMomentumContainer::iterator clusIterator = clusItCont.begin(); clusIterator != clusItCont.end(); ++clusIterator) {
121  clus = static_cast<AliVCluster *>(clusIterator->second);
122 
123  if (!clus->IsEMCAL()) continue;
124 
125  if (fCreateHisto) {
126  fEnergyDistBefore->Fill(clus->E());
127  fEnergyTimeHistBefore->Fill(clus->E(), clus->GetTOF());
128  }
129 
130  if (fRecoUtils) {
131  if (fRecoUtils->GetNonLinearityFunction() != AliEMCALRecoUtils::kNoCorrection) {
132  Double_t energy = fRecoUtils->CorrectClusterEnergyLinearity(clus);
133  clus->SetNonLinCorrEnergy(energy);
134  }
135  }
136 
137  // Fill histograms only if cluster is not exotic, as in ClusterMaker (the clusters are flagged, not removed)
138  if (fCreateHisto && !clus->GetIsExotic()) {
139  fEnergyDistAfter->Fill(clus->GetNonLinCorrEnergy());
140  fEnergyTimeHistAfter->Fill(clus->GetNonLinCorrEnergy(), clus->GetTOF());
141  }
142  }
143  }
144 
145  return kTRUE;
146 }
double Double_t
Definition: External.C:58
Definition: External.C:236
static const std::map< std::string, AliEMCALRecoUtils::NonlinearityFunctions > fgkNonlinearityFunctionMap
Relates string to the non-linearity function enumeration for YAML configuration.
bidirectional stl iterator over the EMCAL iterable container
energy
Definition: HFPtSpectrum.C:44
TH2F * fEnergyTimeHistBefore
!energy/time distribution before
AliEMCALRecoUtils * fRecoUtils
Pointer to RecoUtils.
TObjArray fClusterCollArray
Cluster collection array.
unsigned int UInt_t
Definition: External.C:33
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
Cluster energy non-linearity correction component in the EMCal correction framework.
bool Bool_t
Definition: External.C:53
Container structure for EMCAL clusters.
static RegisterCorrectionComponent< AliEmcalCorrectionClusterNonLinearity > reg
TH2F * fEnergyTimeHistAfter
!energy/time distribution after
bool GetProperty(std::string propertyName, T &property, bool requiredProperty=true, std::string correctionName="")
Retrieve property.