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