AliPhysics  764b6ea (764b6ea)

Calculate the weight to the event to be applied when filling histograms. More...

#include <AliAnaWeights.h>

Inheritance diagram for AliAnaWeights:

Public Member Functions

 AliAnaWeights ()
 Constructor. More...
 
virtual ~AliAnaWeights ()
 Destructor. More...
 
TListGetCreateOutputHistograms ()
 
Int_t GetDebug () const
 
void SetDebug (Int_t d)
 
void InitCentralityWeightsHistogram (Int_t nbins=100, Int_t minCen=0, Int_t maxCen=100)
 
TH1F * GetCentralityWeightsHistogram ()
 
void SetCentralityWeightsHistogram (TH1F *h)
 
void SetCentrality (Float_t cen)
 
Bool_t IsCentralityWeightOn () const
 
void SwitchOnCentralityWeight ()
 
void SwitchOffCentralityWeight ()
 
virtual Double_t GetPythiaCrossSection ()
 
virtual Double_t GetWeight ()
 
Bool_t IsWeightSettingOn () const
 
Bool_t IsMCCrossSectionCalculationOn () const
 
Bool_t IsMCCrossSectionJustHistoFillOn () const
 
void SetPythiaEventHeader (AliGenPythiaEventHeader *py)
 
void SwitchOnMCCrossSectionCalculation ()
 
void SwitchOffMCCrossSectionCalculation ()
 
void SwitchOnMCCrossSectionHistoFill ()
 
void SwitchOnMCCrossSectionFromEventHeader ()
 
void SwitchOffMCCrossSectionFromEventHeader ()
 
Double_t GetParticlePtWeight (Float_t pt, Int_t pdg, TString genName, Int_t igen) const
 
void SetEtaFunction (TF1 *fun)
 
void SetPi0Function (TF1 *fun)
 
void SwitchOnMCParticlePtWeights ()
 
void SwitchOffMCParticlePtWeights ()
 

Static Public Member Functions

static Bool_t GetPythiaInfoFromFile (TString currFile, Float_t &xsec, Float_t &trials)
 

Private Member Functions

 AliAnaWeights (const AliAnaWeights &)
 Copy constructor not implemented. More...
 
AliAnaWeightsoperator= (const AliAnaWeights &)
 Assignment operator not implemented. More...
 

Private Attributes

Int_t fDebug
 Debug level. More...
 
TH1F * fhCentralityWeight
 Container of centrality weights. More...
 
Float_t fCentrality
 Container of centrality percentile. More...
 
Bool_t fUseCentralityWeight
 Return the centratlity weight. More...
 
Bool_t fDoMCParticlePtWeights
 activate the generation of a pT weight depending on MC particle pdg and generator More...
 
TF1 * fEtaFunction
 ! eta spectrum parametrization More...
 
TF1 * fPi0Function
 ! pi0 spectrum parametrization More...
 
Double_t fMCWeight
 pT-hard bin MC weight. It is used only internally. More...
 
TString fCurrFileName
 Current file path name. More...
 
Bool_t fCheckMCCrossSection
 Retrieve from the pyxsec.root file the cross section, only if requested. More...
 
Bool_t fJustFillCrossSecHist
 Do not provide a weight, just fill cross section histograms. More...
 
TH1F * fhXsec
 ! Cross section in PYTHIA. More...
 
TH1F * fhTrials
 ! Number of event trials in PYTHIA. More...
 
AliGenPythiaEventHeader * fPyEventHeader
 ! Pythia event header, needed to retrieve cross section, only in recent MC More...
 
Bool_t fCheckPythiaEventHeader
 Get cross section from pythia event header. More...
 

Detailed Description

Calculate the weight to the event to be applied when filling histograms.

This class is intended to get/calculate the global event weight needed to fill histograms. Right now it only gets the weight corresponding to:

  • the cross section of MC PYTHIA productions produced in pT-hard bins.
  • in case of centrality dependent weight.

This task is called by the class AliCaloTrackReader in method AliCaloTrackReader::FillInputEvent(). Also, in case of MC pT-hard bins analysis, the cross section and trials used are stored in dedicated histograms via AliAnaCaloTrackCorrMaker, being initialized in AliAnaCaloTrackCorrMaker::GetOutputContainer()

The different analysis subtasks should recover the weight stored in the reader, and apply it when filling the histogram.

Author
Gustavo Conesa Balbastre Gusta.nosp@m.vo.C.nosp@m.onesa.nosp@m..Bal.nosp@m.bastr.nosp@m.e@ce.nosp@m.rn.ch, LPSC-IN2P3-CNRS

Definition at line 34 of file AliAnaWeights.h.

Constructor & Destructor Documentation

AliAnaWeights::AliAnaWeights ( )

Constructor.

Definition at line 39 of file AliAnaWeights.cxx.

AliAnaWeights::~AliAnaWeights ( )
virtual

Destructor.

Definition at line 60 of file AliAnaWeights.cxx.

AliAnaWeights::AliAnaWeights ( const AliAnaWeights )
private

Copy constructor not implemented.

Member Function Documentation

TH1F * AliAnaWeights::GetCentralityWeightsHistogram ( )
Returns
histogram with centrality weights.

Definition at line 355 of file AliAnaWeights.cxx.

Referenced by SetDebug().

TList * AliAnaWeights::GetCreateOutputHistograms ( )

Init histogram pointers.

Returns
TList containing histograms.

Definition at line 73 of file AliAnaWeights.cxx.

Referenced by AliAnaCaloTrackCorrMaker::GetOutputContainer().

Int_t AliAnaWeights::GetDebug ( ) const
inline

Definition at line 48 of file AliAnaWeights.h.

Referenced by AliCaloTrackReader::Init().

Double_t AliAnaWeights::GetParticlePtWeight ( Float_t  pt,
Int_t  pdg,
TString  genName,
Int_t  igen 
) const
Returns
the particle pT weight depending on
Parameters
ptinput particle transverse momentum
pdginput particle type pdg code
genNameTString with generator name
igenindex of generator (not used yet).

currently, only pi0 and eta cases are considered. The parametrizations are passed via SetEtaFunction(TF1* fun) and SetPi0Function(TF1* fun) ex param: TF1* PowerEta0 = new TF1("PowerEta0","[0]*pow([1]/(([1]*exp(-[3]*x)+x)),[2])",4,25); Only active when SwitchOnMCParticlePtWeights()

Definition at line 140 of file AliAnaWeights.cxx.

Referenced by SwitchOffMCCrossSectionFromEventHeader().

Double_t AliAnaWeights::GetPythiaCrossSection ( )
virtual

Read the cross sections and number of trials from pyxsec.root (ESD) or pysec_hists.root (AODs), values stored in specific histograms-trees. If no file available, get information from Pythia event header

Fill the control histograms on number of trials and xsection optionally, with fJustFillCrossSecHist, just do that and do not return a weight. For cross section obtained from pythia event header, this is already the case.

Definition at line 169 of file AliAnaWeights.cxx.

Referenced by GetWeight(), and SwitchOffCentralityWeight().

Bool_t AliAnaWeights::GetPythiaInfoFromFile ( TString  file,
Float_t xsec,
Float_t trials 
)
static

This method gets and returns the

Parameters
file: either pyxsec.root (ESDs) or pysec_hists.root (AODs) files
xsec: cross section
trials: number of event trials that should be located where the main data file is. This is called in Notify and should provide the path to the AOD/ESD file

Definition at line 261 of file AliAnaWeights.cxx.

Referenced by GetPythiaCrossSection(), and SwitchOffCentralityWeight().

Double_t AliAnaWeights::GetWeight ( )
virtual
Returns
weight to be applied to the event. The weight can be right now:
  • pT-hard bin in PYTHIA MC events
  • centrality dependent weight

Definition at line 99 of file AliAnaWeights.cxx.

Referenced by AliCaloTrackReader::FillInputEvent(), and SwitchOffCentralityWeight().

void AliAnaWeights::InitCentralityWeightsHistogram ( Int_t  nbins = 100,
Int_t  minCen = 0,
Int_t  maxCen = 100 
)

This method intializes the histogram containing the centrality weights. Use in case of different binning than default.

Parameters
nbins: number of bins in histogram, default 100
minCen: minimum value of centrality, default 0
maxCen: maximum value of centrality, default 100

Definition at line 342 of file AliAnaWeights.cxx.

Referenced by GetCentralityWeightsHistogram(), SetDebug(), and SwitchOnCentralityWeight().

Bool_t AliAnaWeights::IsCentralityWeightOn ( ) const
inline

Definition at line 63 of file AliAnaWeights.h.

Referenced by AliAnaCaloTrackCorrMaker::FillControlHistograms().

Bool_t AliAnaWeights::IsMCCrossSectionCalculationOn ( ) const
inline
Bool_t AliAnaWeights::IsMCCrossSectionJustHistoFillOn ( ) const
inline
Bool_t AliAnaWeights::IsWeightSettingOn ( ) const
inline

Definition at line 78 of file AliAnaWeights.h.

Referenced by AliCaloTrackReader::FillInputEvent().

AliAnaWeights& AliAnaWeights::operator= ( const AliAnaWeights )
private

Assignment operator not implemented.

void AliAnaWeights::SetCentrality ( Float_t  cen)
inline

Definition at line 61 of file AliAnaWeights.h.

Referenced by AliCaloTrackReader::FillInputEvent().

void AliAnaWeights::SetCentralityWeightsHistogram ( TH1F *  h)
inline

Definition at line 59 of file AliAnaWeights.h.

void AliAnaWeights::SetDebug ( Int_t  d)
inline

Definition at line 50 of file AliAnaWeights.h.

void AliAnaWeights::SetEtaFunction ( TF1 *  fun)
inline

Definition at line 98 of file AliAnaWeights.h.

void AliAnaWeights::SetPi0Function ( TF1 *  fun)
inline

Definition at line 99 of file AliAnaWeights.h.

void AliAnaWeights::SetPythiaEventHeader ( AliGenPythiaEventHeader *  py)
inline

Definition at line 83 of file AliAnaWeights.h.

Referenced by AliCaloTrackReader::FillInputEvent().

void AliAnaWeights::SwitchOffCentralityWeight ( )
inline

Definition at line 67 of file AliAnaWeights.h.

void AliAnaWeights::SwitchOffMCCrossSectionCalculation ( )
inline

Definition at line 86 of file AliAnaWeights.h.

void AliAnaWeights::SwitchOffMCCrossSectionFromEventHeader ( )
inline

Definition at line 91 of file AliAnaWeights.h.

void AliAnaWeights::SwitchOffMCParticlePtWeights ( )
inline

Definition at line 102 of file AliAnaWeights.h.

void AliAnaWeights::SwitchOnCentralityWeight ( )
inline

Definition at line 65 of file AliAnaWeights.h.

void AliAnaWeights::SwitchOnMCCrossSectionCalculation ( )
inline

Definition at line 85 of file AliAnaWeights.h.

void AliAnaWeights::SwitchOnMCCrossSectionFromEventHeader ( )
inline

Definition at line 90 of file AliAnaWeights.h.

Referenced by AddTaskClusterShape().

void AliAnaWeights::SwitchOnMCCrossSectionHistoFill ( )
inline
void AliAnaWeights::SwitchOnMCParticlePtWeights ( )
inline

Definition at line 101 of file AliAnaWeights.h.

Member Data Documentation

Float_t AliAnaWeights::fCentrality
private

Container of centrality percentile.

Definition at line 114 of file AliAnaWeights.h.

Referenced by GetWeight(), and SetCentrality().

Bool_t AliAnaWeights::fCheckMCCrossSection
private
Bool_t AliAnaWeights::fCheckPythiaEventHeader
private

Get cross section from pythia event header.

Definition at line 146 of file AliAnaWeights.h.

Referenced by GetPythiaCrossSection(), SwitchOffMCCrossSectionFromEventHeader(), and SwitchOnMCCrossSectionFromEventHeader().

TString AliAnaWeights::fCurrFileName
private

Current file path name.

Definition at line 134 of file AliAnaWeights.h.

Referenced by GetPythiaCrossSection().

Int_t AliAnaWeights::fDebug
private

Debug level.

Definition at line 106 of file AliAnaWeights.h.

Referenced by GetDebug(), GetPythiaCrossSection(), and SetDebug().

Bool_t AliAnaWeights::fDoMCParticlePtWeights
private

activate the generation of a pT weight depending on MC particle pdg and generator

Definition at line 122 of file AliAnaWeights.h.

Referenced by GetParticlePtWeight(), SwitchOffMCParticlePtWeights(), and SwitchOnMCParticlePtWeights().

TF1* AliAnaWeights::fEtaFunction
private

! eta spectrum parametrization

Definition at line 124 of file AliAnaWeights.h.

Referenced by GetParticlePtWeight(), SetEtaFunction(), and ~AliAnaWeights().

TH1F* AliAnaWeights::fhCentralityWeight
private
TH1F* AliAnaWeights::fhTrials
private

! Number of event trials in PYTHIA.

Definition at line 142 of file AliAnaWeights.h.

Referenced by GetCreateOutputHistograms(), and GetPythiaCrossSection().

TH1F* AliAnaWeights::fhXsec
private

! Cross section in PYTHIA.

Definition at line 140 of file AliAnaWeights.h.

Referenced by GetCreateOutputHistograms(), and GetPythiaCrossSection().

Bool_t AliAnaWeights::fJustFillCrossSecHist
private

Do not provide a weight, just fill cross section histograms.

Definition at line 138 of file AliAnaWeights.h.

Referenced by GetPythiaCrossSection(), IsMCCrossSectionJustHistoFillOn(), and SwitchOnMCCrossSectionHistoFill().

Double_t AliAnaWeights::fMCWeight
private

pT-hard bin MC weight. It is used only internally.

Definition at line 132 of file AliAnaWeights.h.

Referenced by GetPythiaCrossSection().

TF1* AliAnaWeights::fPi0Function
private

! pi0 spectrum parametrization

Definition at line 126 of file AliAnaWeights.h.

Referenced by GetParticlePtWeight(), SetPi0Function(), and ~AliAnaWeights().

AliGenPythiaEventHeader* AliAnaWeights::fPyEventHeader
private

! Pythia event header, needed to retrieve cross section, only in recent MC

Definition at line 144 of file AliAnaWeights.h.

Referenced by GetPythiaCrossSection(), and SetPythiaEventHeader().

Bool_t AliAnaWeights::fUseCentralityWeight
private

Return the centratlity weight.

Definition at line 116 of file AliAnaWeights.h.

Referenced by GetWeight(), IsCentralityWeightOn(), IsWeightSettingOn(), SwitchOffCentralityWeight(), and SwitchOnCentralityWeight().


The documentation for this class was generated from the following files: