CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HFJetShowerShape Class Reference

#include <HFJetShowerShape/HFJetShowerShape/plugins/HFJetShowerShape.cc>

Inheritance diagram for HFJetShowerShape:
edm::stream::EDProducer<>

Public Member Functions

 HFJetShowerShape (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 
template<typename T >
void putInEvent (const std::string &, const edm::Handle< edm::View< reco::Jet >> &, std::vector< T >, edm::Event &) const
 Function to put product into event. More...
 

Private Attributes

const double hfTowerEtaWidth_
 
const double hfTowerPhiWidth_
 
const double jetEtaThreshold_
 
const double jetPtThreshold_
 
const double jetReferenceRadius_
 
const edm::EDGetTokenT< edm::View< reco::Jet > > jets_token_
 
const double offsetPerPU_
 
const double stripPtThreshold_
 
const double vertexRecoEffcy_
 
const edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_token_
 
const double widthPtThreshold_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 43 of file HFJetShowerShape.cc.

Constructor & Destructor Documentation

◆ HFJetShowerShape()

HFJetShowerShape::HFJetShowerShape ( const edm::ParameterSet iConfig)
explicit

Definition at line 78 of file HFJetShowerShape.cc.

79  : jets_token_(consumes<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("jets"))),
80  vertices_token_(consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("vertices"))),
81  jetPtThreshold_(iConfig.getParameter<double>("jetPtThreshold")),
82  jetEtaThreshold_(iConfig.getParameter<double>("jetEtaThreshold")),
83  hfTowerEtaWidth_(iConfig.getParameter<double>("hfTowerEtaWidth")),
84  hfTowerPhiWidth_(iConfig.getParameter<double>("hfTowerPhiWidth")),
85  vertexRecoEffcy_(iConfig.getParameter<double>("vertexRecoEffcy")),
86  offsetPerPU_(iConfig.getParameter<double>("offsetPerPU")),
87  jetReferenceRadius_(iConfig.getParameter<double>("jetReferenceRadius")),
88  stripPtThreshold_(iConfig.getParameter<double>("stripPtThreshold")),
89  widthPtThreshold_(iConfig.getParameter<double>("widthPtThreshold")) {
90  produces<edm::ValueMap<float>>("sigmaEtaEta");
91  produces<edm::ValueMap<float>>("sigmaPhiPhi");
92  produces<edm::ValueMap<int>>("centralEtaStripSize");
93  produces<edm::ValueMap<int>>("adjacentEtaStripsSize");
94 }
const double jetReferenceRadius_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const double vertexRecoEffcy_
const double widthPtThreshold_
const double hfTowerPhiWidth_
const double jetPtThreshold_
const edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_token_
const edm::EDGetTokenT< edm::View< reco::Jet > > jets_token_
const double offsetPerPU_
const double stripPtThreshold_
const double jetEtaThreshold_
const double hfTowerEtaWidth_

Member Function Documentation

◆ fillDescriptions()

void HFJetShowerShape::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 212 of file HFJetShowerShape.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, and HLT_2022v15_cff::InputTag.

212  {
214  desc.add<edm::InputTag>("jets", edm::InputTag("ak4PFJetsCHS"));
215  desc.add<edm::InputTag>("vertices", edm::InputTag("offlineSlimmedPrimaryVertices"));
216  desc.add<double>("jetPtThreshold", 25.);
217  desc.add<double>("jetEtaThreshold", 2.9);
218  desc.add<double>("hfTowerEtaWidth", 0.175);
219  desc.add<double>("hfTowerPhiWidth", 0.175);
220  desc.add<double>("vertexRecoEffcy", 0.7);
221  desc.add<double>("offsetPerPU", 0.4);
222  desc.add<double>("jetReferenceRadius", 0.4);
223  desc.add<double>("stripPtThreshold", 10.);
224  desc.add<double>("widthPtThreshold", 3.);
225  descriptions.add("hfJetShowerShape", desc);
226 }
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

void HFJetShowerShape::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 101 of file HFJetShowerShape.cc.

References funct::abs(), SiPixelRawToDigiRegional_cfi::deltaPhi, reco::Candidate::eta(), hfTowerEtaWidth_, hfTowerPhiWidth_, mps_fire::i, iEvent, metsig::jet, jetEtaThreshold_, jetPtThreshold_, jetReferenceRadius_, jets_token_, M_PI, nPV, offsetPerPU_, reco::Candidate::pdgId(), reco::Candidate::phi(), reco::Candidate::pt(), putInEvent(), mathSSE::sqrt(), stripPtThreshold_, vertexRecoEffcy_, vertices_token_, mps_merge::weight, and widthPtThreshold_.

101  {
102  using namespace edm;
103 
104  //Jets
105  auto theJets = iEvent.getHandle(jets_token_);
106 
107  //Vertices
108  int nPV = iEvent.get(vertices_token_).size();
109 
110  //Products
111  std::vector<float> v_sigmaEtaEta, v_sigmaPhiPhi;
112  v_sigmaEtaEta.reserve(theJets->size());
113  v_sigmaPhiPhi.reserve(theJets->size());
114  std::vector<int> v_size_CentralEtaStrip, v_size_AdjacentEtaStrips;
115  v_size_CentralEtaStrip.reserve(theJets->size());
116  v_size_AdjacentEtaStrips.reserve(theJets->size());
117 
118  //Et offset for HF PF candidates
121 
122  for (auto const& jet : *theJets) {
123  double pt_jet = jet.pt();
124  double eta_jet = jet.eta();
125 
126  //If central or low pt jets, fill with dummy variables
127  if (pt_jet <= jetPtThreshold_ || std::abs(eta_jet) <= jetEtaThreshold_) {
128  v_sigmaEtaEta.push_back(-1.);
129  v_sigmaPhiPhi.push_back(-1.);
130  v_size_CentralEtaStrip.push_back(0);
131  v_size_AdjacentEtaStrips.push_back(0);
132  } else {
133  //First loop over PF candidates to compute some global variables needed for shower shape calculations
134  double sumptPFcands = 0.;
135 
136  for (unsigned i = 0; i < jet.numberOfSourceCandidatePtrs(); ++i) {
137  const reco::Candidate* icand = jet.sourceCandidatePtr(i).get();
138  //Only look at pdgId =1,2 (HF PF cands)
139  if (std::abs(icand->pdgId()) > 2)
140  continue;
141  double pt_PUsub = icand->pt() - puoffset;
142  if (pt_PUsub < widthPtThreshold_)
143  continue;
144  sumptPFcands += pt_PUsub;
145  }
146 
147  //Second loop to compute the various shower shape variables
148  int size_CentralEtaStrip(0.), size_AdjacentEtaStrips(0.);
149  double sigmaEtaEtaSq(0.), sigmaPhiPhiSq(0.);
150  double sumweightsPFcands = 0;
151  for (unsigned i = 0; i < jet.numberOfSourceCandidatePtrs(); ++i) {
152  const reco::Candidate* icand = jet.sourceCandidatePtr(i).get();
153  //Only look at pdgId =1,2 (HF PF cands)
154  if (std::abs(icand->pdgId()) > 2)
155  continue;
156 
157  double deta = std::abs(icand->eta() - jet.eta());
158  double dphi = std::abs(deltaPhi(icand->phi(), jet.phi()));
159  double pt_PUsub = icand->pt() - puoffset;
160 
161  //This is simply the size of the central eta strip and the adjacent strips
162  if (pt_PUsub >= stripPtThreshold_) {
163  if (dphi <= hfTowerPhiWidth_ * 0.5)
164  size_CentralEtaStrip++;
165  else if (dphi <= hfTowerPhiWidth_ * 1.5)
166  size_AdjacentEtaStrips++;
167  }
168 
169  //Now computing sigmaEtaEta/PhiPhi
170  if (pt_PUsub >= widthPtThreshold_ && sumptPFcands > 0) {
171  double weight = pt_PUsub / sumptPFcands;
172  sigmaEtaEtaSq += deta * deta * weight;
173  sigmaPhiPhiSq += dphi * dphi * weight;
174  sumweightsPFcands += weight;
175  }
176  }
177 
178  v_size_CentralEtaStrip.push_back(size_CentralEtaStrip);
179  v_size_AdjacentEtaStrips.push_back(size_AdjacentEtaStrips);
180 
181  if (sumweightsPFcands > 0 && sigmaEtaEtaSq > 0 && sigmaPhiPhiSq > 0) {
182  v_sigmaEtaEta.push_back(sqrt(sigmaEtaEtaSq / sumweightsPFcands));
183  v_sigmaPhiPhi.push_back(sqrt(sigmaPhiPhiSq / sumweightsPFcands));
184  } else {
185  v_sigmaEtaEta.push_back(-1.);
186  v_sigmaPhiPhi.push_back(-1.);
187  }
188 
189  } //End loop over jets
190  }
191 
192  putInEvent("sigmaEtaEta", theJets, v_sigmaEtaEta, iEvent);
193  putInEvent("sigmaPhiPhi", theJets, v_sigmaPhiPhi, iEvent);
194  putInEvent("centralEtaStripSize", theJets, v_size_CentralEtaStrip, iEvent);
195  putInEvent("adjacentEtaStripsSize", theJets, v_size_AdjacentEtaStrips, iEvent);
196 }
const double jetReferenceRadius_
const double vertexRecoEffcy_
virtual double pt() const =0
transverse momentum
Definition: weight.py:1
const double widthPtThreshold_
void putInEvent(const std::string &, const edm::Handle< edm::View< reco::Jet >> &, std::vector< T >, edm::Event &) const
Function to put product into event.
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
const double hfTowerPhiWidth_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const double jetPtThreshold_
const edm::EDGetTokenT< std::vector< reco::Vertex > > vertices_token_
const edm::EDGetTokenT< edm::View< reco::Jet > > jets_token_
const double offsetPerPU_
#define M_PI
virtual int pdgId() const =0
PDG identifier.
HLT enums.
const double stripPtThreshold_
const double jetEtaThreshold_
const double hfTowerEtaWidth_
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity

◆ putInEvent()

template<typename T >
void HFJetShowerShape::putInEvent ( const std::string &  name,
const edm::Handle< edm::View< reco::Jet >> &  jets,
std::vector< T product,
edm::Event iEvent 
) const
private

Function to put product into event.

Definition at line 200 of file HFJetShowerShape.cc.

References trigObjTnPSource_cfi::filler, iEvent, PDWG_EXODelayedJetMET_cff::jets, eostools::move(), Skims_PA_cff::name, and MillePedeFileConverter_cfg::out.

Referenced by produce().

203  {
204  auto out = std::make_unique<edm::ValueMap<T>>();
206  filler.insert(jets, product.begin(), product.end());
207  filler.fill();
208  iEvent.put(std::move(out), name);
209 }
int iEvent
Definition: GenABIO.cc:224
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ hfTowerEtaWidth_

const double HFJetShowerShape::hfTowerEtaWidth_
private

Definition at line 60 of file HFJetShowerShape.cc.

Referenced by produce().

◆ hfTowerPhiWidth_

const double HFJetShowerShape::hfTowerPhiWidth_
private

Definition at line 60 of file HFJetShowerShape.cc.

Referenced by produce().

◆ jetEtaThreshold_

const double HFJetShowerShape::jetEtaThreshold_
private

Definition at line 58 of file HFJetShowerShape.cc.

Referenced by produce().

◆ jetPtThreshold_

const double HFJetShowerShape::jetPtThreshold_
private

Definition at line 58 of file HFJetShowerShape.cc.

Referenced by produce().

◆ jetReferenceRadius_

const double HFJetShowerShape::jetReferenceRadius_
private

Definition at line 62 of file HFJetShowerShape.cc.

Referenced by produce().

◆ jets_token_

const edm::EDGetTokenT<edm::View<reco::Jet> > HFJetShowerShape::jets_token_
private

Definition at line 54 of file HFJetShowerShape.cc.

Referenced by produce().

◆ offsetPerPU_

const double HFJetShowerShape::offsetPerPU_
private

Definition at line 62 of file HFJetShowerShape.cc.

Referenced by produce().

◆ stripPtThreshold_

const double HFJetShowerShape::stripPtThreshold_
private

Definition at line 64 of file HFJetShowerShape.cc.

Referenced by produce().

◆ vertexRecoEffcy_

const double HFJetShowerShape::vertexRecoEffcy_
private

Definition at line 62 of file HFJetShowerShape.cc.

Referenced by produce().

◆ vertices_token_

const edm::EDGetTokenT<std::vector<reco::Vertex> > HFJetShowerShape::vertices_token_
private

Definition at line 55 of file HFJetShowerShape.cc.

Referenced by produce().

◆ widthPtThreshold_

const double HFJetShowerShape::widthPtThreshold_
private

Definition at line 64 of file HFJetShowerShape.cc.

Referenced by produce().