CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonTrackingRegionBuilder.cc
Go to the documentation of this file.
1 
11 
12 //-------------------------------
13 // Collaborating Class Headers --
14 //-------------------------------
15 
17 
21 
24 
28 
29 //
30 // constructor
31 //
33 
34  // Adjust errors on Eta, Phi, Z
35  theNsigmaEta = par.getParameter<double>("Rescale_eta");
36  theNsigmaPhi = par.getParameter<double>("Rescale_phi");
37  theNsigmaDz = par.getParameter<double>("Rescale_Dz");
38 
39  // Upper limits parameters
40  theEtaRegionPar1 = par.getParameter<double>("EtaR_UpperLimit_Par1");
41  theEtaRegionPar2 = par.getParameter<double>("EtaR_UpperLimit_Par2");
42  thePhiRegionPar1 = par.getParameter<double>("PhiR_UpperLimit_Par1");
43  thePhiRegionPar2 = par.getParameter<double>("PhiR_UpperLimit_Par2");
44 
45  // Flag to switch to use Vertices instead of BeamSpot
46  useVertex = par.getParameter<bool>("UseVertex");
47 
48  // Flag to use fixed limits for Eta, Phi, Z, pT
49  useFixedZ = par.getParameter<bool>("Z_fixed");
50  useFixedPt = par.getParameter<bool>("Pt_fixed");
51  useFixedPhi = par.getParameter<bool>("Phi_fixed");
52  useFixedEta = par.getParameter<bool>("Eta_fixed");
53 
54  // Minimum value for pT
55  thePtMin = par.getParameter<double>("Pt_min");
56 
57  // Minimum value for Phi
58  thePhiMin = par.getParameter<double>("Phi_min");
59 
60  // Minimum value for Eta
61  theEtaMin = par.getParameter<double>("Eta_min");
62 
63  // The static region size along the Z direction
64  theHalfZ = par.getParameter<double>("DeltaZ");
65 
66  // The transverse distance of the region from the BS/PV
67  theDeltaR = par.getParameter<double>("DeltaR");
68 
69  // The static region size in Eta
70  theDeltaEta = par.getParameter<double>("DeltaEta");
71 
72  // The static region size in Phi
73  theDeltaPhi = par.getParameter<double>("DeltaPhi");
74 
75  // Maximum number of regions to build when looping over Muons
76  theMaxRegions = par.getParameter<int>("maxRegions");
77 
78  // Flag to use precise??
79  thePrecise = par.getParameter<bool>("precise");
80 
81  // perigee reference point ToDo: Check this
85  }
86 
87  // Vertex collection and Beam Spot
90 
91  // Input muon collection
93 }
94 
95 
96 //
97 // Member function to be compatible with TrackingRegionProducerFactory: create many ROI for many tracks
98 //
99 std::vector<std::unique_ptr<TrackingRegion>> MuonTrackingRegionBuilder::regions(const edm::Event& ev, const edm::EventSetup& es) const {
100 
101  std::vector<std::unique_ptr<TrackingRegion>> result;
102 
104  ev.getByToken(inputCollectionToken, tracks);
105 
106  int nRegions = 0;
107  for (auto it = tracks->cbegin(), ed = tracks->cend(); it != ed && nRegions < theMaxRegions; ++it) {
108  result.push_back(region(*it,ev));
109  nRegions++;
110  }
111 
112  return result;
113 
114 }
115 
116 
117 //
118 // Call region on Track from TrackRef
119 //
120 std::unique_ptr<RectangularEtaPhiTrackingRegion> MuonTrackingRegionBuilder::region(const reco::TrackRef& track) const {
121  return region(*track);
122 }
123 
124 
125 //
126 // ToDo: Not sure if this is needed?
127 //
129  theEvent = &event;
130 }
131 
132 
133 //
134 // Main member function called to create the ROI
135 //
136 std::unique_ptr<RectangularEtaPhiTrackingRegion> MuonTrackingRegionBuilder::region(const reco::Track& staTrack, const edm::Event& ev) const {
137 
138  // get track momentum/direction at vertex
139  const math::XYZVector& mom = staTrack.momentum();
140  GlobalVector dirVector(mom.x(),mom.y(),mom.z());
141  double pt = staTrack.pt();
142 
143  // Fix for StandAlone tracks with low momentum
144  const math::XYZVector& innerMomentum = staTrack.innerMomentum();
145  GlobalVector forSmallMomentum(innerMomentum.x(),innerMomentum.y(),innerMomentum.z());
146  if ( staTrack.p() <= 1.5 ) {
147  pt = std::abs(forSmallMomentum.perp());
148  }
149 
150  // initial vertex position - in the following it is replaced with beamspot/vertexing
151  GlobalPoint vertexPos(0.0,0.0,0.0);
152  // standard 15.9, if useVertex than use error from vertex
153  double deltaZ = theHalfZ;
154 
155  // retrieve beam spot information
157  bool bsHandleFlag = ev.getByToken(beamSpotToken, bs);
158 
159  // check the validity, otherwise vertexing
160  if ( bsHandleFlag && bs.isValid() && !useVertex ) {
161  vertexPos = GlobalPoint(bs->x0(), bs->y0(), bs->z0());
162  deltaZ = useFixedZ ? theHalfZ : bs->sigmaZ() * theNsigmaDz;
163  } else {
164  // get originZPos from list of reconstructed vertices (first or all)
166  bool vtxHandleFlag = ev.getByToken(vertexCollectionToken, vertexCollection);
167  // check if there exists at least one reconstructed vertex
168  if ( vtxHandleFlag && !vertexCollection->empty() ) {
169  // use the first vertex in the collection and assume it is the primary event vertex
170  reco::VertexCollection::const_iterator vtx = vertexCollection->begin();
171  if (!vtx->isFake() && vtx->isValid() ) {
172  vertexPos = GlobalPoint(vtx->x(),vtx->y(),vtx->z());
173  deltaZ = useFixedZ ? theHalfZ : vtx->zError() * theNsigmaDz;
174  }
175  }
176  }
177 
178  // inizialize to the maximum possible value
179  double deta = 0.4;
180  double dphi = 0.6;
181 
182  // evaluate the dynamical region if possible
183  deta = theNsigmaEta*(staTrack.etaError());
184  dphi = theNsigmaPhi*(staTrack.phiError());
185 
186  // Region_Parametrizations to take into account possible L2 error matrix inconsistencies
187  double region_dEta = 0;
188  double region_dPhi = 0;
189  double eta = 0;
190  double phi = 0;
191 
192  // eta, pt parametrization from MC study (circa 2009?)
193  if ( pt <= 10. ) {
194  // angular coefficients
195  float acoeff_Phi = (thePhiRegionPar2 - thePhiRegionPar1)/5;
196  float acoeff_Eta = (theEtaRegionPar2 - theEtaRegionPar1)/5;
197 
198  eta = theEtaRegionPar1 + (acoeff_Eta)*(pt-5.);
199  phi = thePhiRegionPar1 + (acoeff_Phi)*(pt-5.) ;
200  }
201  // parametrization 2nd bin in pt from MC study
202  if ( pt > 10. && pt < 100. ) {
203  eta = theEtaRegionPar2;
204  phi = thePhiRegionPar2;
205  }
206  // parametrization 3rd bin in pt from MC study
207  if ( pt >= 100. ) {
208  // angular coefficients
209  float acoeff_Phi = (thePhiRegionPar1 - thePhiRegionPar2)/900;
210  float acoeff_Eta = (theEtaRegionPar1 - theEtaRegionPar2)/900;
211 
212  eta = theEtaRegionPar2 + (acoeff_Eta)*(pt-100.);
213  phi = thePhiRegionPar2 + (acoeff_Phi)*(pt-100.);
214  }
215 
216  double region_dPhi1 = std::min(phi,dphi);
217  double region_dEta1 = std::min(eta,deta);
218 
219  // decide to use either a parametrization or a dynamical region
220  region_dPhi = useFixedPhi ? theDeltaPhi : std::max(thePhiMin,region_dPhi1);
221  region_dEta = useFixedEta ? theDeltaEta : std::max(theEtaMin,region_dEta1);
222 
223  float deltaR = theDeltaR;
224  double minPt = useFixedPt ? thePtMin : std::max(thePtMin,pt*0.6);
225 
226 
231  measurementTracker = hmte.product();
232  }
233 
234  auto region = std::make_unique<RectangularEtaPhiTrackingRegion>(dirVector, vertexPos,
235  minPt, deltaR,
236  deltaZ, region_dEta, region_dPhi,
237  theOnDemand,
238  thePrecise,
240 
241  LogDebug("MuonTrackingRegionBuilder")<<"the region parameters are:\n"
242  <<"\n dirVector: "<<dirVector
243  <<"\n vertexPos: "<<vertexPos
244  <<"\n minPt: "<<minPt
245  <<"\n deltaR:"<<deltaR
246  <<"\n deltaZ:"<<deltaZ
247  <<"\n region_dEta:"<<region_dEta
248  <<"\n region_dPhi:"<<region_dPhi
249  <<"\n on demand parameter: "<<static_cast<int>(theOnDemand);
250 
251  return region;
252 
253 }
254 
256  {
258  desc.add<double>("EtaR_UpperLimit_Par1",0.25);
259  desc.add<double>("DeltaR",0.2);
260  desc.add<edm::InputTag>("beamSpot",edm::InputTag(""));
261  desc.add<int>("OnDemand",-1);
262  desc.add<edm::InputTag>("vertexCollection",edm::InputTag(""));
263  desc.add<double>("Rescale_phi",3.0);
264  desc.add<bool>("Eta_fixed",false);
265  desc.add<double>("Rescale_eta",3.0);
266  desc.add<double>("PhiR_UpperLimit_Par2",0.2);
267  desc.add<double>("Eta_min",0.05);
268  desc.add<bool>("Phi_fixed",false);
269  desc.add<double>("Phi_min",0.05);
270  desc.add<double>("PhiR_UpperLimit_Par1",0.6);
271  desc.add<double>("EtaR_UpperLimit_Par2",0.15);
272  desc.add<edm::InputTag>("MeasurementTrackerName",edm::InputTag(""));
273  desc.add<bool>("UseVertex",false);
274  desc.add<double>("Rescale_Dz",3.0);
275  desc.add<bool>("Pt_fixed",false);
276  desc.add<bool>("Z_fixed",true);
277  desc.add<double>("Pt_min",1.5);
278  desc.add<double>("DeltaZ",15.9);
279  desc.add<double>("DeltaEta",0.2);
280  desc.add<double>("DeltaPhi",0.2);
281  desc.add<int>("maxRegions",1);
282  desc.add<bool>("precise",true);
283  desc.add<edm::InputTag>("input",edm::InputTag(""));
284  descriptions.add("MuonTrackingRegionBuilder",desc);
285  }
286  {
288  desc.add<double>("EtaR_UpperLimit_Par1",0.25);
289  desc.add<double>("DeltaR",0.2);
290  desc.add<edm::InputTag>("beamSpot",edm::InputTag("hltOnlineBeamSpot"));
291  desc.add<int>("OnDemand",-1);
292  desc.add<edm::InputTag>("vertexCollection",edm::InputTag("pixelVertices"));
293  desc.add<double>("Rescale_phi",3.0);
294  desc.add<bool>("Eta_fixed",false);
295  desc.add<double>("Rescale_eta",3.0);
296  desc.add<double>("PhiR_UpperLimit_Par2",0.2);
297  desc.add<double>("Eta_min",0.05);
298  desc.add<bool>("Phi_fixed",false);
299  desc.add<double>("Phi_min",0.05);
300  desc.add<double>("PhiR_UpperLimit_Par1",0.6);
301  desc.add<double>("EtaR_UpperLimit_Par2",0.15);
302  desc.add<edm::InputTag>("MeasurementTrackerName",edm::InputTag("hltESPMeasurementTracker"));
303  desc.add<bool>("UseVertex",false);
304  desc.add<double>("Rescale_Dz",3.0);
305  desc.add<bool>("Pt_fixed",false);
306  desc.add<bool>("Z_fixed",true);
307  desc.add<double>("Pt_min",1.5);
308  desc.add<double>("DeltaZ",15.9);
309  desc.add<double>("DeltaEta",0.2);
310  desc.add<double>("DeltaPhi",0.2);
311  desc.add<int>("maxRegions",1);
312  desc.add<bool>("precise",true);
313  desc.add<edm::InputTag>("input",edm::InputTag("hltL2Muons","UpdatedAtVtx"));
314  descriptions.add("hltMuonTrackingRegionBuilder",desc);
315  }
316  descriptions.setComment("Build a TrackingRegion around a standalone muon. Options to define region around beamspot or primary vertex and dynamic regions are included.");
317 }
#define LogDebug(id)
double p() const
momentum vector magnitude
Definition: TrackBase.h:602
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
static void fillDescriptions(edm::ParameterSetDescription &descriptions)
Add Fill Descriptions.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken
virtual void setEvent(const edm::Event &)
Pass the Event to the algo at each event.
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
double etaError() const
error on eta
Definition: TrackBase.h:771
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool ev
const Vector & momentum() const
track momentum vector
Definition: TrackBase.h:662
tuple vertexCollection
std::unique_ptr< RectangularEtaPhiTrackingRegion > region(const reco::TrackRef &) const
Define tracking region.
void setComment(std::string const &value)
void build(const edm::ParameterSet &, edm::ConsumesCollector &)
double pt() const
track transverse momentum
Definition: TrackBase.h:608
tuple result
Definition: query.py:137
double phiError() const
error on phi
Definition: TrackBase.h:777
RectangularEtaPhiTrackingRegion::UseMeasurementTracker theOnDemand
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:75
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
edm::EDGetTokenT< reco::TrackCollection > inputCollectionToken
virtual std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &, const edm::EventSetup &) const override
Create Region of Interest.
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
edm::EDGetTokenT< reco::VertexCollection > vertexCollectionToken
bool isUninitialized() const
Definition: EDGetToken.h:73
static UseMeasurementTracker intToUseMeasurementTracker(int value)
edm::EDGetTokenT< MeasurementTrackerEvent > theMeasurementTrackerToken