AliPhysics  6cf2591 (6cf2591)
 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 //________________________________________________________________________
19  AliEmcalCorrectionComponent("AliEmcalCorrectionClusterNonLinearity"),
20  fEnergyDistBefore(0),
21  fEnergyTimeHistBefore(0),
22  fEnergyDistAfter(0),
23  fEnergyTimeHistAfter(0)
24 
25 {
26  // Default constructor
27  AliDebug(3, Form("%s", __PRETTY_FUNCTION__));
28 
29 }
30 
31 //________________________________________________________________________
33 {
34  // Destructor
35 }
36 
37 //________________________________________________________________________
39 {
40  // Initialization
41  AliDebug(3, Form("%s", __PRETTY_FUNCTION__));
43  // Do base class initializations and if it fails -> bail out
44  //AliAnalysisTaskEmcal::ExecOnce();
45  //if (!fInitialized) return;
46 
47  GetProperty("createHistos", fCreateHisto);
48 
49  std::string nonLinFunctStr = "";
50  GetProperty("nonLinFunct", nonLinFunctStr);
51  UInt_t nonLinFunct = nonlinearityFunctionMap[nonLinFunctStr];
52 
53  // init reco utils
54  if (!fRecoUtils)
55  fRecoUtils = new AliEMCALRecoUtils;
56  fRecoUtils->SetNonLinearityFunction(nonLinFunct);
57 
58  if (fRecoUtils) {
59  fRecoUtils->InitNonLinearityParam();
60  fRecoUtils->Print("");
61  }
62 
63  return kTRUE;
64 }
65 
66 //________________________________________________________________________
68 {
69  AliDebug(3, Form("%s", __PRETTY_FUNCTION__));
71 
72  // Create my user objects.
73  if (fCreateHisto){
74  fEnergyDistBefore = new TH1F("hEnergyDistBefore","hEnergyDistBefore;E_{clus} (GeV)",1500,0,150);
76  fEnergyTimeHistBefore = new TH2F("hEnergyTimeDistBefore","hEnergyTimeDistBefore;E_{clus} (GeV);time",1500,0,150,500,0,1e-6);
78  fEnergyDistAfter = new TH1F("hEnergyDistAfter","hEnergyDistAfter;E_{clus} (GeV)",1500,0,150);
80  fEnergyTimeHistAfter = new TH2F("hEnergyTimeDistAfter","hEnergyTimeDistAfter;E_{clus} (GeV);time",1500,0,150,500,0,1e-6);
82 
83  // Take ownership of output list
84  fOutput->SetOwner(kTRUE);
85  }
86 }
87 
88 //________________________________________________________________________
90 {
91  // Run
92  AliDebug(3, Form("%s", __PRETTY_FUNCTION__));
94 
95  if (!fClusCont) return kFALSE;
96 
97  // loop over clusters
98  fClusCont->ResetCurrentID();
99  AliVCluster *clus = 0;
100  while ((clus = fClusCont->GetNextCluster())) {
101  if (!clus->IsEMCAL()) continue;
102 
103  if (fCreateHisto) {
104  fEnergyDistBefore->Fill(clus->E());
105  fEnergyTimeHistBefore->Fill(clus->E(), clus->GetTOF());
106  }
107 
108  if (fRecoUtils) {
109  if (fRecoUtils->GetNonLinearityFunction() != AliEMCALRecoUtils::kNoCorrection) {
110  Double_t energy = fRecoUtils->CorrectClusterEnergyLinearity(clus);
111  clus->SetNonLinCorrEnergy(energy);
112  }
113  }
114 
115  // Fill histograms only if cluster is not exotic, as in ClusterMaker (the clusters are flagged, not removed)
116  if (fCreateHisto && !clus->GetIsExotic()) {
117  fEnergyDistAfter->Fill(clus->GetNonLinCorrEnergy());
118  fEnergyTimeHistAfter->Fill(clus->GetNonLinCorrEnergy(), clus->GetTOF());
119  }
120  }
121 
122  return kTRUE;
123 }
double Double_t
Definition: External.C:58
Definition: External.C:236
energy
Definition: HFPtSpectrum.C:44
TH2F * fEnergyTimeHistBefore
!energy/time distribution before
AliEMCALRecoUtils * fRecoUtils
Pointer to RecoUtils.
AliClusterContainer * fClusCont
Pointer to the cluster container.
unsigned int UInt_t
Definition: External.C:33
Base class for correction components in the EMCal correction framework.
AliVCluster * GetNextCluster()
TList * fOutput
! List of output histograms
Bool_t fCreateHisto
Flag to make some basic histograms.
Cluster energy non-linearity correction component in the EMCal correction framework.
ClassImp(AliAnalysisTaskCRC) AliAnalysisTaskCRC
std::map< std::string, AliEMCALRecoUtils::NonlinearityFunctions > nonlinearityFunctionMap
bool Bool_t
Definition: External.C:53
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.