CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingRegionsFromBeamSpotAndL2Tau.h
Go to the documentation of this file.
1 #ifndef TrackingRegionsFromBeamSpotAndL2Tau_h
2 #define TrackingRegionsFromBeamSpotAndL2Tau_h
3 
4 //
5 // Class: TrackingRegionsFromBeamSpotAndL2Tau
6 //
7 
8 
12 
21 
26 {
27 public:
28 
31  {
32  edm::LogInfo ("TrackingRegionsFromBeamSpotAndL2Tau") << "Enter the TrackingRegionsFromBeamSpotAndL2Tau";
33 
34  edm::ParameterSet regionPSet = conf.getParameter<edm::ParameterSet>("RegionPSet");
35 
36  m_ptMin = regionPSet.getParameter<double>("ptMin");
37  m_originRadius = regionPSet.getParameter<double>("originRadius");
38  m_originHalfLength = regionPSet.getParameter<double>("originHalfLength");
39  m_deltaEta = regionPSet.getParameter<double>("deltaEta");
40  m_deltaPhi = regionPSet.getParameter<double>("deltaPhi");
41  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
42  m_jetMinPt = regionPSet.getParameter<double>("JetMinPt");
43  m_jetMaxEta = regionPSet.getParameter<double>("JetMaxEta");
44  m_jetMaxN = regionPSet.getParameter<int>("JetMaxN");
45  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
46  m_precise = regionPSet.getParameter<bool>("precise");
47 
48  if (regionPSet.exists("searchOpt")) m_searchOpt = regionPSet.getParameter<bool>("searchOpt");
49  else m_searchOpt = false;
50 
53  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<edm::InputTag>("measurementTrackerName"));
54  }
55  }
56 
58 
59 
60  virtual std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& e, const edm::EventSetup& es) const override
61  {
62  std::vector<std::unique_ptr<TrackingRegion> > result;
63 
64  // use beam spot to pick up the origin
66  e.getByToken( token_beamSpot, bsHandle);
67  if(!bsHandle.isValid()) return result;
68  const reco::BeamSpot & bs = *bsHandle;
69  GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
70 
71  // pick up the candidate objects of interest
73  e.getByToken( token_jet, objects );
74  size_t n_objects = objects->size();
75  if (n_objects == 0) return result;
76 
81  measurementTracker = hmte.product();
82  }
83 
84  // create maximum JetMaxN tracking regions in directions of
85  // highest pt jets that are above threshold and are within allowed eta
86  // (we expect that jet collection was sorted in decreasing pt order)
87  int n_regions = 0;
88  for (size_t i =0; i < n_objects && n_regions < m_jetMaxN; ++i)
89  {
90  const reco::Candidate & jet = (*objects)[i];
91  if ( jet.pt() < m_jetMinPt || std::abs(jet.eta()) > m_jetMaxEta ) continue;
92 
93  GlobalVector direction(jet.momentum().x(), jet.momentum().y(), jet.momentum().z());
94 
95  result.push_back(std::make_unique<RectangularEtaPhiTrackingRegion>(
96  direction,
97  origin,
98  m_ptMin,
101  m_deltaEta,
102  m_deltaPhi,
104  m_precise,
105  measurementTracker,
107  ));
108  ++n_regions;
109  }
110  //std::cout<<"nregions = "<<n_regions<<std::endl;
111  return result;
112  }
113 
114 private:
115 
116  float m_ptMin;
119  float m_deltaEta;
120  float m_deltaPhi;
122  float m_jetMinPt;
123  float m_jetMaxEta;
129  bool m_precise;
130 };
131 
132 #endif
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
int i
Definition: DBlmapReader.cc:9
RectangularEtaPhiTrackingRegion::UseMeasurementTracker m_whereToUseMeasurementTracker
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
virtual double pt() const =0
transverse momentum
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
static UseMeasurementTracker stringToUseMeasurementTracker(const std::string &name)
edm::EDGetTokenT< reco::CandidateView > token_jet
tuple result
Definition: mps_fire.py:84
virtual Vector momentum() const =0
spatial momentum vector
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &e, const edm::EventSetup &es) const override
bool isValid() const
Definition: HandleBase.h:75
T const * product() const
Definition: Handle.h:81
double y0() const
y coordinate
Definition: BeamSpot.h:66
bool isUninitialized() const
Definition: EDGetToken.h:73
TrackingRegionsFromBeamSpotAndL2Tau(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
virtual double eta() const =0
momentum pseudorapidity
double x0() const
x coordinate
Definition: BeamSpot.h:64