Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
JetAnalyzer Class Reference

#include <analysis/blob/master/sPhenixAj/src/JetAnalyzer.hh>

+ Inheritance diagram for JetAnalyzer:
+ Collaboration diagram for JetAnalyzer:

Public Member Functions

 JetAnalyzer (std::vector< fastjet::PseudoJet > &InOrigParticles, fastjet::JetDefinition &JetDef, fastjet::AreaDefinition &AreaDef, fastjet::Selector selector_bkgd=fastjet::SelectorAbsRapMax(0.6)*(!fastjet::SelectorNHardest(2)))
 
 JetAnalyzer (std::vector< fastjet::PseudoJet > &InOrigParticles, fastjet::JetDefinition &JetDef)
 
 ~JetAnalyzer ()
 
fastjet::Subtractor * GetBackgroundSubtractor ()
 
fastjet::JetMedianBackgroundEstimator * GetBackgroundEstimator ()
 

Static Public Member Functions

static double phimod2pi (double phi)
 

Static Public Attributes

static const double pi = 3.141592653589793238462643383279502884
 

Private Attributes

std::vector< fastjet::PseudoJet > & OrigParticles
 
bool CanDoBackground
 
fastjet::JetMedianBackgroundEstimator * bkgd_estimator
 
fastjet::Subtractor * bkgd_subtractor
 
fastjet::Selector selector_bkgd
 
fastjet::JetDefinition * jet_def_bkgd
 
fastjet::AreaDefinition * area_def_bkgd
 

Detailed Description

The main class.

Definition at line 74 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 74 of file JetAnalyzer.hh

Constructor & Destructor Documentation

JetAnalyzer::JetAnalyzer ( std::vector< fastjet::PseudoJet > &  InOrigParticles,
fastjet::JetDefinition &  JetDef,
fastjet::AreaDefinition &  AreaDef,
fastjet::Selector  selector_bkgd = fastjet::SelectorAbsRapMax( 0.6 ) * (!fastjet::SelectorNHardest(2)) 
)

Standard constructor. Pass through to Fastjet and set up internal background estimator with default values. Note that initialization as ClusterSequenceArea( pHi, jet_def, 0 ) works fine and skips the area computation. Instead, we provide a second ctor and use this information to determine whether to provide background capability.

Parameters
InOrigParticlesis the full set of constituent candidates. Passed by reference!
JetDefis a fastjet::JetDefinition for the clustering. Passed by reference!
AreaDefis a fastjet::AreaDefinition for the clustering
selector_bkgdis a fastjet::Selector for background subtraction

Definition at line 18 of file JetAnalyzer.cxx.

View newest version in sPHENIX GitHub at line 18 of file JetAnalyzer.cxx

References area_def_bkgd, bkgd_estimator, bkgd_subtractor, CanDoBackground, and jet_def_bkgd.

JetAnalyzer::JetAnalyzer ( std::vector< fastjet::PseudoJet > &  InOrigParticles,
fastjet::JetDefinition &  JetDef 
)

Use as ClusterSequence, without area computation and without BG options

Parameters
InOrigParticlesis the full set of constituent candidates. handed by reference!
JetDefis a fastjet::JetDefinition for the clustering

Definition at line 45 of file JetAnalyzer.cxx.

View newest version in sPHENIX GitHub at line 45 of file JetAnalyzer.cxx

References area_def_bkgd, bkgd_estimator, bkgd_subtractor, CanDoBackground, and jet_def_bkgd.

JetAnalyzer::~JetAnalyzer ( )
inline

Destructor. Take care of all the objects created with new.

Definition at line 126 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 126 of file JetAnalyzer.hh

References area_def_bkgd, bkgd_estimator, bkgd_subtractor, and jet_def_bkgd.

Member Function Documentation

fastjet::JetMedianBackgroundEstimator* JetAnalyzer::GetBackgroundEstimator ( )
inline

Handle to BackgroundEstimator()

Definition at line 162 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 162 of file JetAnalyzer.hh

fastjet::Subtractor * JetAnalyzer::GetBackgroundSubtractor ( )

Background functionality. Currently, the jet definition is hard-coded to fastjet::kt_algorithm, jet_def().R(), and the area definition is computed internally. Expand and modify as needed.

Definition at line 68 of file JetAnalyzer.cxx.

View newest version in sPHENIX GitHub at line 68 of file JetAnalyzer.cxx

References area_def_bkgd, bkgd_estimator, bkgd_subtractor, CanDoBackground, jet_def_bkgd, kt_algorithm, OrigParticles, Acts::IntegrationTest::R, and selector_bkgd.

Referenced by PHAJMaker::process_event().

+ Here is the caller graph for this function:

double JetAnalyzer::phimod2pi ( double  phi)
static

Returns an angle between -pi and pi

Definition at line 185 of file JetAnalyzer.cxx.

View newest version in sPHENIX GitHub at line 185 of file JetAnalyzer.cxx

References ActsTests::PropagationDatasets::phi, and pi.

Referenced by SelectorDijetWorker::terminator().

+ Here is the caller graph for this function:

Member Data Documentation

fastjet::AreaDefinition* JetAnalyzer::area_def_bkgd
private

Background area definiton

Definition at line 99 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 99 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor(), JetAnalyzer(), and ~JetAnalyzer().

fastjet::JetMedianBackgroundEstimator* JetAnalyzer::bkgd_estimator
private

For background subtraction. Pointer for now, so that it can be 0 if unused. Default specs and areas will be supplied.

Definition at line 90 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 90 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor(), JetAnalyzer(), and ~JetAnalyzer().

fastjet::Subtractor* JetAnalyzer::bkgd_subtractor
private

Subtractor

Definition at line 92 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 92 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor(), JetAnalyzer(), and ~JetAnalyzer().

bool JetAnalyzer::CanDoBackground
private

Determined by whether we have an area definition

Definition at line 84 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 84 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor(), and JetAnalyzer().

fastjet::JetDefinition* JetAnalyzer::jet_def_bkgd
private

Background jet definiton

Definition at line 97 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 97 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor(), JetAnalyzer(), and ~JetAnalyzer().

std::vector<fastjet::PseudoJet>& JetAnalyzer::OrigParticles
private

Keep a copy of the original constituents. In principle, they should be accessible via this->jets(), but there is extra space allocated in that vector and I'd rather not muck around with it.

Definition at line 81 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 81 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor().

const double JetAnalyzer::pi = 3.141592653589793238462643383279502884
static

The nth version of this important constant.

Definition at line 162 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 162 of file JetAnalyzer.hh

Referenced by phimod2pi(), and SelectorDijetWorker::terminator().

fastjet::Selector JetAnalyzer::selector_bkgd
private

Background jet cut

Definition at line 95 of file JetAnalyzer.hh.

View newest version in sPHENIX GitHub at line 95 of file JetAnalyzer.hh

Referenced by GetBackgroundSubtractor().


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