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