CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
TrackingRegionsFromBeamSpotAndL2Tau Class Reference

#include <TrackingRegionsFromBeamSpotAndL2Tau.h>

Inheritance diagram for TrackingRegionsFromBeamSpotAndL2Tau:
TrackingRegionProducer

Public Member Functions

virtual std::vector
< TrackingRegion * > 
regions (const edm::Event &e, const edm::EventSetup &es) const
 
 TrackingRegionsFromBeamSpotAndL2Tau (const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
 
virtual ~TrackingRegionsFromBeamSpotAndL2Tau ()
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Private Attributes

float m_deltaEta
 
float m_deltaPhi
 
float m_jetMaxEta
 
int m_jetMaxN
 
float m_jetMinPt
 
std::string m_measurementTrackerName
 
float m_originHalfLength
 
float m_originRadius
 
bool m_precise
 
float m_ptMin
 
bool m_searchOpt
 
float m_whereToUseMeasurementTracker
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT
< reco::CandidateView
token_jet
 

Detailed Description

class TrackingRegionsFromBeamSpotAndL2Tau plugin for creating eta-phi TrackingRegions in directions of L2 taus

Definition at line 24 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

Constructor & Destructor Documentation

TrackingRegionsFromBeamSpotAndL2Tau::TrackingRegionsFromBeamSpotAndL2Tau ( const edm::ParameterSet conf,
edm::ConsumesCollector &&  iC 
)
inlineexplicit

Definition at line 28 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), m_deltaEta, m_deltaPhi, m_jetMaxEta, m_jetMaxN, m_jetMinPt, m_measurementTrackerName, m_originHalfLength, m_originRadius, m_precise, m_ptMin, m_searchOpt, m_whereToUseMeasurementTracker, AlCaHLTBitMon_QueryRunRegistry::string, token_beamSpot, and token_jet.

30  {
31  edm::LogInfo ("TrackingRegionsFromBeamSpotAndL2Tau") << "Enter the TrackingRegionsFromBeamSpotAndL2Tau";
32 
33  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
34 
35  m_ptMin = regionPSet.getParameter<double>("ptMin");
36  m_originRadius = regionPSet.getParameter<double>("originRadius");
37  m_originHalfLength = regionPSet.getParameter<double>("originHalfLength");
38  m_deltaEta = regionPSet.getParameter<double>("deltaEta");
39  m_deltaPhi = regionPSet.getParameter<double>("deltaPhi");
41  m_jetMinPt = regionPSet.getParameter<double>("JetMinPt");
42  m_jetMaxEta = regionPSet.getParameter<double>("JetMaxEta");
43  m_jetMaxN = regionPSet.getParameter<int>("JetMaxN");
44  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
45  m_precise = regionPSet.getParameter<bool>("precise");
46 
47  if (regionPSet.exists("searchOpt")) m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
48  else m_searchOpt = false;
49 
52  if (regionPSet.exists("measurementTrackerName"))
53  {
54  m_measurementTrackerName = regionPSet.getParameter<std::string>("measurementTrackerName");
55  if (regionPSet.exists("whereToUseMeasurementTracker"))
56  m_whereToUseMeasurementTracker = regionPSet.getParameter<double>("whereToUseMeasurementTracker");
57  }
58  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::CandidateView > token_jet
virtual TrackingRegionsFromBeamSpotAndL2Tau::~TrackingRegionsFromBeamSpotAndL2Tau ( )
inlinevirtual

Definition at line 60 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

60 {}

Member Function Documentation

virtual std::vector<TrackingRegion* > TrackingRegionsFromBeamSpotAndL2Tau::regions ( const edm::Event e,
const edm::EventSetup es 
) const
inlinevirtual

Implements TrackingRegionProducer.

Definition at line 63 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

References funct::abs(), reco::Candidate::eta(), edm::Event::getByToken(), i, edm::HandleBase::isValid(), metsig::jet, m_deltaEta, m_deltaPhi, m_jetMaxEta, m_jetMaxN, m_jetMinPt, m_measurementTrackerName, m_originHalfLength, m_originRadius, m_precise, m_ptMin, m_searchOpt, m_whereToUseMeasurementTracker, reco::Candidate::momentum(), reco::Candidate::pt(), query::result, token_beamSpot, token_jet, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

64  {
65  std::vector<TrackingRegion* > result;
66 
67  // use beam spot to pick up the origin
69  e.getByToken( token_beamSpot, bsHandle);
70  if(!bsHandle.isValid()) return result;
71  const reco::BeamSpot & bs = *bsHandle;
72  GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
73 
74  // pick up the candidate objects of interest
76  e.getByToken( token_jet, objects );
77  size_t n_objects = objects->size();
78  if (n_objects == 0) return result;
79 
80  // create maximum JetMaxN tracking regions in directions of
81  // highest pt jets that are above threshold and are within allowed eta
82  // (we expect that jet collection was sorted in decreasing pt order)
83  int n_regions = 0;
84  for (size_t i =0; i < n_objects && n_regions < m_jetMaxN; ++i)
85  {
86  const reco::Candidate & jet = (*objects)[i];
87  if ( jet.pt() < m_jetMinPt || std::abs(jet.eta()) > m_jetMaxEta ) continue;
88 
89  GlobalVector direction(jet.momentum().x(), jet.momentum().y(), jet.momentum().z());
90 
92  direction,
93  origin,
94  m_ptMin,
97  m_deltaEta,
98  m_deltaPhi,
100  m_precise,
103  );
104  result.push_back(etaphiRegion);
105  ++n_regions;
106  }
107  //std::cout<<"nregions = "<<n_regions<<std::endl;
108  return result;
109  }
double z0() const
z coordinate
Definition: BeamSpot.h:68
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual float eta() const =0
momentum pseudorapidity
edm::EDGetTokenT< reco::CandidateView > token_jet
virtual Vector momentum() const =0
spatial momentum vector
virtual float pt() const =0
transverse momentum
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:76
double y0() const
y coordinate
Definition: BeamSpot.h:66
double x0() const
x coordinate
Definition: BeamSpot.h:64

Member Data Documentation

float TrackingRegionsFromBeamSpotAndL2Tau::m_deltaEta
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_deltaPhi
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_jetMaxEta
private
int TrackingRegionsFromBeamSpotAndL2Tau::m_jetMaxN
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_jetMinPt
private
std::string TrackingRegionsFromBeamSpotAndL2Tau::m_measurementTrackerName
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_originHalfLength
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_originRadius
private
bool TrackingRegionsFromBeamSpotAndL2Tau::m_precise
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_ptMin
private
bool TrackingRegionsFromBeamSpotAndL2Tau::m_searchOpt
private
float TrackingRegionsFromBeamSpotAndL2Tau::m_whereToUseMeasurementTracker
private
edm::EDGetTokenT<reco::BeamSpot> TrackingRegionsFromBeamSpotAndL2Tau::token_beamSpot
private
edm::EDGetTokenT<reco::CandidateView> TrackingRegionsFromBeamSpotAndL2Tau::token_jet
private