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)
 
virtual ~TrackingRegionsFromBeamSpotAndL2Tau ()
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Private Attributes

edm::InputTag m_beamSpotTag
 
float m_deltaEta
 
float m_deltaPhi
 
float m_jetMaxEta
 
int m_jetMaxN
 
float m_jetMinPt
 
edm::InputTag m_jetSrc
 
std::string m_measurementTrackerName
 
float m_originHalfLength
 
float m_originRadius
 
bool m_precise
 
float m_ptMin
 
bool m_searchOpt
 
float m_whereToUseMeasurementTracker
 

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)
inlineexplicit

Definition at line 28 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), m_beamSpotTag, m_deltaEta, m_deltaPhi, m_jetMaxEta, m_jetMaxN, m_jetMinPt, m_jetSrc, m_measurementTrackerName, m_originHalfLength, m_originRadius, m_precise, m_ptMin, m_searchOpt, and m_whereToUseMeasurementTracker.

29  {
30  edm::LogInfo ("TrackingRegionsFromBeamSpotAndL2Tau") << "Enter the TrackingRegionsFromBeamSpotAndL2Tau";
31 
32  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
33 
34  m_ptMin = regionPSet.getParameter<double>("ptMin");
35  m_originRadius = regionPSet.getParameter<double>("originRadius");
36  m_originHalfLength = regionPSet.getParameter<double>("originHalfLength");
37  m_deltaEta = regionPSet.getParameter<double>("deltaEta");
38  m_deltaPhi = regionPSet.getParameter<double>("deltaPhi");
39  m_jetSrc = regionPSet.getParameter<edm::InputTag>("JetSrc");
40  m_jetMinPt = regionPSet.getParameter<double>("JetMinPt");
41  m_jetMaxEta = regionPSet.getParameter<double>("JetMaxEta");
42  m_jetMaxN = regionPSet.getParameter<int>("JetMaxN");
43  m_beamSpotTag = regionPSet.getParameter<edm::InputTag>("beamSpot");
44  m_precise = regionPSet.getParameter<bool>("precise");
45 
46  if (regionPSet.exists("searchOpt")) m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
47  else m_searchOpt = false;
48 
51  if (regionPSet.exists("measurementTrackerName"))
52  {
53  m_measurementTrackerName = regionPSet.getParameter<std::string>("measurementTrackerName");
54  if (regionPSet.exists("whereToUseMeasurementTracker"))
55  m_whereToUseMeasurementTracker = regionPSet.getParameter<double>("whereToUseMeasurementTracker");
56  }
57  }
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual TrackingRegionsFromBeamSpotAndL2Tau::~TrackingRegionsFromBeamSpotAndL2Tau ( )
inlinevirtual

Definition at line 59 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

59 {}

Member Function Documentation

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

Implements TrackingRegionProducer.

Definition at line 62 of file TrackingRegionsFromBeamSpotAndL2Tau.h.

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

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

Member Data Documentation

edm::InputTag TrackingRegionsFromBeamSpotAndL2Tau::m_beamSpotTag
private
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
edm::InputTag TrackingRegionsFromBeamSpotAndL2Tau::m_jetSrc
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