CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTTimingExtractor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MuonIdentification
4 // Class: DTTimingExtractor
5 //
11 //
12 // Original Author: Traczyk Piotr
13 // Created: Thu Oct 11 15:01:28 CEST 2007
14 // $Id: DTTimingExtractor.cc,v 1.15 2011/02/24 15:41:53 farrell3 Exp $
15 //
16 //
17 
19 
20 
21 // user include files
24 
27 
29 
31 
40 
44 
52 
57 
58 
59 // system include files
60 #include <memory>
61 #include <vector>
62 #include <string>
63 #include <iostream>
64 
65 namespace edm {
66  class ParameterSet;
67  class EventSetup;
68  class InputTag;
69 }
70 
71 class MuonServiceProxy;
72 
73 //
74 // constructors and destructor
75 //
77  :
78  DTSegmentTags_(iConfig.getParameter<edm::InputTag>("DTsegments")),
79  theHitsMin_(iConfig.getParameter<int>("HitsMin")),
80  thePruneCut_(iConfig.getParameter<double>("PruneCut")),
81  theTimeOffset_(iConfig.getParameter<double>("DTTimeOffset")),
82  theError_(iConfig.getParameter<double>("HitError")),
83  useSegmentT0_(iConfig.getParameter<bool>("UseSegmentT0")),
84  doWireCorr_(iConfig.getParameter<bool>("DoWireCorr")),
85  dropTheta_(iConfig.getParameter<bool>("DropTheta")),
86  requireBothProjections_(iConfig.getParameter<bool>("RequireBothProjections")),
87  debug(iConfig.getParameter<bool>("debug"))
88 {
89  edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
90  theService = new MuonServiceProxy(serviceParameters);
91 
92  edm::ParameterSet matchParameters = iConfig.getParameter<edm::ParameterSet>("MatchParameters");
93 
94  theMatcher = new MuonSegmentMatcher(matchParameters, theService);
95 }
96 
97 
99 {
100  if (theService) delete theService;
101  if (theMatcher) delete theMatcher;
102 }
103 
104 
105 //
106 // member functions
107 //
108 
109 // ------------ method called to produce the data ------------
110 void
112 {
113 
114 // using reco::TrackCollection;
115 
116  if (debug)
117  std::cout << " *** Muon Timimng Extractor ***" << std::endl;
118 
119  theService->update(iSetup);
120 
122 
123  // get the DT geometry
124  edm::ESHandle<DTGeometry> theDTGeom;
125  iSetup.get<MuonGeometryRecord>().get(theDTGeom);
126 
128  iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", propagator);
129  const Propagator *propag = propagator.product();
130 
131  double invbeta=0;
132  double invbetaerr=0;
133  double totalWeightInvbeta=0;
134  double totalWeightVertex=0;
135  std::vector<TimeMeasurement> tms;
136 
137  math::XYZPoint pos=muonTrack->innerPosition();
138  math::XYZVector mom=muonTrack->innerMomentum();
139  GlobalPoint posp(pos.x(), pos.y(), pos.z());
140  GlobalVector momv(mom.x(), mom.y(), mom.z());
141  FreeTrajectoryState muonFTS(posp, momv, (TrackCharge)muonTrack->charge(), theService->magneticField().product());
142 
143  // get the DT segments that were used to construct the muon
144  std::vector<const DTRecSegment4D*> range = theMatcher->matchDT(*muonTrack,iEvent);
145 
146  // create a collection on TimeMeasurements for the track
147  for (std::vector<const DTRecSegment4D*>::iterator rechit = range.begin(); rechit!=range.end();++rechit) {
148 
149  // Create the ChamberId
150  DetId id = (*rechit)->geographicalId();
151  DTChamberId chamberId(id.rawId());
152  int station = chamberId.station();
153 
154  // use only segments with both phi and theta projections present (optional)
155  bool bothProjections = ( ((*rechit)->hasPhi()) && ((*rechit)->hasZed()) );
156 
157  if (requireBothProjections_ && !bothProjections) continue;
158 
159  // loop over (theta, phi) segments
160  for (int phi=0; phi<2; phi++) {
161 
162  if (dropTheta_ && !phi) continue;
163 
164  const DTRecSegment2D* segm;
165  if (phi) segm = dynamic_cast<const DTRecSegment2D*>((*rechit)->phiSegment());
166  else segm = dynamic_cast<const DTRecSegment2D*>((*rechit)->zSegment());
167 
168  if(segm == 0) continue;
169  if (!segm->specificRecHits().size()) continue;
170 
171  const GeomDet* geomDet = theTrackingGeometry->idToDet(segm->geographicalId());
172  const std::vector<DTRecHit1D> hits1d = segm->specificRecHits();
173 
174  // store all the hits from the segment
175  for (std::vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
176 
177  const GeomDet* dtcell = theTrackingGeometry->idToDet(hiti->geographicalId());
178  TimeMeasurement thisHit;
179 
180  std::pair< TrajectoryStateOnSurface, double> tsos;
181  tsos=propag->propagateWithPath(muonFTS,dtcell->surface());
182 
183  double dist;
184  double dist_straight = dtcell->toGlobal(hiti->localPosition()).mag();
185  if (tsos.first.isValid()) {
186  dist = tsos.second+posp.mag();
187 // std::cout << "Propagate distance: " << dist << " ( innermost: " << posp.mag() << ")" << std::endl;
188  } else {
189  dist = dist_straight;
190 // std::cout << "Geom distance: " << dist << std::endl;
191  }
192 
193  thisHit.driftCell = hiti->geographicalId();
194  if (hiti->lrSide()==DTEnums::Left) thisHit.isLeft=true; else thisHit.isLeft=false;
195  thisHit.isPhi = phi;
196  thisHit.posInLayer = geomDet->toLocal(dtcell->toGlobal(hiti->localPosition())).x();
197  thisHit.distIP = dist;
198  thisHit.station = station;
199  if (useSegmentT0_ && segm->ist0Valid()) thisHit.timeCorr=segm->t0();
200  else thisHit.timeCorr=0.;
201  thisHit.timeCorr += theTimeOffset_;
202 
203  // signal propagation along the wire correction for unmached theta or phi segment hits
204  if (doWireCorr_ && !bothProjections && tsos.first.isValid()) {
205  const DTLayer* layer = theDTGeom->layer(hiti->wireId());
206  float propgL = layer->toLocal( tsos.first.globalPosition() ).y();
207  float wirePropCorr = propgL/24.4*0.00543; // signal propagation speed along the wire
208  if (thisHit.isLeft) wirePropCorr=-wirePropCorr;
209  thisHit.posInLayer += wirePropCorr;
210  const DTSuperLayer *sl = layer->superLayer();
211  float tofCorr = sl->surface().position().mag()-tsos.first.globalPosition().mag();
212  tofCorr = (tofCorr/29.979)*0.00543;
213  if (thisHit.isLeft) tofCorr=-tofCorr;
214  thisHit.posInLayer += tofCorr;
215  } else {
216  // non straight-line path correction
217  float slCorr = (dist_straight-dist)/29.979*0.00543;
218  if (thisHit.isLeft) slCorr=-slCorr;
219  thisHit.posInLayer += slCorr;
220  }
221 
222  tms.push_back(thisHit);
223  }
224  } // phi = (0,1)
225  } // rechit
226 
227  bool modified = false;
228  std::vector <double> dstnc, dsegm, dtraj, hitWeightVertex, hitWeightInvbeta, left;
229 
230  // Now loop over the measurements, calculate 1/beta and cut away outliers
231  do {
232 
233  modified = false;
234  dstnc.clear();
235  dsegm.clear();
236  dtraj.clear();
237  hitWeightVertex.clear();
238  hitWeightInvbeta.clear();
239  left.clear();
240 
241  std::vector <int> hit_idx;
242  totalWeightInvbeta=0;
243  totalWeightVertex=0;
244 
245  // Rebuild segments
246  for (int sta=1;sta<5;sta++)
247  for (int phi=0;phi<2;phi++) {
248  std::vector <TimeMeasurement> seg;
249  std::vector <int> seg_idx;
250  int tmpos=0;
251  for (std::vector<TimeMeasurement>::iterator tm=tms.begin(); tm!=tms.end(); ++tm) {
252  if ((tm->station==sta) && (tm->isPhi==phi)) {
253  seg.push_back(*tm);
254  seg_idx.push_back(tmpos);
255  }
256  tmpos++;
257  }
258 
259  unsigned int segsize = seg.size();
260  if (segsize<theHitsMin_) continue;
261 
262  double a=0, b=0;
263  std::vector <double> hitxl,hitxr,hityl,hityr;
264 
265  for (std::vector<TimeMeasurement>::iterator tm=seg.begin(); tm!=seg.end(); ++tm) {
266 
267  DetId id = tm->driftCell;
268  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
269  DTChamberId chamberId(id.rawId());
270  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
271 
272  double celly=dtcham->toLocal(dtcell->position()).z();
273 
274  if (tm->isLeft) {
275  hitxl.push_back(celly);
276  hityl.push_back(tm->posInLayer);
277  } else {
278  hitxr.push_back(celly);
279  hityr.push_back(tm->posInLayer);
280  }
281  }
282 
283  if (!fitT0(a,b,hitxl,hityl,hitxr,hityr)) {
284  if (debug)
285  std::cout << " t0 = zero, Left hits: " << hitxl.size() << " Right hits: " << hitxr.size() << std::endl;
286  continue;
287  }
288 
289  // a segment must have at least one left and one right hit
290  if ((!hitxl.size()) || (!hityl.size())) continue;
291 
292  int segidx=0;
293  for (std::vector<TimeMeasurement>::const_iterator tm=seg.begin(); tm!=seg.end(); ++tm) {
294 
295  DetId id = tm->driftCell;
296  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
297  DTChamberId chamberId(id.rawId());
298  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
299 
300  double layerZ = dtcham->toLocal(dtcell->position()).z();
301  double segmLocalPos = b+layerZ*a;
302  double hitLocalPos = tm->posInLayer;
303  int hitSide = -tm->isLeft*2+1;
304  double t0_segm = (-(hitSide*segmLocalPos)+(hitSide*hitLocalPos))/0.00543+tm->timeCorr;
305 
306  dstnc.push_back(tm->distIP);
307  dsegm.push_back(t0_segm);
308  left.push_back(hitSide);
309  hitWeightInvbeta.push_back(((double)seg.size()-2.)*tm->distIP*tm->distIP/((double)seg.size()*30.*30.*theError_*theError_));
310  hitWeightVertex.push_back(((double)seg.size()-2.)/((double)seg.size()*theError_*theError_));
311  hit_idx.push_back(seg_idx.at(segidx));
312  segidx++;
313  totalWeightInvbeta+=((double)seg.size()-2.)*tm->distIP*tm->distIP/((double)seg.size()*30.*30.*theError_*theError_);
314  totalWeightVertex+=((double)seg.size()-2.)/((double)seg.size()*theError_*theError_);
315  }
316  }
317 
318  if (totalWeightInvbeta==0) break;
319 
320  // calculate the value and error of 1/beta from the complete set of 1D hits
321  if (debug)
322  std::cout << " Points for global fit: " << dstnc.size() << std::endl;
323 
324  // inverse beta - weighted average of the contributions from individual hits
325  invbeta=0;
326  for (unsigned int i=0;i<dstnc.size();i++)
327  invbeta+=(1.+dsegm.at(i)/dstnc.at(i)*30.)*hitWeightInvbeta.at(i)/totalWeightInvbeta;
328 
329  double chimax=0.;
330  std::vector<TimeMeasurement>::iterator tmmax;
331 
332  // the dispersion of inverse beta
333  double diff;
334  for (unsigned int i=0;i<dstnc.size();i++) {
335  diff=(1.+dsegm.at(i)/dstnc.at(i)*30.)-invbeta;
336  diff=diff*diff*hitWeightInvbeta.at(i);
337  invbetaerr+=diff;
338  if (diff>chimax) {
339  tmmax=tms.begin()+hit_idx.at(i);
340  chimax=diff;
341  }
342  }
343 
344  invbetaerr=sqrt(invbetaerr/totalWeightInvbeta);
345 
346  // cut away the outliers
347  if (chimax>thePruneCut_) {
348  tms.erase(tmmax);
349  modified=true;
350  }
351 
352  if (debug)
353  std::cout << " Measured 1/beta: " << invbeta << " +/- " << invbetaerr << std::endl;
354 
355  } while (modified);
356 
357  for (unsigned int i=0;i<dstnc.size();i++) {
358  tmSequence.dstnc.push_back(dstnc.at(i));
359  tmSequence.local_t0.push_back(dsegm.at(i));
360  tmSequence.weightInvbeta.push_back(hitWeightInvbeta.at(i));
361  tmSequence.weightVertex.push_back(hitWeightVertex.at(i));
362  }
363 
364  tmSequence.totalWeightInvbeta=totalWeightInvbeta;
365  tmSequence.totalWeightVertex=totalWeightVertex;
366 
367 }
368 
369 double
370 DTTimingExtractor::fitT0(double &a, double &b, std::vector<double> xl, std::vector<double> yl, std::vector<double> xr, std::vector<double> yr ) {
371 
372  double sx=0,sy=0,sxy=0,sxx=0,ssx=0,ssy=0,s=0,ss=0;
373 
374  for (unsigned int i=0; i<xl.size(); i++) {
375  sx+=xl[i];
376  sy+=yl[i];
377  sxy+=xl[i]*yl[i];
378  sxx+=xl[i]*xl[i];
379  s++;
380  ssx+=xl[i];
381  ssy+=yl[i];
382  ss++;
383  }
384 
385  for (unsigned int i=0; i<xr.size(); i++) {
386  sx+=xr[i];
387  sy+=yr[i];
388  sxy+=xr[i]*yr[i];
389  sxx+=xr[i]*xr[i];
390  s++;
391  ssx-=xr[i];
392  ssy-=yr[i];
393  ss--;
394  }
395 
396  double delta = ss*ss*sxx+s*sx*sx+s*ssx*ssx-s*s*sxx-2*ss*sx*ssx;
397 
398  double t0_corr=0.;
399 
400  if (delta) {
401  a=(ssy*s*ssx+sxy*ss*ss+sy*sx*s-sy*ss*ssx-ssy*sx*ss-sxy*s*s)/delta;
402  b=(ssx*sy*ssx+sxx*ssy*ss+sx*sxy*s-sxx*sy*s-ssx*sxy*ss-sx*ssy*ssx)/delta;
403  t0_corr=(ssx*s*sxy+sxx*ss*sy+sx*sx*ssy-sxx*s*ssy-sx*ss*sxy-ssx*sx*sy)/delta;
404  }
405 
406  // convert drift distance to time
407  t0_corr/=-0.00543;
408 
409  return t0_corr;
410 }
411 
412 
413 //define this as a plug-in
414 //DEFINE_FWK_MODULE(DTTimingExtractor);
void update(const edm::EventSetup &setup)
update the services each event
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
DTTimingExtractor(const edm::ParameterSet &)
Constructor.
std::vector< double > local_t0
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:47
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:62
void fillTiming(TimeMeasurementSequence &tmSequence, reco::TrackRef muonTrack, const edm::Event &iEvent, const edm::EventSetup &iSetup)
edm::ESHandle< MagneticField > magneticField() const
get the magnetic field
double double double z
int TrackCharge
Definition: TrackCharge.h:4
virtual const GeomDet * idToDet(DetId) const
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:77
int iEvent
Definition: GenABIO.cc:243
T mag() const
Definition: PV3DBase.h:66
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:41
T sqrt(T t)
Definition: SSEVec.h:46
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
std::vector< double > weightInvbeta
~DTTimingExtractor()
Destructor.
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool ist0Valid() const
unsigned int theHitsMin_
Definition: DetId.h:20
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
double b
Definition: hdecay.h:120
MuonSegmentMatcher * theMatcher
double fitT0(double &a, double &b, std::vector< double > xl, std::vector< double > yl, std::vector< double > xr, std::vector< double > yr)
const DTSuperLayer * superLayer() const
Definition: DTLayer.cc:56
edm::ESHandle< GlobalTrackingGeometry > trackingGeometry() const
get the tracking geometry
const BoundPlane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
double a
Definition: hdecay.h:121
std::vector< const DTRecSegment4D * > matchDT(const reco::Track &muon, const edm::Event &event)
perform the matching
tuple cout
Definition: gather_cfg.py:121
DetId geographicalId() const
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
int station() const
Return the station number.
Definition: DTChamberId.h:53
std::vector< double > dstnc
#define debug
Definition: MEtoEDMFormat.h:34
T x() const
Definition: PV3DBase.h:61
const PositionType & position() const
std::vector< double > weightVertex
MuonServiceProxy * theService
Definition: DDAxes.h:10