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 
11 
20 
25 {
26 public:
27 
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");
40  token_jet = iC.consumes<reco::CandidateView>(regionPSet.getParameter<edm::InputTag>("JetSrc"));
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  }
59 
61 
62 
63  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
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  }
110 
111 private:
112 
113  float m_ptMin;
116  float m_deltaEta;
117  float m_deltaPhi;
119  float m_jetMinPt;
120  float m_jetMaxEta;
126  bool m_precise;
127 };
128 
129 #endif
T getParameter(std::string const &) const
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
bool exists(std::string const &parameterName) const
checks if a parameter exists
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
virtual std::vector< TrackingRegion * > regions(const edm::Event &e, const edm::EventSetup &es) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:76
tuple conf
Definition: dbtoconf.py:185
double y0() const
y coordinate
Definition: BeamSpot.h:66
TrackingRegionsFromBeamSpotAndL2Tau(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
double x0() const
x coordinate
Definition: BeamSpot.h:64