CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoJets/JetAlgorithms/interface/CATopJetAlgorithm.h

Go to the documentation of this file.
00001 #ifndef RecoJets_JetAlgorithms_CATopJetAlgorithm_h
00002 #define RecoJets_JetAlgorithms_CATopJetAlgorithm_h
00003 
00004 
00005 
00006 /* *********************************************************
00007  * \class CATopJetAlgorithm
00008  * Jet producer to produce top jets using the C-A algorithm to break
00009  * jets into subjets as described here:
00010  * "Top-tagging: A Method for Identifying Boosted Hadronic Tops"
00011  * David E. Kaplan, Keith Rehermann, Matthew D. Schwartz, Brock Tweedie
00012  * arXiv:0806.0848v1 [hep-ph] 
00013  *
00014  ************************************************************/
00015 
00016 #include <vector>
00017 #include <list>
00018 #include <functional>
00019 #include <TMath.h>
00020 #include <iostream>
00021 
00022 #include "DataFormats/Candidate/interface/Candidate.h"
00023 #include "DataFormats/CaloTowers/interface/CaloTower.h"
00024 #include "DataFormats/CaloTowers/interface/CaloTowerFwd.h"
00025 #include "DataFormats/JetReco/interface/CaloJet.h"
00026 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00027 #include "DataFormats/JetReco/interface/BasicJet.h"
00028 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
00029 #include "DataFormats/Math/interface/deltaR.h"
00030 #include "RecoJets/JetAlgorithms/interface/JetAlgoHelper.h"
00031 #include "RecoJets/JetAlgorithms/interface/CompoundPseudoJet.h"
00032 #include "DataFormats/Candidate/interface/LeafCandidate.h"
00033 #include "FWCore/Framework/interface/Event.h"
00034 
00035 #include "RecoJets/JetAlgorithms/interface/CompoundPseudoJet.h"
00036 
00037 
00038 #include <fastjet/JetDefinition.hh>
00039 #include <fastjet/PseudoJet.hh>
00040 #include <fastjet/ClusterSequence.hh>
00041 #include <fastjet/GhostedAreaSpec.hh>
00042 #include <fastjet/ClusterSequenceArea.hh>
00043 
00044 class CATopJetAlgorithm{
00045  public:
00048   CATopJetAlgorithm(edm::InputTag mSrc,
00049                     bool   verbose,
00050                     int algorithm,
00051                     int useAdjacency,
00052                     double centralEtaCut,            
00053                     double ptMin,                    
00054                     const std::vector<double> & sumEtBins,      
00055                     const std::vector<double> & rBins,       
00056                     const std::vector<double> & ptFracBins,  
00057                     const std::vector<double> & deltarBins,
00058                     const std::vector<double> & nCellBins,
00059                     double seedThreshold,            
00060                     bool   useMaxTower,
00061                     double sumEtEtaCut,
00062                     double etFrac) :
00063     mSrc_          (mSrc          ),
00064     verbose_       (verbose       ),
00065     algorithm_     (algorithm     ),
00066     useAdjacency_  (useAdjacency  ),
00067     centralEtaCut_ (centralEtaCut ), 
00068     ptMin_         (ptMin         ),         
00069     sumEtBins_     (sumEtBins     ),        
00070     rBins_         (rBins         ),         
00071     ptFracBins_    (ptFracBins    ),  
00072     deltarBins_    (deltarBins    ),
00073     nCellBins_     (nCellBins     ),
00074     seedThreshold_ (seedThreshold ), 
00075     useMaxTower_   (useMaxTower   ),
00076     sumEtEtaCut_   (sumEtEtaCut   ),   
00077     etFrac_        (etFrac        )
00078 
00079       { }
00080 
00082     void run( const std::vector<fastjet::PseudoJet> & cell_particles, 
00083               std::vector<fastjet::PseudoJet> & hardjetsOutput ,
00084               boost::shared_ptr<fastjet::ClusterSequence> & fjClusterSeq
00085               );
00086 
00087  private:
00088 
00089   edm::InputTag       mSrc_;                            //<! calo tower input source
00090   bool                            verbose_;                 //<!
00091   int                 algorithm_;                       //<! 0 = KT, 1 = CA, 2 = anti-KT
00092   int                 useAdjacency_;                    //<! choose adjacency requirement:
00093                                                                                                 //<!    0 = no adjacency
00094                                                                                                 //<!    1 = deltar adjacency 
00095                                                                                                 //<!    2 = modified adjacency
00096                                                                                                 //<!    3 = calotower neirest neigbor based adjacency (untested)
00097   double              centralEtaCut_;                   //<! eta for defining "central" jets                                     
00098   double              ptMin_;                           //<! lower pt cut on which jets to reco                                  
00099   std::vector<double> sumEtBins_;                       //<! sumEt bins over which cuts vary. vector={bin 0 lower bound, bin 1 lower bound, ...}   
00100   std::vector<double> rBins_;                           //<! Jet distance paramter R. R values depend on sumEt bins.         
00101   std::vector<double> ptFracBins_;                      //<! deltap = fraction of full jet pt for a subjet to be consider "hard". 
00102   std::vector<double> deltarBins_;              //<! Applicable only if useAdjacency=1. deltar adjacency values for each sumEtBin
00103   std::vector<double> nCellBins_;                       //<! Applicable only if useAdjacency=3. number of cells apart for two subjets to be considered "independent"
00104   // NOT USED:
00105   double              seedThreshold_;                   //<! calo tower seed threshold - NOT USED 
00106   bool                useMaxTower_;                     //<! use max tower for jet adjacency criterion, false is to use the centroid - NOT USED
00107   double              sumEtEtaCut_;                     //<! eta for event SumEt - NOT USED                                 
00108   double              etFrac_;                          //<! fraction of event sumEt / 2 for a jet to be considered "hard" - NOT USED 
00109   std::string         jetType_;                         //<! CaloJets or GenJets - NOT USED
00110 
00111 
00112   // Decide if the two jets are in adjacent cells    
00113   bool adjacentCells(const fastjet::PseudoJet & jet1, const fastjet::PseudoJet & jet2, 
00114                      const std::vector<fastjet::PseudoJet> & cell_particles,
00115                      const fastjet::ClusterSequence & theClusterSequence,
00116                      double nCellMin ) const;
00117 
00118 
00119   // Attempt to break up one "hard" jet into two "soft" jets
00120 
00121   bool decomposeJet(const fastjet::PseudoJet & theJet, 
00122                     const fastjet::ClusterSequence & theClusterSequence, 
00123                     const std::vector<fastjet::PseudoJet> & cell_particles,
00124                     double ptHard, double nCellMin, double deltarcut,
00125                     fastjet::PseudoJet & ja, fastjet::PseudoJet & jb, 
00126                     std::vector<fastjet::PseudoJet> & leftovers) const;
00127 
00128 };
00129 
00130 
00131 
00132 #endif