AliPhysics  master (3d17d9d)
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  { "kPCMv1", AliEMCALRecoUtils::kPCMv1 },
33  { "kPCMplusBTCv1", AliEMCALRecoUtils::kPCMplusBTCv1 },
34  { "kPCMsysv1", AliEMCALRecoUtils::kPCMsysv1 },
35  { "kBeamTestCorrectedv4", AliEMCALRecoUtils::kBeamTestCorrectedv4 },
36  { "kBeamTestNS", AliEMCALRecoUtils::kBeamTestNS },
37  { "kPi0MCNS", AliEMCALRecoUtils::kPi0MCNS }
38 };
39 
44  AliEmcalCorrectionComponent("AliEmcalCorrectionClusterNonLinearity"),
45  fEnergyDistBefore(0),
46  fEnergyTimeHistBefore(0),
47  fEnergyDistAfter(0),
48  fEnergyTimeHistAfter(0),
49  fSetForceClusterE(kFALSE)
50 {
51 }
52 
57 {
58 }
59 
64 {
65  // Initialization
67 
68  std::string nonLinFunctStr = "";
69  GetProperty("nonLinFunct", nonLinFunctStr);
70  UInt_t nonLinFunct = fgkNonlinearityFunctionMap.at(nonLinFunctStr);
71 
72  GetProperty("setForceClusterE", fSetForceClusterE);
73 
74  // init reco utils
75  if (!fRecoUtils)
77  fRecoUtils->SetNonLinearityFunction(nonLinFunct);
78 
79  if (fRecoUtils) {
81  fRecoUtils->Print("");
82  }
83 
84  return kTRUE;
85 }
86 
91 {
93 
94  // Create my user objects.
95  if (fCreateHisto){
96  fEnergyDistBefore = new TH1F("hEnergyDistBefore","hEnergyDistBefore;E_{clus} (GeV)",1500,0,150);
98  fEnergyTimeHistBefore = new TH2F("hEnergyTimeDistBefore","hEnergyTimeDistBefore;E_{clus} (GeV);time (s)",1500,0,150,500,-1e-6,1e-6);
100  fEnergyDistAfter = new TH1F("hEnergyDistAfter","hEnergyDistAfter;E_{clus} (GeV)",1500,0,150);
102  fEnergyTimeHistAfter = new TH2F("hEnergyTimeDistAfter","hEnergyTimeDistAfter;E_{clus} (GeV);time (s)",1500,0,150,500,-1e-6,1e-6);
104 
105  // Take ownership of output list
106  fOutput->SetOwner(kTRUE);
107  }
108 }
109 
114 {
116 
117  // loop over clusters
118  AliVCluster *clus = 0;
119  AliClusterContainer * clusCont = 0;
120  TIter nextClusCont(&fClusterCollArray);
121  while ((clusCont = static_cast<AliClusterContainer*>(nextClusCont()))) {
122 
123  if (!clusCont) return kFALSE;
124  auto clusItCont = clusCont->all_momentum();
125 
126  for (AliClusterIterableMomentumContainer::iterator clusIterator = clusItCont.begin(); clusIterator != clusItCont.end(); ++clusIterator) {
127  clus = static_cast<AliVCluster *>(clusIterator->second);
128 
129  if (!clus->IsEMCAL()) continue;
130 
131  if (fCreateHisto) {
132  fEnergyDistBefore->Fill(clus->E());
133  fEnergyTimeHistBefore->Fill(clus->E(), clus->GetTOF());
134  }
135 
136  if (fRecoUtils) {
139  clus->SetNonLinCorrEnergy(energy);
140  if ( fSetForceClusterE ) clus->SetE(energy);
141  }
142  }
143 
144  // Fill histograms only if cluster is not exotic, as in ClusterMaker (the clusters are flagged, not removed)
145  if (fCreateHisto && !clus->GetIsExotic()) {
146  Float_t energy = clus->GetNonLinCorrEnergy();
147  if(fSetForceClusterE) energy = clus->E();
148 
149  fEnergyDistAfter->Fill(energy);
150  fEnergyTimeHistAfter->Fill(energy, clus->GetTOF());
151  }
152  }
153  }
154 
155  return kTRUE;
156 }
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:45
void SetNonLinearityFunction(Int_t fun)
TH2F * fEnergyTimeHistBefore
!energy/time distribution before
Bool_t fSetForceClusterE
Only for backwards compatibility, force cluster->E() to be set to the cluster non-linearity corrected...
Some utilities for cluster and cell treatment.
AliEMCALRecoUtils * fRecoUtils
Pointer to RecoUtils.
TObjArray fClusterCollArray
Cluster collection array.
void Print(const Option_t *) const
unsigned int UInt_t
Definition: External.C:33
float Float_t
Definition: External.C:68
Base class for correction components in the EMCal correction framework.
Int_t GetNonLinearityFunction() const
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
Float_t CorrectClusterEnergyLinearity(AliVCluster *clu)
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.