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