baseAnalysis Class Reference

Base class for all HepMCanalysis classes to implement common functionality and interfaces. More...

#include <baseAnalysis.h>

Inheritance diagram for baseAnalysis:

bbbarAnalysis ElasScatAnalysis EtMissAnalysis JetAnalysis JetFinder ttbarAnalysis UEAnalysis UserAnalysis WplusJetAnalysis WtaunuAnalysis ZAnalysis ZtautauAnalysis List of all members.

Public Member Functions

 baseAnalysis ()
virtual ~baseAnalysis ()
virtual int Finalize ()
virtual int Finalize (TFile *output)
virtual TH1D * popHisto ()
virtual int Init (double maxeta=2.5, double minpt=0.5)
virtual int Process (HepMC::GenEvent *hepmcevt)
virtual std::vector< TH1D * > averagedHistograms ()
virtual std::vector< TH1D * > & Histograms ()
int InitJetFinder (double coneRadius, double overlapThreshold, double jet_ptmin, double lepton_ptmin, double DeltaR_lepton_track, int jetalgorithm=2)
int FindJetableParticles (HepMC::GenEvent *hepmcevt)
int FindJet (HepMC::GenEvent *hepmcevt)
int DeleteJetObject (HepMC::GenEvent *hepmcevt)
std::vector< fastjet::PseudoJet > GetJet (HepMC::GenEvent *hepmcevt)
virtual void SetJet (std::vector< fastjet::PseudoJet > *inclusive_jets)
int FindMissingEt (HepMC::GenEvent *hepmcevt)
int ClearMissingEt (HepMC::GenEvent *hepmcevt)
int ClearEvent (HepMC::GenEvent *hepmcevt)
int setOutputFileName (const std::string &filename)
int setOutputRootDir (const std::string &dirname)
virtual int IsFinalStateParticle (HepMC::GenParticle *p)
int setMaxEtaCut (double maxeta)
int setMinPtCut (double minpt)
TH1D * initHist (const std::string &name, const std::string &title, const std::string &xlabel, int nrBins=100, double xmin=0.0, double xmax=100.0)
TH1D * initHistVariableBin (const std::string &name, const std::string &title, const std::string &xlabel, int nbin, double nbinRange[])
bool checkDaughter (HepMC::GenParticle *mother, HepMC::GenParticle *daughter, int maxGenerations=-1)
bool trackfromPID (int pid, HepMC::GenParticle *track, int maxGenerations=-1)

Static Public Member Functions

static bool IsNeutrino (int pid)
static bool IsGamma (int pid)
static bool IsNeutron (int pid)
static bool IsK0 (int pid)
static bool IsPi0 (int pid)
static bool IsElectron (int pid)
static bool IsMuon (int pid)
static bool IsGluon (int pid)
static bool chargedParticle (int pid)
static double getRapidity (const HepMC::GenEvent::particle_const_iterator &p)

Protected Attributes

fastjet::JetDefinition::Plugin * m_plugin
fastjet::JetDefinition * m_jetDef
std::vector< fastjet::PseudoJet > m_input_particles
std::vector< fastjet::PseudoJet > m_inclusive_jets
fastjet::ClusterSequence * m_clust_seq
std::vector< TH1D * > m_histVector
std::string m_outputFileName
std::string m_outputRootDir
double m_coneRadius
double m_overlapThreshold
double m_jet_ptmin
double m_lepton_ptmin
double m_DeltaR_lepton_track
bool m_Jetfinder_enabled
double m_max_eta
double m_min_pt
double exMissTruth
double eyMissTruth
double etMissTruth
double etsumMissTruth
size_t m_NanNum

Detailed Description

Base class for all HepMCanalysis classes to implement common functionality and interfaces.

Author:
Cano Ay Dec 2008

Definition at line 33 of file baseAnalysis.h.


Constructor & Destructor Documentation

baseAnalysis::baseAnalysis (  ) 

Default constructor.

Definition at line 35 of file baseAnalysis.cc.

References etMissTruth, etsumMissTruth, exMissTruth, eyMissTruth, m_clust_seq, m_coneRadius, m_DeltaR_lepton_track, m_jet_ptmin, m_jetDef, m_Jetfinder_enabled, m_lepton_ptmin, m_max_eta, m_min_pt, m_NanNum, m_overlapThreshold, and m_plugin.

baseAnalysis::~baseAnalysis (  )  [virtual]

Virtual destructor.

Definition at line 65 of file baseAnalysis.cc.

References m_NanNum.


Member Function Documentation

vector< TH1D * > baseAnalysis::averagedHistograms (  )  [virtual]

Reimplemented in UEAnalysis.

Definition at line 119 of file baseAnalysis.cc.

References m_histVector.

bool baseAnalysis::chargedParticle ( int  pid  )  [static]

Check if neutral particle

Definition at line 520 of file baseAnalysis.cc.

References IsGamma(), IsGluon(), IsK0(), IsNeutrino(), IsNeutron(), and IsPi0().

Referenced by bbbarAnalysis::getBbbar(), ZtautauAnalysis::Process(), ZAnalysis::Process(), WtaunuAnalysis::Process(), WplusJetAnalysis::Process(), UEAnalysis::Process(), ttbarAnalysis::Process(), and JetAnalysis::Process().

bool baseAnalysis::checkDaughter ( HepMC::GenParticle *  mother,
HepMC::GenParticle *  daughter,
int  maxGenerations = -1 
)

check if mother decayed into daugther

Definition at line 669 of file baseAnalysis.cc.

int baseAnalysis::ClearEvent ( HepMC::GenEvent *  hepmcevt  ) 

Clear some values from the event

Definition at line 437 of file baseAnalysis.cc.

References ClearMissingEt(), and DeleteJetObject().

Referenced by main().

int baseAnalysis::ClearMissingEt ( HepMC::GenEvent *  hepmcevt  ) 

ClearMissingEt: setting the variables to calculate missing energy to zero.

Definition at line 404 of file baseAnalysis.cc.

References etMissTruth, etsumMissTruth, exMissTruth, and eyMissTruth.

Referenced by ClearEvent().

int baseAnalysis::DeleteJetObject ( HepMC::GenEvent *  hepmcevt  ) 

DeleteJetObject: delete all the jet objects.

Definition at line 349 of file baseAnalysis.cc.

References m_clust_seq, m_inclusive_jets, and m_input_particles.

Referenced by ClearEvent().

int baseAnalysis::Finalize ( TFile *  output  )  [virtual]

In the final step all the histogramms are stored in a rootfile. The name of the rootfile can be set with the function setOutputFileName(const string &fileName).

Reimplemented in UEAnalysis.

Definition at line 563 of file baseAnalysis.cc.

References m_histVector, and m_outputRootDir.

int baseAnalysis::Finalize (  )  [virtual]

Finalize: In the final step all the histogramms are stored in a rootfile. The name of the rootfile can be set with the function setOutputFileName(const string& fileName).

Reimplemented in UEAnalysis.

Definition at line 549 of file baseAnalysis.cc.

Referenced by UEAnalysis::Finalize().

int baseAnalysis::FindJet ( HepMC::GenEvent *  hepmcevt  ) 

Run JetFinder.

Definition at line 300 of file baseAnalysis.cc.

References FindJetableParticles(), m_clust_seq, m_inclusive_jets, m_input_particles, m_jet_ptmin, m_jetDef, and m_Jetfinder_enabled.

Referenced by GetJet().

int baseAnalysis::FindJetableParticles ( HepMC::GenEvent *  hepmcevt  ) 

Find Jetable Particles.

Definition at line 185 of file baseAnalysis.cc.

References IsElectron(), IsMuon(), IsNeutrino(), IsNeutron(), m_DeltaR_lepton_track, m_input_particles, m_jetDef, m_Jetfinder_enabled, m_lepton_ptmin, m_max_eta, m_min_pt, m_NanNum, and q.

Referenced by FindJet().

int baseAnalysis::FindMissingEt ( HepMC::GenEvent *  hepmcevt  ) 

Calculate Missing Et

Definition at line 363 of file baseAnalysis.cc.

References etMissTruth, etsumMissTruth, exMissTruth, eyMissTruth, and IsNeutrino().

Referenced by EtMissAnalysis::Process().

vector< fastjet::PseudoJet > baseAnalysis::GetJet ( HepMC::GenEvent *  hepmcevt  ) 

Definition at line 331 of file baseAnalysis.cc.

References FindJet(), and m_inclusive_jets.

Referenced by main().

double baseAnalysis::getRapidity ( const HepMC::GenEvent::particle_const_iterator &  p  )  [static]

calculate the rapidity of a particle

Definition at line 739 of file baseAnalysis.cc.

References log.

Referenced by WtaunuAnalysis::Process(), and WplusJetAnalysis::Process().

vector< TH1D * > & baseAnalysis::Histograms (  )  [virtual]

Definition at line 125 of file baseAnalysis.cc.

References m_histVector.

int baseAnalysis::Init ( double  maxeta = 2.5,
double  minpt = 0.5 
) [virtual]

Reimplemented in bbbarAnalysis, ElasScatAnalysis, EtMissAnalysis, JetAnalysis, JetFinder, ttbarAnalysis, UEAnalysis, UserAnalysis, WplusJetAnalysis, WtaunuAnalysis, ZAnalysis, and ZtautauAnalysis.

Definition at line 103 of file baseAnalysis.cc.

TH1D* baseAnalysis::initHist ( const std::string &  name,
const std::string &  title,
const std::string &  xlabel,
int  nrBins = 100,
double  xmin = 0.0,
double  xmax = 100.0 
)

Initialization of histograms

Referenced by ZtautauAnalysis::Init(), ZAnalysis::Init(), WtaunuAnalysis::Init(), WplusJetAnalysis::Init(), UserAnalysis::Init(), UEAnalysis::Init(), ttbarAnalysis::Init(), JetAnalysis::Init(), EtMissAnalysis::Init(), ElasScatAnalysis::Init(), and bbbarAnalysis::Init().

TH1D* baseAnalysis::initHistVariableBin ( const std::string &  name,
const std::string &  title,
const std::string &  xlabel,
int  nbin,
double  nbinRange[] 
)

Referenced by UEAnalysis::Init().

int baseAnalysis::InitJetFinder ( double  coneRadius,
double  overlapThreshold,
double  jet_ptmin,
double  lepton_ptmin,
double  DeltaR_lepton_track,
int  jetalgorithm = 2 
)

Initialize FastJet Algorithm

Definition at line 133 of file baseAnalysis.cc.

References m_coneRadius, m_DeltaR_lepton_track, m_jet_ptmin, m_jetDef, m_Jetfinder_enabled, m_lepton_ptmin, m_overlapThreshold, and m_plugin.

Referenced by main().

bool baseAnalysis::IsElectron ( int  pid  )  [static]

Definition at line 499 of file baseAnalysis.cc.

Referenced by FindJetableParticles().

int baseAnalysis::IsFinalStateParticle ( HepMC::GenParticle *  p  )  [virtual]

Check if final state particle

Definition at line 417 of file baseAnalysis.cc.

Referenced by ZAnalysis::FSVector(), ZtautauAnalysis::Process(), ZAnalysis::Process(), WtaunuAnalysis::Process(), WplusJetAnalysis::Process(), UEAnalysis::Process(), ttbarAnalysis::Process(), and JetAnalysis::Process().

bool baseAnalysis::IsGamma ( int  pid  )  [static]

Definition at line 469 of file baseAnalysis.cc.

Referenced by chargedParticle().

bool baseAnalysis::IsGluon ( int  pid  )  [static]

Definition at line 513 of file baseAnalysis.cc.

Referenced by chargedParticle().

bool baseAnalysis::IsK0 ( int  pid  )  [static]

Definition at line 485 of file baseAnalysis.cc.

Referenced by chargedParticle().

bool baseAnalysis::IsMuon ( int  pid  )  [static]

Definition at line 506 of file baseAnalysis.cc.

Referenced by FindJetableParticles().

bool baseAnalysis::IsNeutrino ( int  pid  )  [static]

Check some special neutral Particles

Definition at line 460 of file baseAnalysis.cc.

Referenced by chargedParticle(), FindJetableParticles(), and FindMissingEt().

bool baseAnalysis::IsNeutron ( int  pid  )  [static]

Definition at line 478 of file baseAnalysis.cc.

Referenced by chargedParticle(), and FindJetableParticles().

bool baseAnalysis::IsPi0 ( int  pid  )  [static]

Definition at line 492 of file baseAnalysis.cc.

Referenced by chargedParticle().

TH1D * baseAnalysis::popHisto (  )  [virtual]

popHisto: return the last entry of the histogram vector, needed for ATHENA (ATLAS software) implementation.

Definition at line 589 of file baseAnalysis.cc.

References m_histVector.

int baseAnalysis::Process ( HepMC::GenEvent *  hepmcevt  )  [virtual]

Reimplemented in bbbarAnalysis, ElasScatAnalysis, EtMissAnalysis, JetAnalysis, ttbarAnalysis, UEAnalysis, UserAnalysis, WplusJetAnalysis, WtaunuAnalysis, ZAnalysis, and ZtautauAnalysis.

Definition at line 111 of file baseAnalysis.cc.

virtual void baseAnalysis::SetJet ( std::vector< fastjet::PseudoJet > *  inclusive_jets  )  [virtual]

int baseAnalysis::setMaxEtaCut ( double  maxeta  )  [inline]

Set the maximum allowed eta range

Definition at line 532 of file baseAnalysis.cc.

References m_max_eta.

int baseAnalysis::setMinPtCut ( double  minpt  )  [inline]

Set maximum pt of tracks

Definition at line 539 of file baseAnalysis.cc.

References m_min_pt.

int baseAnalysis::setOutputFileName ( const std::string &  filename  )  [inline]

Set the Output filename

int baseAnalysis::setOutputRootDir ( const std::string &  dirname  )  [inline]

Set the directory name in Output root file

bool baseAnalysis::trackfromPID ( int  pid,
HepMC::GenParticle *  track,
int  maxGenerations = -1 
)

check if the track comes from a specific particle e.g. pid(W-Boson)=24

Definition at line 702 of file baseAnalysis.cc.


Member Data Documentation

double baseAnalysis::etMissTruth [protected]

Definition at line 144 of file baseAnalysis.h.

Referenced by baseAnalysis(), ClearMissingEt(), FindMissingEt(), and EtMissAnalysis::Process().

double baseAnalysis::etsumMissTruth [protected]

Definition at line 145 of file baseAnalysis.h.

Referenced by baseAnalysis(), ClearMissingEt(), FindMissingEt(), and EtMissAnalysis::Process().

double baseAnalysis::exMissTruth [protected]

Definition at line 142 of file baseAnalysis.h.

Referenced by baseAnalysis(), ClearMissingEt(), FindMissingEt(), and EtMissAnalysis::Process().

double baseAnalysis::eyMissTruth [protected]

Definition at line 143 of file baseAnalysis.h.

Referenced by baseAnalysis(), ClearMissingEt(), FindMissingEt(), and EtMissAnalysis::Process().

fastjet::ClusterSequence* baseAnalysis::m_clust_seq [protected]

Definition at line 121 of file baseAnalysis.h.

Referenced by baseAnalysis(), DeleteJetObject(), and FindJet().

double baseAnalysis::m_coneRadius [protected]

Reimplemented in WplusJetAnalysis.

Definition at line 130 of file baseAnalysis.h.

Referenced by baseAnalysis(), and InitJetFinder().

double baseAnalysis::m_DeltaR_lepton_track [protected]

Reimplemented in WplusJetAnalysis.

Definition at line 134 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJetableParticles(), and InitJetFinder().

std::vector<TH1D*> baseAnalysis::m_histVector [protected]

Definition at line 124 of file baseAnalysis.h.

Referenced by UEAnalysis::averagedHistograms(), averagedHistograms(), Finalize(), Histograms(), WplusJetAnalysis::Init(), and popHisto().

std::vector< fastjet::PseudoJet > baseAnalysis::m_inclusive_jets [protected]

Definition at line 120 of file baseAnalysis.h.

Referenced by DeleteJetObject(), FindJet(), GetJet(), ZtautauAnalysis::Process(), ZAnalysis::Process(), WplusJetAnalysis::Process(), UserAnalysis::Process(), UEAnalysis::Process(), ttbarAnalysis::Process(), and JetAnalysis::Process().

std::vector< fastjet::PseudoJet > baseAnalysis::m_input_particles [protected]

Definition at line 119 of file baseAnalysis.h.

Referenced by DeleteJetObject(), FindJet(), and FindJetableParticles().

double baseAnalysis::m_jet_ptmin [protected]

Reimplemented in WplusJetAnalysis.

Definition at line 132 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJet(), and InitJetFinder().

fastjet::JetDefinition* baseAnalysis::m_jetDef [protected]

Definition at line 117 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJet(), FindJetableParticles(), and InitJetFinder().

bool baseAnalysis::m_Jetfinder_enabled [protected]

Definition at line 135 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJet(), FindJetableParticles(), and InitJetFinder().

double baseAnalysis::m_lepton_ptmin [protected]

Reimplemented in WplusJetAnalysis.

Definition at line 133 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJetableParticles(), and InitJetFinder().

double baseAnalysis::m_max_eta [protected]

Definition at line 138 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJetableParticles(), ZtautauAnalysis::Init(), ZAnalysis::Init(), WtaunuAnalysis::Init(), WplusJetAnalysis::Init(), UserAnalysis::Init(), UEAnalysis::Init(), ttbarAnalysis::Init(), JetFinder::Init(), JetAnalysis::Init(), EtMissAnalysis::Init(), ElasScatAnalysis::Init(), bbbarAnalysis::Init(), and setMaxEtaCut().

double baseAnalysis::m_min_pt [protected]

Definition at line 139 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJetableParticles(), ZtautauAnalysis::Init(), ZAnalysis::Init(), WtaunuAnalysis::Init(), WplusJetAnalysis::Init(), UserAnalysis::Init(), UEAnalysis::Init(), ttbarAnalysis::Init(), JetFinder::Init(), JetAnalysis::Init(), EtMissAnalysis::Init(), ElasScatAnalysis::Init(), bbbarAnalysis::Init(), and setMinPtCut().

size_t baseAnalysis::m_NanNum [protected]

Definition at line 147 of file baseAnalysis.h.

Referenced by baseAnalysis(), FindJetableParticles(), and ~baseAnalysis().

std::string baseAnalysis::m_outputFileName [protected]

Definition at line 126 of file baseAnalysis.h.

Referenced by ZtautauAnalysis::Init(), ZAnalysis::Init(), WtaunuAnalysis::Init(), WplusJetAnalysis::Init(), UserAnalysis::Init(), UEAnalysis::Init(), ttbarAnalysis::Init(), JetAnalysis::Init(), EtMissAnalysis::Init(), ElasScatAnalysis::Init(), and bbbarAnalysis::Init().

std::string baseAnalysis::m_outputRootDir [protected]

Definition at line 127 of file baseAnalysis.h.

Referenced by Finalize(), ZtautauAnalysis::Init(), ZAnalysis::Init(), WtaunuAnalysis::Init(), WplusJetAnalysis::Init(), UserAnalysis::Init(), UEAnalysis::Init(), ttbarAnalysis::Init(), JetAnalysis::Init(), EtMissAnalysis::Init(), ElasScatAnalysis::Init(), and bbbarAnalysis::Init().

double baseAnalysis::m_overlapThreshold [protected]

Reimplemented in WplusJetAnalysis.

Definition at line 131 of file baseAnalysis.h.

Referenced by baseAnalysis(), and InitJetFinder().

fastjet::JetDefinition::Plugin* baseAnalysis::m_plugin [protected]

jet finding

Definition at line 116 of file baseAnalysis.h.

Referenced by baseAnalysis(), and InitJetFinder().


The documentation for this class was generated from the following files:
Generated on Wed Aug 31 09:44:59 2011 for HepMCAnalysis by  doxygen 1.4.7