AliPhysics  80ccde44 (80ccde44)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliAnalysisTaskDmesonJets::AnalysisEngine Class Reference

Struct that encapsulates analysis parameters. More...

#include <AliAnalysisTaskDmesonJets.h>

Inheritance diagram for AliAnalysisTaskDmesonJets::AnalysisEngine:

Public Types

enum  EFindParticleOriginMode_t { kFindFirst, kFindLast }
 
typedef std::pair< AliJetInfo
*, Double_t
jet_distance_pair
 

Public Member Functions

 AnalysisEngine ()
 This is the default constructor, used for ROOT I/O purposes. More...
 
 AnalysisEngine (ECandidateType_t type, EMCMode_t MCmode, AliRDHFCuts *cuts=0, Int_t nBins=80, Double_t range=0.50)
 
 AnalysisEngine (const AnalysisEngine &source)
 
AnalysisEngineoperator= (const AnalysisEngine &source)
 
virtual ~AnalysisEngine ()
 
void SetCandidateType (ECandidateType_t t)
 
void SetMCMode (EMCMode_t m)
 
void SetNMassBins (Int_t n)
 
void SetMassRange (Double_t min, Double_t max)
 
void AdoptRDHFCuts (AliRDHFCuts *cuts)
 
void SetRDHFCuts (AliRDHFCuts *cuts)
 
void SetRejectedOriginMap (UInt_t m)
 
void SetAcceptedDecayMap (UInt_t m)
 
void SetRejectISR (Bool_t b)
 
const char * GetCandidateName () const
 
const char * GetName () const
 
const char * GetName (const AliHFJetDefinition &jetDef) const
 
EMCMode_t GetMCMode () const
 
ECandidateType_t GetCandidateType () const
 
AliHFJetDefinitionAddJetDefinition (EJetType_t type, Double_t r, EJetAlgo_t algo, ERecoScheme_t reco)
 
AliHFJetDefinitionAddJetDefinition (const AliHFJetDefinition &def)
 
std::vector
< AliHFJetDefinition >
::iterator 
FindJetDefinition (const AliHFJetDefinition &eng)
 
std::vector
< AliAnalysisTaskDmesonJets::AliHFJetDefinition > & 
GetJetDefinitions ()
 
Bool_t IsAnyJetInAcceptance (const AliDmesonJetInfo &dMesonJet) const
 
std::map< int, AliDmesonJetInfo > & GetDmesons ()
 
void Init (const AliEMCALGeometry *const geom, Int_t runNumber)
 Initialize the analysis engine. More...
 
TTreeBuildTree (const char *taskName)
 
TTreeGetTree () const
 
Bool_t FillTree (Bool_t applyKinCuts)
 
void SetTrackEfficiency (Double_t t)
 
void AssignDataSlot (Int_t n)
 
Int_t GetDataSlotNumber () const
 
void BuildHnSparse (UInt_t enabledAxis)
 
Bool_t FillHnSparse (Bool_t applyKinCuts)
 
Bool_t FillHnSparse (THnSparse *h, const AliDmesonJetInfo &DmesonJet, std::string n)
 
Bool_t FillQA (Bool_t applyKinCuts)
 
Bool_t IsInhibit () const
 

Static Public Member Functions

static AliAODMCParticle * FindParticleOrigin (const AliAODMCParticle *part, TClonesArray *mcArray, EFindParticleOriginMode_t mode, const std::set< UInt_t > &pdgSet)
 
static AliAODMCParticle * FindParticleOrigin (const AliAODMCParticle *part, TClonesArray *mcArray, EFindParticleOriginMode_t mode)
 
static std::pair
< AliAnalysisTaskDmesonJets::EMesonOrigin_t,
AliAODMCParticle * > 
IsPromptCharm (const AliAODMCParticle *part, TClonesArray *mcArray)
 
static EMesonDecayChannel_t CheckDecayChannel (const AliAODMCParticle *part, TClonesArray *mcArray)
 

Public Attributes

std::map< AliAODMCParticle
*, Short_t
fPartons
 ! set of the partons in the shower that produced each D meson More...
 

Protected Member Functions

void RunAnalysis ()
 Run the analysis. More...
 

Protected Attributes

ECandidateType_t fCandidateType
 Candidate type. More...
 
TString fCandidateName
 Candidate name. More...
 
UInt_t fCandidatePDG
 Candidate PDG. More...
 
UChar_t fNDaughters
 Number of daughters. More...
 
TArrayI fPDGdaughters
 List of the PDG code of the daughters. More...
 
TString fBranchName
 AOD branch where the D meson candidate are found. More...
 
EMCMode_t fMCMode
 MC mode: No MC (data and MC detector level), background-only (MC), signal-only (MC), MC truth (particle level) More...
 
Int_t fNMassBins
 Mass number of bins. More...
 
Double_t fMinMass
 Min mass in histogram axis. More...
 
Double_t fMaxMass
 Max mass in histogram axis. More...
 
AliRDHFCutsfRDHFCuts
 D meson candidates cuts. More...
 
UInt_t fRejectedOrigin
 Bit mask with D meson origins that are rejected (used for MC analysis, i.e. signal-only, background-only and particle-level) More...
 
UInt_t fAcceptedDecay
 Bit mask with D meson decays that are accepted (only used for particle-level analysis) More...
 
Bool_t fInhibit
 Inhibit the task. More...
 
std::vector< AliHFJetDefinitionfJetDefinitions
 Jet definitions. More...
 
Float_t fPtBinWidth
 Histogram pt bin width. More...
 
Float_t fMaxPt
 Histogram pt limit. More...
 
TRandom * fRandomGen
 ! Random number generator More...
 
Double_t fTrackEfficiency
 ! Artificial tracking inefficiency (0...1) -> set automatically at ExecOnce by AliAnalysisTaskDmesonJets More...
 
Bool_t fRejectISR
 ! Reject initial state radiation More...
 
Int_t fDataSlotNumber
 ! Data slot where the tree output is posted More...
 
TTreefTree
 ! Output tree More...
 
AliDmesonInfoSummaryfCurrentDmesonJetInfo
 ! Current D meson jet info More...
 
AliJetInfoSummary ** fCurrentJetInfo
 ! Current jet info More...
 
std::map< int, AliDmesonJetInfofDmesonJets
 ! Array containing the D meson jets More...
 
TClonesArray * fCandidateArray
 ! D meson candidate array More...
 
AliHFAODMCParticleContainerfMCContainer
 ! MC particle container More...
 
std::vector< AliTrackContainer * > fTrackContainers
 ! Track containers More...
 
std::vector
< AliClusterContainer * > 
fClusterContainers
 ! Cluster containers More...
 
AliAODEventfAodEvent
 ! AOD event More...
 
AliFJWrapperfFastJetWrapper
 ! Fastjet wrapper More...
 
THistManagerfHistManager
 ! Histograms More...
 
Double_t fCent
 ! Event centrality More...
 

Private Member Functions

void AddInputVectors (AliEmcalContainer *cont, Int_t offset, TH2 *rejectHist=0, Double_t eff=0.)
 
void SetCandidateProperties (Double_t range)
 
AliAODMCParticle * MatchToMC () const
 
void RunDetectorLevelAnalysis ()
 Run a detector level analysis. More...
 
void RunParticleLevelAnalysis ()
 Run a particle level analysis. More...
 
Bool_t ExtractRecoDecayAttributes (const AliAODRecoDecayHF2Prong *Dcand, AliDmesonJetInfo &DmesonJet, UInt_t i)
 
Bool_t ExtractD0Attributes (const AliAODRecoDecayHF2Prong *Dcand, AliDmesonJetInfo &DmesonJet, UInt_t i)
 
Bool_t ExtractDstarAttributes (const AliAODRecoCascadeHF *DstarCand, AliDmesonJetInfo &DmesonJet, UInt_t i)
 
Bool_t FindJet (AliAODRecoDecayHF2Prong *Dcand, AliDmesonJetInfo &DmesonJet, AliHFJetDefinition &jetDef)
 

Friends

class AliAnalysisTaskDmesonJets
 
bool operator< (const AnalysisEngine &lhs, const AnalysisEngine &rhs)
 
bool operator> (const AnalysisEngine &lhs, const AnalysisEngine &rhs)
 
bool operator<= (const AnalysisEngine &lhs, const AnalysisEngine &rhs)
 
bool operator>= (const AnalysisEngine &lhs, const AnalysisEngine &rhs)
 
bool operator== (const AnalysisEngine &lhs, const AnalysisEngine &rhs)
 
bool operator!= (const AnalysisEngine &lhs, const AnalysisEngine &rhs)
 

Detailed Description

Struct that encapsulates analysis parameters.

This struct encapsulates analysis parameters for the D meson jet analysis.

Definition at line 404 of file AliAnalysisTaskDmesonJets.h.

Member Typedef Documentation

Member Enumeration Documentation

Enumerator
kFindFirst 
kFindLast 

Look for the very first particle in the fragmentation tree.

Look for the last particle in the fragmentation tree (closest to the hadron)

Definition at line 408 of file AliAnalysisTaskDmesonJets.h.

Constructor & Destructor Documentation

AliAnalysisTaskDmesonJets::AnalysisEngine::AnalysisEngine ( )

This is the default constructor, used for ROOT I/O purposes.

Definition at line 760 of file AliAnalysisTaskDmesonJets.cxx.

AliAnalysisTaskDmesonJets::AnalysisEngine::AnalysisEngine ( ECandidateType_t  type,
EMCMode_t  MCmode,
AliRDHFCuts cuts = 0,
Int_t  nMassBins = 80,
Double_t  range = 0.50 
)

This is the standard constructor.

Parameters
typeOne of the enum constants of ECandidateType_t
bkgModeOne of the enum constants of EMCMode_t
cutsD meson cuts (if null, it will use standard cuts)
nMassBinsNumber of bins in the mass axis
rangeRange of the mass axis (will be centered around the PDG mass)

Definition at line 805 of file AliAnalysisTaskDmesonJets.cxx.

AliAnalysisTaskDmesonJets::AnalysisEngine::AnalysisEngine ( const AnalysisEngine source)

Copy constructor

Parameters
sourceReference to a valid AnalysisEngine to copy from.

Definition at line 846 of file AliAnalysisTaskDmesonJets.cxx.

AliAnalysisTaskDmesonJets::AnalysisEngine::~AnalysisEngine ( )
virtual

Definition at line 885 of file AliAnalysisTaskDmesonJets.cxx.

Member Function Documentation

void AliAnalysisTaskDmesonJets::AnalysisEngine::AddInputVectors ( AliEmcalContainer *  cont,
Int_t  offset,
TH2 rejectHist = 0,
Double_t  eff = 0. 
)
private

Adds all the particles contained in the container into the fastjet wrapper

Parameters
contPointer to a valid AliEmcalContainer object

Definition at line 1638 of file AliAnalysisTaskDmesonJets.cxx.

AliAnalysisTaskDmesonJets::AliHFJetDefinition * AliAnalysisTaskDmesonJets::AnalysisEngine::AddJetDefinition ( EJetType_t  type,
Double_t  r,
EJetAlgo_t  algo,
ERecoScheme_t  reco 
)

Add a new jet definition If the jet definition is already present, it does nothing.

Parameters
typeJet type
rJet radius
algoJet algorithm
recoRecombination scheme
Returns
Pointer to the new jet definition (or to the one that was already present)

Definition at line 1056 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AliAnalysisTaskDmesonJets::AddAnalysisEngine().

AliAnalysisTaskDmesonJets::AliHFJetDefinition * AliAnalysisTaskDmesonJets::AnalysisEngine::AddJetDefinition ( const AliHFJetDefinition def)

Add a new jet definition If the jet definition is already present, it does nothing.

Parameters
defReference to a AliJetDefinition object
Returns
Pointer to the new jet definition (or to the one that was already present)

Definition at line 1027 of file AliAnalysisTaskDmesonJets.cxx.

void AliAnalysisTaskDmesonJets::AnalysisEngine::AdoptRDHFCuts ( AliRDHFCuts cuts)

Adopt the cuts (this class owns the cuts object, which will be destroyed when needed).

Parameters
Pointerto a AliRDHFCuts object.

Definition at line 965 of file AliAnalysisTaskDmesonJets.cxx.

void AliAnalysisTaskDmesonJets::AnalysisEngine::AssignDataSlot ( Int_t  n)
inline

Definition at line 459 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::BuildHnSparse ( UInt_t  enabledAxis)

Allocate a THnSparse histogram

Parameters
paramAnalysis parameters used to properly set some of the axis

Definition at line 1847 of file AliAnalysisTaskDmesonJets.cxx.

TTree * AliAnalysisTaskDmesonJets::AnalysisEngine::BuildTree ( const char *  taskName)

Builds the tree where the output will be posted

Returns
Pointer to the new tree

Definition at line 1806 of file AliAnalysisTaskDmesonJets.cxx.

AliAnalysisTaskDmesonJets::EMesonDecayChannel_t AliAnalysisTaskDmesonJets::AnalysisEngine::CheckDecayChannel ( const AliAODMCParticle *  part,
TClonesArray *  mcArray 
)
static

Checks the decay channel of a D meson

Parameters
partPointer to an AliAODMCParticle object for which decay channel is requested
mcArrayPointer to a TClonesArray object where to look for particles
Returns
One of the enum constants of AliAnalysisTaskDmesonJets::EMesonDecayChannel_t (D0->Kpi or D*->D0pi->Kpipi)

Definition at line 1295 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AliHFAODMCParticleContainer::IsSpecialPDG().

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::ExtractD0Attributes ( const AliAODRecoDecayHF2Prong Dcand,
AliDmesonJetInfo DmesonJet,
UInt_t  i 
)
private

Extract attributes of the D0 meson candidate.

Parameters
DcandPointer to a AliAODRecoDecayHF2Prong representing the D0 meson candidate
DmesonJetReference to an AliDmesonJetInfo object where the D0 meson candidate information will be copied
iEither 0 or 1, for the two possible mass hypothesis assignments
Returns
kTRUE on success

Definition at line 1129 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::ExtractDstarAttributes ( const AliAODRecoCascadeHF DstarCand,
AliDmesonJetInfo DmesonJet,
UInt_t  i 
)
private

Extract attributes of the D* meson candidate.

Parameters
DstarCandPointer to a AliAODRecoCascadeHF representing the D* meson candidate
DmesonJetReference to an AliDmesonJetInfo object where the D* meson candidate information will be copied
iEither 0 or 1, for the two possible mass hypothesis assignments (since there is only one mass hypothesis for D*, returns kFALSE for i > 0)
Returns
kTRUE on success

Definition at line 1239 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::ExtractRecoDecayAttributes ( const AliAODRecoDecayHF2Prong Dcand,
AliDmesonJetInfo DmesonJet,
UInt_t  i 
)
private

Extract attributes of the D meson candidate.

Parameters
DcandPointer to a AliAODRecoDecayHF2Prong representing the D meson candidate
DmesonJetReference to an AliDmesonJetInfo object where the D meson candidate information will be copied
iEither 0 or 1, for the two possible mass hypothesis assignments
Returns
kTRUE on success

Definition at line 1109 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::FillHnSparse ( Bool_t  applyKinCuts)

Post the output with D meson jets found in the current event

Returns
kTRUE on success

Definition at line 2246 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::FillHnSparse ( THnSparse *  h,
const AliDmesonJetInfo DmesonJet,
std::string  n 
)

Fill a THnSparse using information from a AliDmesonJetInfo object

Parameters
hValid pointer to a THnSparse object
DmesonJetConst reference to an AliDmesonJetInfo object
nJet name

Definition at line 2290 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::FillQA ( Bool_t  applyKinCuts)

Fills QA histograms. This method is not used by the AliAnalysisTaskDmesonJets task, but can be used by derived tasks that have a custom implementation to fill the output objects.

Returns
Always kTRUE

Definition at line 2125 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::FillTree ( Bool_t  applyKinCuts)

Post the output with D meson jets found in the current event

Returns
kTRUE on success

Definition at line 1998 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::FindJet ( AliAODRecoDecayHF2Prong Dcand,
AliDmesonJetInfo DmesonJet,
AliHFJetDefinition jetDef 
)
private

Find the jet that contains a D meson candidate. The jet finding algorithm is always anti-kt Tracks and clusters are accessed through fTrackContainer and fClusterContainer

Parameters
DcandValid pointer to a D meson candidate object
DmesonJetReference to a AliDmesonJetInfo object where the result will be stored
rJet radius
Returns
kTRUE on success, kFALSE otherwise

Definition at line 1549 of file AliAnalysisTaskDmesonJets.cxx.

std::vector< AliAnalysisTaskDmesonJets::AliHFJetDefinition >::iterator AliAnalysisTaskDmesonJets::AnalysisEngine::FindJetDefinition ( const AliHFJetDefinition eng)

Look for a jet definition that is equal

Parameters
defReference to a jet definition object
Returns
An iterator to the jet definition object, if it is found. An iterator to the end if not found.

Definition at line 1068 of file AliAnalysisTaskDmesonJets.cxx.

AliAODMCParticle * AliAnalysisTaskDmesonJets::AnalysisEngine::FindParticleOrigin ( const AliAODMCParticle *  part,
TClonesArray *  mcArray,
EFindParticleOriginMode_t  mode,
const std::set< UInt_t > &  pdgSet 
)
static

Finds a particle in the fragmentation tree of a final state particle

Parameters
partPointer to an AliAODMCParticle object for which originating quark is required
mcArrayPointer to a TClonesArray object where to look for particles
modeSee documentation of the enum type EFindParticleOriginMode_t
pdgSetA set of PDG codes that are being searched
Returns
A pointer to the MC particle found in the fragmentation tree

Definition at line 1392 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AliHFAODMCParticleContainer::IsSpecialPDG().

AliAODMCParticle * AliAnalysisTaskDmesonJets::AnalysisEngine::FindParticleOrigin ( const AliAODMCParticle *  part,
TClonesArray *  mcArray,
EFindParticleOriginMode_t  mode 
)
static

Finds a particle in the fragmentation tree of a final state particle

Parameters
partPointer to an AliAODMCParticle object for which originating quark is required
mcArrayPointer to a TClonesArray object where to look for particles
modeSee documentation of the enum type EFindParticleOriginMode_t
Returns
A pointer to the MC particle found in the fragmentation tree

Definition at line 1376 of file AliAnalysisTaskDmesonJets.cxx.

const char* AliAnalysisTaskDmesonJets::AnalysisEngine::GetCandidateName ( ) const
inline

Definition at line 435 of file AliAnalysisTaskDmesonJets.h.

ECandidateType_t AliAnalysisTaskDmesonJets::AnalysisEngine::GetCandidateType ( ) const
inline

Definition at line 440 of file AliAnalysisTaskDmesonJets.h.

Int_t AliAnalysisTaskDmesonJets::AnalysisEngine::GetDataSlotNumber ( ) const
inline
std::map<int, AliDmesonJetInfo>& AliAnalysisTaskDmesonJets::AnalysisEngine::GetDmesons ( )
inline

Definition at line 449 of file AliAnalysisTaskDmesonJets.h.

std::vector<AliAnalysisTaskDmesonJets::AliHFJetDefinition>& AliAnalysisTaskDmesonJets::AnalysisEngine::GetJetDefinitions ( )
inline

Definition at line 445 of file AliAnalysisTaskDmesonJets.h.

EMCMode_t AliAnalysisTaskDmesonJets::AnalysisEngine::GetMCMode ( ) const
inline

Definition at line 439 of file AliAnalysisTaskDmesonJets.h.

const char * AliAnalysisTaskDmesonJets::AnalysisEngine::GetName ( ) const

Generate a name for this analysis parameter set

Parameters
iIndex of the jet radius array.

Definition at line 996 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AliAnalysisTaskDmesonJets::AddAnalysisEngine().

const char * AliAnalysisTaskDmesonJets::AnalysisEngine::GetName ( const AliHFJetDefinition jetDef) const

Generate a name for this analysis parameter set

Parameters
iIndex of the jet radius array.

Definition at line 984 of file AliAnalysisTaskDmesonJets.cxx.

TTree* AliAnalysisTaskDmesonJets::AnalysisEngine::GetTree ( ) const
inline
void AliAnalysisTaskDmesonJets::AnalysisEngine::Init ( const AliEMCALGeometry *const  geom,
Int_t  runNumber 
)

Initialize the analysis engine.

Definition at line 912 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::IsAnyJetInAcceptance ( const AliDmesonJetInfo dMesonJet) const

Checks whether any of the D meson jets is in the acceptance

Parameters
Constreference to a valid AliDmesonJetInfo object

Definition at line 902 of file AliAnalysisTaskDmesonJets.cxx.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::IsInhibit ( ) const
inline

Definition at line 468 of file AliAnalysisTaskDmesonJets.h.

std::pair< AliAnalysisTaskDmesonJets::EMesonOrigin_t, AliAODMCParticle * > AliAnalysisTaskDmesonJets::AnalysisEngine::IsPromptCharm ( const AliAODMCParticle *  part,
TClonesArray *  mcArray 
)
static

Checks whether a particle is the result of the hadronization of a charm quark or a bottom quark

Parameters
partPointer to an AliAODMCParticle object for which originating quark is required
mcArrayPointer to a TClonesArray object where to look for particles
Returns
A pair: first is either kFromCharm or kFromBottom; second is the pointer to the quark

Definition at line 1348 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AliHFAODMCParticleContainer::IsSpecialPDG().

AliAODMCParticle* AliAnalysisTaskDmesonJets::AnalysisEngine::MatchToMC ( ) const
private
AliAnalysisTaskDmesonJets::AnalysisEngine & AliAnalysisTaskDmesonJets::AnalysisEngine::operator= ( const AnalysisEngine source)

Assignement operator

Parameters
sourceReference to a valid AnalysisEngine to copy from.

Definition at line 893 of file AliAnalysisTaskDmesonJets.cxx.

void AliAnalysisTaskDmesonJets::AnalysisEngine::RunAnalysis ( )
protected

Run the analysis.

Definition at line 1420 of file AliAnalysisTaskDmesonJets.cxx.

void AliAnalysisTaskDmesonJets::AnalysisEngine::RunDetectorLevelAnalysis ( )
private

Run a detector level analysis.

Definition at line 1435 of file AliAnalysisTaskDmesonJets.cxx.

void AliAnalysisTaskDmesonJets::AnalysisEngine::RunParticleLevelAnalysis ( )
private

Run a particle level analysis.

Definition at line 1660 of file AliAnalysisTaskDmesonJets.cxx.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetAcceptedDecayMap ( UInt_t  m)
inline

Definition at line 432 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetCandidateProperties ( Double_t  range)
private

Sets the D meson candidate properties.

Parameters
rangeRange of the mass axis (will be centered around the PDG mass)

Definition at line 919 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AnalysisEngine().

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetCandidateType ( ECandidateType_t  t)
inline

Definition at line 425 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetMassRange ( Double_t  min,
Double_t  max 
)
inline

Definition at line 428 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetMCMode ( EMCMode_t  m)
inline

Definition at line 426 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetNMassBins ( Int_t  n)
inline

Definition at line 427 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetRDHFCuts ( AliRDHFCuts cuts)

Set the cuts (creates a copy, so the original object is not owned by this class).

Parameters
Pointerto a AliRDHFCuts object.

Definition at line 974 of file AliAnalysisTaskDmesonJets.cxx.

Referenced by AliAnalysisTaskDmesonJets::AddAnalysisEngine(), and AnalysisEngine().

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetRejectedOriginMap ( UInt_t  m)
inline

Definition at line 431 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetRejectISR ( Bool_t  b)
inline

Definition at line 433 of file AliAnalysisTaskDmesonJets.h.

void AliAnalysisTaskDmesonJets::AnalysisEngine::SetTrackEfficiency ( Double_t  t)
inline

Definition at line 458 of file AliAnalysisTaskDmesonJets.h.

Friends And Related Function Documentation

friend class AliAnalysisTaskDmesonJets
friend

Definition at line 517 of file AliAnalysisTaskDmesonJets.h.

bool operator!= ( const AnalysisEngine lhs,
const AnalysisEngine rhs 
)
friend

Definition at line 476 of file AliAnalysisTaskDmesonJets.h.

bool operator< ( const AnalysisEngine lhs,
const AnalysisEngine rhs 
)
friend

Compares 2 analysis engines. The ordering is based on the candidate type first and then on the MC mode.

Parameters
lhsReference to the first AnalysisEngine object
rhsReference to the second AnalysisEngine object

Definition at line 1080 of file AliAnalysisTaskDmesonJets.cxx.

bool operator<= ( const AnalysisEngine lhs,
const AnalysisEngine rhs 
)
friend

Definition at line 472 of file AliAnalysisTaskDmesonJets.h.

bool operator== ( const AnalysisEngine lhs,
const AnalysisEngine rhs 
)
friend

Compares 2 analysis engines. Two analysis engines are considerate equal if they have both the same candidate type and MC mode.

Parameters
lhsReference to the first AnalysisEngine object
rhsReference to the second AnalysisEngine object

Definition at line 1095 of file AliAnalysisTaskDmesonJets.cxx.

bool operator> ( const AnalysisEngine lhs,
const AnalysisEngine rhs 
)
friend

Definition at line 471 of file AliAnalysisTaskDmesonJets.h.

bool operator>= ( const AnalysisEngine lhs,
const AnalysisEngine rhs 
)
friend

Definition at line 473 of file AliAnalysisTaskDmesonJets.h.

Member Data Documentation

UInt_t AliAnalysisTaskDmesonJets::AnalysisEngine::fAcceptedDecay
protected

Bit mask with D meson decays that are accepted (only used for particle-level analysis)

Definition at line 495 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetAcceptedDecayMap().

AliAODEvent* AliAnalysisTaskDmesonJets::AnalysisEngine::fAodEvent
protected

! AOD event

Definition at line 512 of file AliAnalysisTaskDmesonJets.h.

TString AliAnalysisTaskDmesonJets::AnalysisEngine::fBranchName
protected

AOD branch where the D meson candidate are found.

Definition at line 488 of file AliAnalysisTaskDmesonJets.h.

TClonesArray* AliAnalysisTaskDmesonJets::AnalysisEngine::fCandidateArray
protected

! D meson candidate array

Definition at line 508 of file AliAnalysisTaskDmesonJets.h.

TString AliAnalysisTaskDmesonJets::AnalysisEngine::fCandidateName
protected

Candidate name.

Definition at line 484 of file AliAnalysisTaskDmesonJets.h.

Referenced by GetCandidateName().

UInt_t AliAnalysisTaskDmesonJets::AnalysisEngine::fCandidatePDG
protected

Candidate PDG.

Definition at line 485 of file AliAnalysisTaskDmesonJets.h.

ECandidateType_t AliAnalysisTaskDmesonJets::AnalysisEngine::fCandidateType
protected

Candidate type.

Definition at line 483 of file AliAnalysisTaskDmesonJets.h.

Referenced by GetCandidateType(), operator<(), operator==(), and SetCandidateType().

Double_t AliAnalysisTaskDmesonJets::AnalysisEngine::fCent
protected

! Event centrality

Definition at line 515 of file AliAnalysisTaskDmesonJets.h.

std::vector<AliClusterContainer*> AliAnalysisTaskDmesonJets::AnalysisEngine::fClusterContainers
protected

! Cluster containers

Definition at line 511 of file AliAnalysisTaskDmesonJets.h.

AliDmesonInfoSummary* AliAnalysisTaskDmesonJets::AnalysisEngine::fCurrentDmesonJetInfo
protected

! Current D meson jet info

Definition at line 505 of file AliAnalysisTaskDmesonJets.h.

AliJetInfoSummary** AliAnalysisTaskDmesonJets::AnalysisEngine::fCurrentJetInfo
protected

! Current jet info

Definition at line 506 of file AliAnalysisTaskDmesonJets.h.

Int_t AliAnalysisTaskDmesonJets::AnalysisEngine::fDataSlotNumber
protected

! Data slot where the tree output is posted

Definition at line 503 of file AliAnalysisTaskDmesonJets.h.

Referenced by AssignDataSlot(), and GetDataSlotNumber().

std::map<int, AliDmesonJetInfo> AliAnalysisTaskDmesonJets::AnalysisEngine::fDmesonJets
protected

! Array containing the D meson jets

Definition at line 507 of file AliAnalysisTaskDmesonJets.h.

Referenced by GetDmesons().

AliFJWrapper* AliAnalysisTaskDmesonJets::AnalysisEngine::fFastJetWrapper
protected

! Fastjet wrapper

Definition at line 513 of file AliAnalysisTaskDmesonJets.h.

THistManager* AliAnalysisTaskDmesonJets::AnalysisEngine::fHistManager
protected

! Histograms

Definition at line 514 of file AliAnalysisTaskDmesonJets.h.

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::fInhibit
protected

Inhibit the task.

Definition at line 496 of file AliAnalysisTaskDmesonJets.h.

Referenced by IsInhibit().

std::vector<AliHFJetDefinition> AliAnalysisTaskDmesonJets::AnalysisEngine::fJetDefinitions
protected

Jet definitions.

Definition at line 497 of file AliAnalysisTaskDmesonJets.h.

Referenced by AliAnalysisTaskDmesonJets::AddAnalysisEngine(), and GetJetDefinitions().

Double_t AliAnalysisTaskDmesonJets::AnalysisEngine::fMaxMass
protected

Max mass in histogram axis.

Definition at line 492 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetMassRange().

Float_t AliAnalysisTaskDmesonJets::AnalysisEngine::fMaxPt
protected

Histogram pt limit.

Definition at line 499 of file AliAnalysisTaskDmesonJets.h.

AliHFAODMCParticleContainer* AliAnalysisTaskDmesonJets::AnalysisEngine::fMCContainer
protected

! MC particle container

Definition at line 509 of file AliAnalysisTaskDmesonJets.h.

EMCMode_t AliAnalysisTaskDmesonJets::AnalysisEngine::fMCMode
protected

MC mode: No MC (data and MC detector level), background-only (MC), signal-only (MC), MC truth (particle level)

Definition at line 489 of file AliAnalysisTaskDmesonJets.h.

Referenced by GetMCMode(), operator<(), operator==(), and SetMCMode().

Double_t AliAnalysisTaskDmesonJets::AnalysisEngine::fMinMass
protected

Min mass in histogram axis.

Definition at line 491 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetMassRange().

UChar_t AliAnalysisTaskDmesonJets::AnalysisEngine::fNDaughters
protected

Number of daughters.

Definition at line 486 of file AliAnalysisTaskDmesonJets.h.

Int_t AliAnalysisTaskDmesonJets::AnalysisEngine::fNMassBins
protected

Mass number of bins.

Definition at line 490 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetNMassBins().

std::map<AliAODMCParticle*, Short_t> AliAnalysisTaskDmesonJets::AnalysisEngine::fPartons

! set of the partons in the shower that produced each D meson

Definition at line 478 of file AliAnalysisTaskDmesonJets.h.

TArrayI AliAnalysisTaskDmesonJets::AnalysisEngine::fPDGdaughters
protected

List of the PDG code of the daughters.

Definition at line 487 of file AliAnalysisTaskDmesonJets.h.

Float_t AliAnalysisTaskDmesonJets::AnalysisEngine::fPtBinWidth
protected

Histogram pt bin width.

Definition at line 498 of file AliAnalysisTaskDmesonJets.h.

TRandom* AliAnalysisTaskDmesonJets::AnalysisEngine::fRandomGen
protected

! Random number generator

Definition at line 500 of file AliAnalysisTaskDmesonJets.h.

AliRDHFCuts* AliAnalysisTaskDmesonJets::AnalysisEngine::fRDHFCuts
protected

D meson candidates cuts.

Definition at line 493 of file AliAnalysisTaskDmesonJets.h.

Referenced by AliAnalysisTaskDmesonJets::AddAnalysisEngine(), and AnalysisEngine().

UInt_t AliAnalysisTaskDmesonJets::AnalysisEngine::fRejectedOrigin
protected

Bit mask with D meson origins that are rejected (used for MC analysis, i.e. signal-only, background-only and particle-level)

Definition at line 494 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetRejectedOriginMap().

Bool_t AliAnalysisTaskDmesonJets::AnalysisEngine::fRejectISR
protected

! Reject initial state radiation

Definition at line 502 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetRejectISR().

std::vector<AliTrackContainer*> AliAnalysisTaskDmesonJets::AnalysisEngine::fTrackContainers
protected

! Track containers

Definition at line 510 of file AliAnalysisTaskDmesonJets.h.

Double_t AliAnalysisTaskDmesonJets::AnalysisEngine::fTrackEfficiency
protected

! Artificial tracking inefficiency (0...1) -> set automatically at ExecOnce by AliAnalysisTaskDmesonJets

Definition at line 501 of file AliAnalysisTaskDmesonJets.h.

Referenced by SetTrackEfficiency().

TTree* AliAnalysisTaskDmesonJets::AnalysisEngine::fTree
protected

! Output tree

Definition at line 504 of file AliAnalysisTaskDmesonJets.h.

Referenced by GetTree().


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