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 
12 
13 //---------------
14 // C++ Headers --
15 //---------------
16 
17 
18 //-------------------------------
19 // Collaborating Class Headers --
20 //-------------------------------
21 
32 
33 using namespace std;
34 
35 //
36 // constructor
37 //
38 
39 void MuonTrackingRegionBuilder::init(const MuonServiceProxy* service) { theService= service;}
41 {
42  build(par, iC);
43 }
45  // vertex Collection and Beam Spot
46  theBeamSpotTag = par.getParameter<edm::InputTag>("beamSpot");
47  theVertexCollTag = par.getParameter<edm::InputTag>("vertexCollection");
48 
49  // parmeters
50  theNsigmaEta = par.getParameter<double>("Rescale_eta");
51  theNsigmaPhi = par.getParameter<double>("Rescale_phi");
52  theNsigmaDz = par.getParameter<double>("Rescale_Dz");
53  theTkEscapePt = par.getParameter<double>("EscapePt");
54 
55  // upper limits parameters
56  theEtaRegionPar1 = par.getParameter<double>("EtaR_UpperLimit_Par1");
57  theEtaRegionPar2 = par.getParameter<double>("EtaR_UpperLimit_Par2");
58  thePhiRegionPar1 = par.getParameter<double>("PhiR_UpperLimit_Par1");
59  thePhiRegionPar2 = par.getParameter<double>("PhiR_UpperLimit_Par2");
60 
61  useVertex = par.getParameter<bool>("UseVertex");
62  useFixedRegion = par.getParameter<bool>("UseFixedRegion");
63 
64  // fixed limits
65  thePhiFixed = par.getParameter<double>("Phi_fixed");
66  theEtaFixed = par.getParameter<double>("Eta_fixed");
67 
68  thePhiMin = par.getParameter<double>("Phi_min");
69  theEtaMin = par.getParameter<double>("Eta_min");
70  theHalfZ = par.getParameter<double>("DeltaZ_Region");
71  theDeltaR = par.getParameter<double>("DeltaR");
72 
73  // perigee reference point
74  // FIXME: when next time altering the configuration of this
75  // class, please change the types of the following parameters:
76  // - OnDemand to at least int32 or to a string
77  // corresponding to the UseMeasurementTracker enumeration
78  // - MeasurementTrackerName to InputTag
81  theMeasurementTrackerToken = iC.consumes<MeasurementTrackerEvent>(par.getParameter<std::string>("MeasurementTrackerName"));
82  }
83 
84  bsHandleToken=iC.consumes<reco::BeamSpot>(theBeamSpotTag);
85  vertexCollectionToken=iC.consumes<reco::VertexCollection>(theVertexCollTag);
86 }
87 
88 
89 //
90 //
91 //
94 
95  return region(*track);
96 
97 }
98 
99 
100 //
101 //
102 //
104 
105  theEvent = &event;
106 
107 }
108 
109 
110 //
111 //
112 //
115 
116  // get the free trajectory state of the muon updated at vertex
117  TSCPBuilderNoMaterial tscpBuilder;
118 
119  FreeTrajectoryState muFTS = trajectoryStateTransform::initialFreeState(staTrack,&*theService->magneticField());
120 
121  LogDebug("MuonTrackingRegionBuilder")<<"from state: "<<muFTS;
122 
123  // get track direction at vertex
124  GlobalVector dirVector(muFTS.momentum());
125 
126  // get track momentum
127  const math::XYZVector& mo = staTrack.innerMomentum();
128  GlobalVector mom(mo.x(),mo.y(),mo.z());
129  if ( staTrack.p() > 1.5 ) {
130  mom = dirVector;
131  }
132 
133  // initial vertex position - in the following it is replaced with beam spot/vertexing
134  GlobalPoint vertexPos(0.0,0.0,0.0);
135  double deltaZatVTX = 0.0;
136 
137  // retrieve beam spot information
139  bool bsHandleFlag = theEvent->getByToken(bsHandleToken, bsHandle);
140  // check the validity, otherwise vertexing
141  // inizialization of BS
142 
143  if ( bsHandleFlag && !useVertex ) {
144  const reco::BeamSpot& bs = *bsHandle;
145  vertexPos = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
146  } else {
147  // get originZPos from list of reconstructed vertices (first or all)
149  bool vtxHandleFlag = theEvent->getByToken(vertexCollectionToken, vertexCollection);
150  // check if there exists at least one reconstructed vertex
151  if ( vtxHandleFlag && !vertexCollection->empty() ) {
152  // use the first vertex in the collection and assume it is the primary event vertex
153  reco::VertexCollection::const_iterator vtx = vertexCollection->begin();
154  vertexPos = GlobalPoint(vtx->x(),vtx->y(),vtx->z());
155  // delta Z from vertex error
156  deltaZatVTX = vtx->zError() * theNsigmaDz;
157  }
158  }
159 
160 
161  // inizialize to the maximum possible value to avoit
162  // problems with TSCBL
163 
164  double deta = 0.4;
165  double dphi = 0.6;
166 
167  // take into account the correct beanspot rotation
168  if ( bsHandleFlag ) {
169 
170  const reco::BeamSpot& bs = *bsHandle;
171 
172  TSCBLBuilderNoMaterial tscblBuilder;
173  TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(muFTS,bs);
174 
175 
176  // evaluate the dynamical region if possible
177  if(tscbl.isValid()){
178 
180  GlobalVector pTrack = tscbl.trackStateAtPCA().momentum();
181 
182  // calculate deltaEta from deltaTheta
183  double deltaTheta = trackPerigeeErrors.thetaError();
184  double theta = pTrack.theta();
185  double sin_theta = sin(theta);
186 
187  // get dEta and dPhi
188  deta = theNsigmaEta*(1/fabs(sin_theta))*deltaTheta;
189  dphi = theNsigmaPhi*(trackPerigeeErrors.phiError());
190 
191  }
192  }
193 
194  /* Region_Parametrizations to take into account possible
195  L2 error matrix inconsistencies. Detailed Explanation in TWIKI
196  page.
197  */
198  double region_dEta = 0;
199  double region_dPhi = 0;
200  double eta=0; double phi=0;
201 
202  // eta, pt parametrization from MC study
203  float pt = fabs(mom.perp());
204 
205  if ( pt <= 10. ) {
206  // angular coefficients
207  float acoeff_Phi = (thePhiRegionPar2 - thePhiRegionPar1)/5;
208  float acoeff_Eta = (theEtaRegionPar2 - theEtaRegionPar1)/5;
209 
210  eta = theEtaRegionPar1 + (acoeff_Eta)*(mom.perp()-5.);
211  phi = thePhiRegionPar1 + (acoeff_Phi)*(mom.perp()-5.) ;
212  }
213  // parametrization 2nd bin in pt from MC study
214  if ( pt > 10. && pt < 100. ) {
215  eta = theEtaRegionPar2;
216  phi = thePhiRegionPar2;
217  }
218  // parametrization 3rd bin in pt from MC study
219  if ( pt >= 100. ) {
220  // angular coefficients
221  float acoeff_Phi = (thePhiRegionPar1 - thePhiRegionPar2)/900;
222  float acoeff_Eta = (theEtaRegionPar1 - theEtaRegionPar2)/900;
223 
224  eta = theEtaRegionPar2 + (acoeff_Eta)*(mom.perp()-100.);
225  phi = thePhiRegionPar2 + (acoeff_Phi)*(mom.perp()-100.);
226  }
227 
228  // decide to use either a parametrization or a dynamical region
229  double region_dPhi1 = min(phi,dphi);
230  double region_dEta1 = min(eta,deta);
231 
232  // minimum size
233  region_dPhi = max(thePhiMin,region_dPhi1);
234  region_dEta = max(theEtaMin,region_dEta1);
235 
236  float deltaZ = 0.0;
237  // standard 15.9 is useVertex than region from vertexing
238  if ( useVertex ) {
239  deltaZ = max(theHalfZ,deltaZatVTX);
240  } else {
241  deltaZ = theHalfZ;
242  }
243 
244  float deltaR = theDeltaR;
245  double minPt = max(theTkEscapePt,mom.perp()*0.6);
246 
247  RectangularEtaPhiTrackingRegion* region = 0;
248 
249  if (useFixedRegion) {
250  region_dEta = theEtaFixed;
251  region_dPhi = thePhiFixed;
252  }
253 
255  if(!theMeasurementTrackerToken.isUninitialized()) {
257  theEvent->getByToken(theMeasurementTrackerToken, hmte);
258  measurementTracker = hmte.product();
259  }
260 
261  region = new RectangularEtaPhiTrackingRegion(dirVector, vertexPos,
262  minPt, deltaR,
263  deltaZ, region_dEta, region_dPhi,
264  theOnDemand,
265  true,/*default in the header*/
266  measurementTracker);
267 
268  LogDebug("MuonTrackingRegionBuilder")<<"the region parameters are:\n"
269  <<"\n dirVector: "<<dirVector
270  <<"\n vertexPos: "<<vertexPos
271  <<"\n minPt: "<<minPt
272  <<"\n deltaR:"<<deltaR
273  <<"\n deltaZ:"<<deltaZ
274  <<"\n region_dEta:"<<region_dEta
275  <<"\n region_dPhi:"<<region_dPhi
276  <<"\n on demand parameter: "<<static_cast<int>(theOnDemand);
277 
278 
279  return region;
280 
281 }
#define LogDebug(id)
RectangularEtaPhiTrackingRegion * region(const reco::TrackRef &) const
define tracking region
double p() const
momentum vector magnitude
Definition: TrackBase.h:663
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
MuonTrackingRegionBuilder(const edm::ParameterSet &, edm::ConsumesCollector &)
constructor
T getParameter(std::string const &) const
double z0() const
z coordinate
Definition: BeamSpot.h:68
virtual void setEvent(const edm::Event &)
pass the Event to the algo at each event
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Geom::Theta< T > theta() const
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
T eta() const
tuple vertexCollection
void build(const edm::ParameterSet &, edm::ConsumesCollector &)
T min(T a, T b)
Definition: MathUtil.h:58
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
GlobalVector momentum() const
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
static UseMeasurementTracker doubleToUseMeasurementTracker(double value)
void init(const MuonServiceProxy *)
T const * product() const
Definition: Handle.h:81
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
double y0() const
y coordinate
Definition: BeamSpot.h:66
PerigeeTrajectoryError ftsToPerigeeError(const FTS &originalFTS)
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field, bool withErr=true)
double x0() const
x coordinate
Definition: BeamSpot.h:64
Definition: DDAxes.h:10