Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:24

0001 #ifndef RecoJets_JetAlgorithms_CATopJetAlgorithm_h
0002 #define RecoJets_JetAlgorithms_CATopJetAlgorithm_h
0003 
0004 /* *********************************************************
0005  * \class CATopJetAlgorithm
0006  * Jet producer to produce top jets using the C-A algorithm to break
0007  * jets into subjets as described here:
0008  * "Top-tagging: A Method for Identifying Boosted Hadronic Tops"
0009  * David E. Kaplan, Keith Rehermann, Matthew D. Schwartz, Brock Tweedie
0010  * arXiv:0806.0848v1 [hep-ph] 
0011  *
0012  ************************************************************/
0013 
0014 #include <vector>
0015 #include <list>
0016 #include <functional>
0017 #include <TMath.h>
0018 #include <iostream>
0019 
0020 #include "DataFormats/Candidate/interface/Candidate.h"
0021 #include "DataFormats/CaloTowers/interface/CaloTower.h"
0022 #include "DataFormats/CaloTowers/interface/CaloTowerDefs.h"
0023 #include "DataFormats/JetReco/interface/CaloJet.h"
0024 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0025 #include "DataFormats/JetReco/interface/BasicJet.h"
0026 #include "DataFormats/JetReco/interface/BasicJetCollection.h"
0027 #include "DataFormats/Math/interface/deltaR.h"
0028 #include "RecoJets/JetAlgorithms/interface/JetAlgoHelper.h"
0029 #include "RecoJets/JetAlgorithms/interface/CompoundPseudoJet.h"
0030 #include "DataFormats/Candidate/interface/LeafCandidate.h"
0031 #include "FWCore/Framework/interface/Event.h"
0032 
0033 #include "RecoJets/JetAlgorithms/interface/CompoundPseudoJet.h"
0034 
0035 #include <fastjet/JetDefinition.hh>
0036 #include <fastjet/PseudoJet.hh>
0037 #include <fastjet/ClusterSequence.hh>
0038 #include <fastjet/GhostedAreaSpec.hh>
0039 #include <fastjet/ClusterSequenceArea.hh>
0040 
0041 class CATopJetAlgorithm {
0042 public:
0043   /** Constructor
0044   */
0045   CATopJetAlgorithm(const edm::InputTag& mSrc,
0046                     bool verbose,
0047                     int algorithm,
0048                     int useAdjacency,
0049                     double centralEtaCut,
0050                     double ptMin,
0051                     const std::vector<double>& sumEtBins,
0052                     const std::vector<double>& rBins,
0053                     const std::vector<double>& ptFracBins,
0054                     const std::vector<double>& deltarBins,
0055                     const std::vector<double>& nCellBins,
0056                     double seedThreshold,
0057                     bool useMaxTower,
0058                     double sumEtEtaCut,
0059                     double etFrac)
0060       : mSrc_(mSrc),
0061         verbose_(verbose),
0062         algorithm_(algorithm),
0063         useAdjacency_(useAdjacency),
0064         centralEtaCut_(centralEtaCut),
0065         ptMin_(ptMin),
0066         sumEtBins_(sumEtBins),
0067         rBins_(rBins),
0068         ptFracBins_(ptFracBins),
0069         deltarBins_(deltarBins),
0070         nCellBins_(nCellBins),
0071         seedThreshold_(seedThreshold),
0072         useMaxTower_(useMaxTower),
0073         sumEtEtaCut_(sumEtEtaCut),
0074         etFrac_(etFrac)
0075 
0076   {}
0077 
0078   /// Find the ProtoJets from the collection of input Candidates.
0079   void run(const std::vector<fastjet::PseudoJet>& cell_particles,
0080            std::vector<fastjet::PseudoJet>& hardjetsOutput,
0081            std::shared_ptr<fastjet::ClusterSequence>& fjClusterSeq);
0082 
0083 private:
0084   edm::InputTag mSrc_;    //<! calo tower input source
0085   bool verbose_;          //<!
0086   int algorithm_;         //<! 0 = KT, 1 = CA, 2 = anti-KT
0087   int useAdjacency_;      //<! choose adjacency requirement:
0088                           //<!  0 = no adjacency
0089                           //<!  1 = deltar adjacency
0090                           //<!  2 = modified adjacency
0091                           //<!  3 = calotower neirest neigbor based adjacency (untested)
0092   double centralEtaCut_;  //<! eta for defining "central" jets
0093   double ptMin_;          //<! lower pt cut on which jets to reco
0094   std::vector<double>
0095       sumEtBins_;              //<! sumEt bins over which cuts vary. vector={bin 0 lower bound, bin 1 lower bound, ...}
0096   std::vector<double> rBins_;  //<! Jet distance paramter R. R values depend on sumEt bins.
0097   std::vector<double> ptFracBins_;  //<! deltap = fraction of full jet pt for a subjet to be consider "hard".
0098   std::vector<double> deltarBins_;  //<! Applicable only if useAdjacency=1. deltar adjacency values for each sumEtBin
0099   std::vector<double>
0100       nCellBins_;  //<! Applicable only if useAdjacency=3. number of cells apart for two subjets to be considered "independent"
0101   // NOT USED:
0102   double seedThreshold_;  //<! calo tower seed threshold - NOT USED
0103   bool useMaxTower_;      //<! use max tower for jet adjacency criterion, false is to use the centroid - NOT USED
0104   double sumEtEtaCut_;    //<! eta for event SumEt - NOT USED
0105   double etFrac_;         //<! fraction of event sumEt / 2 for a jet to be considered "hard" - NOT USED
0106   std::string jetType_;   //<! CaloJets or GenJets - NOT USED
0107 
0108   // Decide if the two jets are in adjacent cells
0109   bool adjacentCells(const fastjet::PseudoJet& jet1,
0110                      const fastjet::PseudoJet& jet2,
0111                      const std::vector<fastjet::PseudoJet>& cell_particles,
0112                      const fastjet::ClusterSequence& theClusterSequence,
0113                      double nCellMin) const;
0114 
0115   // Attempt to break up one "hard" jet into two "soft" jets
0116 
0117   bool decomposeJet(const fastjet::PseudoJet& theJet,
0118                     const fastjet::ClusterSequence& theClusterSequence,
0119                     const std::vector<fastjet::PseudoJet>& cell_particles,
0120                     double ptHard,
0121                     double nCellMin,
0122                     double deltarcut,
0123                     fastjet::PseudoJet& ja,
0124                     fastjet::PseudoJet& jb,
0125                     std::vector<fastjet::PseudoJet>& leftovers) const;
0126 };
0127 
0128 #endif