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 
52  if (regionPSet.exists("measurementTrackerName"))
53  {
54  // FIXME: when next time altering the configuration of this
55  // class, please change the types of the following parameters:
56  // - whereToUseMeasurementTracker to at least int32 or to a string
57  // corresponding to the UseMeasurementTracker enumeration
58  // - measurementTrackerName to InputTag
59  if (regionPSet.exists("whereToUseMeasurementTracker"))
62  token_measurementTracker = iC.consumes<MeasurementTrackerEvent>(regionPSet.getParameter<std::string>("measurementTrackerName"));
63  }
64  }
65 
67 
68 
69  virtual std::vector<TrackingRegion* > regions(const edm::Event& e, const edm::EventSetup& es) const
70  {
71  std::vector<TrackingRegion* > result;
72 
73  // use beam spot to pick up the origin
75  e.getByToken( token_beamSpot, bsHandle);
76  if(!bsHandle.isValid()) return result;
77  const reco::BeamSpot & bs = *bsHandle;
78  GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
79 
80  // pick up the candidate objects of interest
82  e.getByToken( token_jet, objects );
83  size_t n_objects = objects->size();
84  if (n_objects == 0) return result;
85 
86  const MeasurementTrackerEvent *measurementTracker = nullptr;
90  measurementTracker = hmte.product();
91  }
92 
93  // create maximum JetMaxN tracking regions in directions of
94  // highest pt jets that are above threshold and are within allowed eta
95  // (we expect that jet collection was sorted in decreasing pt order)
96  int n_regions = 0;
97  for (size_t i =0; i < n_objects && n_regions < m_jetMaxN; ++i)
98  {
99  const reco::Candidate & jet = (*objects)[i];
100  if ( jet.pt() < m_jetMinPt || std::abs(jet.eta()) > m_jetMaxEta ) continue;
101 
102  GlobalVector direction(jet.momentum().x(), jet.momentum().y(), jet.momentum().z());
103 
105  direction,
106  origin,
107  m_ptMin,
110  m_deltaEta,
111  m_deltaPhi,
113  m_precise,
114  measurementTracker,
116  );
117  result.push_back(etaphiRegion);
118  ++n_regions;
119  }
120  //std::cout<<"nregions = "<<n_regions<<std::endl;
121  return result;
122  }
123 
124 private:
125 
126  float m_ptMin;
129  float m_deltaEta;
130  float m_deltaPhi;
132  float m_jetMinPt;
133  float m_jetMaxEta;
139  bool m_precise;
140 };
141 
142 #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:446
virtual float eta() const =0
momentum pseudorapidity
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< MeasurementTrackerEvent > token_measurementTracker
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
static UseMeasurementTracker doubleToUseMeasurementTracker(double value)
T const * product() const
Definition: Handle.h:81
double y0() const
y coordinate
Definition: BeamSpot.h:66
bool isUninitialized() const
Definition: EDGetToken.h:71
TrackingRegionsFromBeamSpotAndL2Tau(const edm::ParameterSet &conf, edm::ConsumesCollector &&iC)
double x0() const
x coordinate
Definition: BeamSpot.h:64