CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

CATopJetAlgorithm Class Reference

#include <CATopJetAlgorithm.h>

List of all members.

Public Member Functions

 CATopJetAlgorithm (edm::InputTag mSrc, bool verbose, int algorithm, int useAdjacency, double centralEtaCut, double ptMin, const std::vector< double > &sumEtBins, const std::vector< double > &rBins, const std::vector< double > &ptFracBins, const std::vector< double > &deltarBins, const std::vector< double > &nCellBins, double seedThreshold, bool useMaxTower, double sumEtEtaCut, double etFrac)
void run (const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< fastjet::PseudoJet > &hardjetsOutput, boost::shared_ptr< fastjet::ClusterSequence > &fjClusterSeq)
 Find the ProtoJets from the collection of input Candidates.

Private Member Functions

bool adjacentCells (const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2, const std::vector< fastjet::PseudoJet > &cell_particles, const fastjet::ClusterSequence &theClusterSequence, double nCellMin) const
bool decomposeJet (const fastjet::PseudoJet &theJet, const fastjet::ClusterSequence &theClusterSequence, const std::vector< fastjet::PseudoJet > &cell_particles, double ptHard, double nCellMin, double deltarcut, fastjet::PseudoJet &ja, fastjet::PseudoJet &jb, std::vector< fastjet::PseudoJet > &leftovers) const

Private Attributes

int algorithm_
double centralEtaCut_
std::vector< double > deltarBins_
double etFrac_
std::string jetType_
edm::InputTag mSrc_
std::vector< double > nCellBins_
std::vector< double > ptFracBins_
double ptMin_
std::vector< double > rBins_
double seedThreshold_
std::vector< double > sumEtBins_
double sumEtEtaCut_
int useAdjacency_
bool useMaxTower_
bool verbose_

Detailed Description

Definition at line 44 of file CATopJetAlgorithm.h.


Constructor & Destructor Documentation

CATopJetAlgorithm::CATopJetAlgorithm ( edm::InputTag  mSrc,
bool  verbose,
int  algorithm,
int  useAdjacency,
double  centralEtaCut,
double  ptMin,
const std::vector< double > &  sumEtBins,
const std::vector< double > &  rBins,
const std::vector< double > &  ptFracBins,
const std::vector< double > &  deltarBins,
const std::vector< double > &  nCellBins,
double  seedThreshold,
bool  useMaxTower,
double  sumEtEtaCut,
double  etFrac 
) [inline]

Constructor

Definition at line 48 of file CATopJetAlgorithm.h.

                                   :
    mSrc_          (mSrc          ),
    verbose_       (verbose       ),
    algorithm_     (algorithm     ),
    useAdjacency_  (useAdjacency  ),
    centralEtaCut_ (centralEtaCut ), 
    ptMin_         (ptMin         ),         
    sumEtBins_     (sumEtBins     ),        
    rBins_         (rBins         ),         
    ptFracBins_    (ptFracBins    ),  
    deltarBins_    (deltarBins    ),
    nCellBins_     (nCellBins     ),
    seedThreshold_ (seedThreshold ), 
    useMaxTower_   (useMaxTower   ),
    sumEtEtaCut_   (sumEtEtaCut   ),   
    etFrac_        (etFrac        )

      { }

Member Function Documentation

bool CATopJetAlgorithm::adjacentCells ( const fastjet::PseudoJet &  jet1,
const fastjet::PseudoJet &  jet2,
const std::vector< fastjet::PseudoJet > &  cell_particles,
const fastjet::ClusterSequence &  theClusterSequence,
double  nCellMin 
) const [private]

Definition at line 269 of file CATopJetAlgorithm.cc.

References abs, and reco::deltaPhi().

                                                                                                  {
        
        
        double eta1 = jet1.rapidity();
        double phi1 = jet1.phi();
        double eta2 = jet2.rapidity();
        double phi2 = jet2.phi();
        
        double deta = abs(eta2 - eta1) / 0.087;
        double dphi = fabs( reco::deltaPhi(phi2,phi1) ) / 0.087;
        
        return ( ( deta + dphi ) <= nCellMin );
}
bool CATopJetAlgorithm::decomposeJet ( const fastjet::PseudoJet &  theJet,
const fastjet::ClusterSequence &  theClusterSequence,
const std::vector< fastjet::PseudoJet > &  cell_particles,
double  ptHard,
double  nCellMin,
double  deltarcut,
fastjet::PseudoJet &  ja,
fastjet::PseudoJet &  jb,
std::vector< fastjet::PseudoJet > &  leftovers 
) const [private]

Adjacency Requirement ///

Pt Fraction Requirement ///

Definition at line 291 of file CATopJetAlgorithm.cc.

References gather_cfg::cout, SiPixelRawToDigiRegional_cfi::deltaPhi, deltaR(), and j.

                                                                                                                     {
        
        bool goodBreak;
        fastjet::PseudoJet j = theJet;
        double InputObjectPt = j.perp();
        if ( verbose_ )cout<<"Input Object Pt = "<<InputObjectPt<<endl;
        if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
        leftovers.clear();
        if ( verbose_ )cout<<"start while loop"<<endl;
        
        while (1) {                                                      // watch out for infinite loop!
                goodBreak = theClusterSequence.has_parents(j,ja,jb);
                if (!goodBreak){
                        if ( verbose_ )cout<<"bad break. this is one cell. can't decluster anymore."<<endl;
                        break;         // this is one cell, can't decluster anymore
                }
                
                if ( verbose_ )cout<<"good break. ja Pt = "<<ja.perp()<<" jb Pt = "<<jb.perp()<<endl;
                
                
                // check if clusters are adjacent using a constant deltar adjacency.
                double clusters_deltar=fabs(ja.eta()-jb.eta())+fabs(deltaPhi(ja.phi(),jb.phi()));
                
                if ( verbose_  && useAdjacency_ ==1)cout<<"clusters_deltar = "<<clusters_deltar<<endl;
                if ( verbose_  && useAdjacency_ ==1)cout<<"deltar cut = "<<deltarcut<<endl;
                
                if ( useAdjacency_==1 && clusters_deltar < deltarcut){
                        if ( verbose_ )cout<<"clusters too close. consant adj. break."<<endl;
                        break;
                } 
                
                // Check if clusters are adjacent using a DeltaR adjacency which is a function of pT.
                double clusters_deltaR=deltaR( ja.eta(), ja.phi(), jb.eta(), jb.phi() );
                
                if ( verbose_  && useAdjacency_ ==2)cout<<"clusters_deltaR = "<<clusters_deltaR<<endl;
                if ( verbose_  && useAdjacency_ ==2)cout<<"0.4-0.0004*InputObjectPt = "<<0.4-0.0004*InputObjectPt<<endl;
                
                if ( useAdjacency_==2 && clusters_deltaR < 0.4-0.0004*InputObjectPt)
                {
                        if ( verbose_ )cout<<"clusters too close. modified adj. break."<<endl;
                        break;
                } 

                // Check if clusters are adjacent in the calorimeter. 
                if ( useAdjacency_==3 &&  adjacentCells(ja,jb,cell_particles,theClusterSequence,nCellMin) ){                  
                        if ( verbose_ )cout<<"clusters too close in the calorimeter. calorimeter adj. break."<<endl;
                        break;         // the clusters are "adjacent" in the calorimeter => shouldn't have decomposed
                }
                                
                if ( verbose_ )cout<<"clusters pass distance cut"<<endl;
                
                
                if ( verbose_ )cout<<"ptHard = "<<ptHard<<endl;
                
                if (ja.perp() < ptHard && jb.perp() < ptHard){
                        if ( verbose_ )cout<<"two soft clusters. dead end"<<endl;
                        break;         // broke into two soft clusters, dead end
                }
                
                if (ja.perp() > ptHard && jb.perp() > ptHard){
                        if ( verbose_ )cout<<"two hard clusters. done"<<endl;
                        return true;   // broke into two hard clusters, we're done!
                }
                
                else if (ja.perp() > jb.perp()) {                              // broke into one hard and one soft, ditch the soft one and try again
                        if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = ja"<<endl; 
                        j = ja;
                        vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(jb);
                        leftovers.insert(leftovers.end(),particles.begin(),particles.end());
                }
                else {
                        if ( verbose_ )cout<<"ja hard jb soft. try to split hard. j = jb"<<endl; 
                        j = jb;
                        vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(ja);
                        leftovers.insert(leftovers.end(),particles.begin(),particles.end());
                }
        }
        
        if ( verbose_ )cout<<"did not decluster."<<endl;  // did not decluster into hard subjets
        
        ja.reset(0,0,0,0);
        jb.reset(0,0,0,0);
        leftovers.clear();
        return false;
}
void CATopJetAlgorithm::run ( const std::vector< fastjet::PseudoJet > &  cell_particles,
std::vector< fastjet::PseudoJet > &  hardjetsOutput,
boost::shared_ptr< fastjet::ClusterSequence > &  fjClusterSeq 
)

Find the ProtoJets from the collection of input Candidates.

Referenced by cms::CATopJetProducer::runAlgorithm().


Member Data Documentation

Definition at line 91 of file CATopJetAlgorithm.h.

Definition at line 97 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::deltarBins_ [private]

Definition at line 102 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::etFrac_ [private]

Definition at line 108 of file CATopJetAlgorithm.h.

std::string CATopJetAlgorithm::jetType_ [private]

Definition at line 109 of file CATopJetAlgorithm.h.

Definition at line 89 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::nCellBins_ [private]

Definition at line 103 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::ptFracBins_ [private]

Definition at line 101 of file CATopJetAlgorithm.h.

double CATopJetAlgorithm::ptMin_ [private]

Definition at line 98 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::rBins_ [private]

Definition at line 100 of file CATopJetAlgorithm.h.

Definition at line 105 of file CATopJetAlgorithm.h.

std::vector<double> CATopJetAlgorithm::sumEtBins_ [private]

Definition at line 99 of file CATopJetAlgorithm.h.

Definition at line 107 of file CATopJetAlgorithm.h.

Definition at line 92 of file CATopJetAlgorithm.h.

Definition at line 106 of file CATopJetAlgorithm.h.

Definition at line 90 of file CATopJetAlgorithm.h.