CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Private Member Functions | Private Attributes
DTTimingExtractor Class Reference

#include <RecoMuon/MuonIdentification/src/DTTimingExtractor.cc>

Classes

class  TimeMeasurement
 

Public Member Functions

 DTTimingExtractor (const edm::ParameterSet &, edm::ConsumesCollector &iC)
 Constructor. More...
 
void fillTiming (TimeMeasurementSequence &tmSequence, reco::TrackRef muonTrack, const edm::Event &iEvent, const edm::EventSetup &iSetup)
 
 ~DTTimingExtractor ()
 Destructor. More...
 

Private Member Functions

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)
 

Private Attributes

bool debug
 
bool doWireCorr_
 
bool dropTheta_
 
edm::InputTag DTSegmentTags_
 
bool requireBothProjections_
 
double theError_
 
unsigned int theHitsMin_
 
MuonSegmentMatchertheMatcher
 
double thePruneCut_
 
MuonServiceProxytheService
 
double theTimeOffset_
 
bool useSegmentT0_
 

Detailed Description

Extracts timing information associated to a muon track

Description: Produce timing information for a muon track using 1D DT hits from segments used to build the track

Definition at line 54 of file DTTimingExtractor.h.

Constructor & Destructor Documentation

DTTimingExtractor::DTTimingExtractor ( const edm::ParameterSet iConfig,
edm::ConsumesCollector iC 
)

Constructor.

Definition at line 75 of file DTTimingExtractor.cc.

References edm::ParameterSet::getParameter(), MuonServiceProxy_cff::MuonServiceProxy, theMatcher, and theService.

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 }
T getParameter(std::string const &) const
edm::InputTag DTSegmentTags_
unsigned int theHitsMin_
MuonSegmentMatcher * theMatcher
MuonServiceProxy * theService
DTTimingExtractor::~DTTimingExtractor ( )

Destructor.

Definition at line 97 of file DTTimingExtractor.cc.

References theMatcher, and theService.

98 {
99  if (theService) delete theService;
100  if (theMatcher) delete theMatcher;
101 }
MuonSegmentMatcher * theMatcher
MuonServiceProxy * theService

Member Function Documentation

void DTTimingExtractor::fillTiming ( TimeMeasurementSequence tmSequence,
reco::TrackRef  muonTrack,
const edm::Event iEvent,
const edm::EventSetup iSetup 
)

Definition at line 110 of file DTTimingExtractor.cc.

References a, b, gather_cfg::cout, debug, diffTreeTool::diff, DTTimingExtractor::TimeMeasurement::distIP, doWireCorr_, DTTimingExtractor::TimeMeasurement::driftCell, dropTheta_, TimeMeasurementSequence::dstnc, fitT0(), TrackingRecHit::geographicalId(), edm::EventSetup::get(), i, GlobalTrackingGeometry::idToDet(), if(), DTTimingExtractor::TimeMeasurement::isLeft, DTTimingExtractor::TimeMeasurement::isPhi, DTRecSegment2D::ist0Valid(), DTEnums::Left, TimeMeasurementSequence::local_t0, PV3DBase< T, PVType, FrameType >::mag(), mag(), MuonServiceProxy::magneticField(), MuonSegmentMatcher::matchDT(), phi, DTTimingExtractor::TimeMeasurement::posInLayer, GloballyPositioned< T >::position(), GeomDet::position(), edm::ESHandle< class >::product(), Propagator::propagateWithPath(), MuonServiceProxy::propagator(), requireBothProjections_, DTRecSegment2D::specificRecHits(), mathSSE::sqrt(), DTChamberId::station(), relativeConstraints::station, DTTimingExtractor::TimeMeasurement::station, DTLayer::superLayer(), GeomDet::surface(), DTRecSegment2D::t0(), theError_, theHitsMin_, theMatcher, thePruneCut_, theService, theTimeOffset_, MuonServiceProxy::theTrackingGeometry, DTTimingExtractor::TimeMeasurement::timeCorr, GeomDet::toGlobal(), GeomDet::toLocal(), TimeMeasurementSequence::totalWeightInvbeta, TimeMeasurementSequence::totalWeightVertex, MuonServiceProxy::trackingGeometry(), MuonServiceProxy::update(), useSegmentT0_, TimeMeasurementSequence::weightInvbeta, TimeMeasurementSequence::weightVertex, PV3DBase< T, PVType, FrameType >::x(), detailsBasic3DVector::y, and detailsBasic3DVector::z.

Referenced by MuonTimingFiller::fillTiming().

111 {
112 
113 // using reco::TrackCollection;
114 
115  if (debug)
116  std::cout << " *** Muon Timimng Extractor ***" << std::endl;
117 
118  theService->update(iSetup);
119 
120  const GlobalTrackingGeometry *theTrackingGeometry = &*theService->trackingGeometry();
121 
122  // get the DT geometry
123  edm::ESHandle<DTGeometry> theDTGeom;
124  iSetup.get<MuonGeometryRecord>().get(theDTGeom);
125 
126  edm::ESHandle<Propagator> propagator;
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  if (debug)
146  std::cout << " The muon track matches " << range.size() << " segments." << std::endl;
147 
148 
149  // create a collection on TimeMeasurements for the track
150  for (std::vector<const DTRecSegment4D*>::iterator rechit = range.begin(); rechit!=range.end();++rechit) {
151 
152  // Create the ChamberId
153  DetId id = (*rechit)->geographicalId();
154  DTChamberId chamberId(id.rawId());
155  int station = chamberId.station();
156  if (debug) std::cout << "Matched DT segment in station " << station << std::endl;
157 
158  // use only segments with both phi and theta projections present (optional)
159  bool bothProjections = ( ((*rechit)->hasPhi()) && ((*rechit)->hasZed()) );
160 
161  if (requireBothProjections_ && !bothProjections) continue;
162 
163  // loop over (theta, phi) segments
164  for (int phi=0; phi<2; phi++) {
165 
166  if (dropTheta_ && !phi) continue;
167 
168  const DTRecSegment2D* segm;
169  if (phi) segm = dynamic_cast<const DTRecSegment2D*>((*rechit)->phiSegment());
170  else segm = dynamic_cast<const DTRecSegment2D*>((*rechit)->zSegment());
171 
172  if(segm == 0) continue;
173  if (!segm->specificRecHits().size()) continue;
174 
175  const GeomDet* geomDet = theTrackingGeometry->idToDet(segm->geographicalId());
176  const std::vector<DTRecHit1D> hits1d = segm->specificRecHits();
177 
178  // store all the hits from the segment
179  for (std::vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
180 
181  const GeomDet* dtcell = theTrackingGeometry->idToDet(hiti->geographicalId());
182  TimeMeasurement thisHit;
183 
184  std::pair< TrajectoryStateOnSurface, double> tsos;
185  tsos=propag->propagateWithPath(muonFTS,dtcell->surface());
186 
187  double dist;
188  double dist_straight = dtcell->toGlobal(hiti->localPosition()).mag();
189  if (tsos.first.isValid()) {
190  dist = tsos.second+posp.mag();
191 // std::cout << "Propagate distance: " << dist << " ( innermost: " << posp.mag() << ")" << std::endl;
192  } else {
193  dist = dist_straight;
194 // std::cout << "Geom distance: " << dist << std::endl;
195  }
196 
197  thisHit.driftCell = hiti->geographicalId();
198  if (hiti->lrSide()==DTEnums::Left) thisHit.isLeft=true; else thisHit.isLeft=false;
199  thisHit.isPhi = phi;
200  thisHit.posInLayer = geomDet->toLocal(dtcell->toGlobal(hiti->localPosition())).x();
201  thisHit.distIP = dist;
202  thisHit.station = station;
203  if (useSegmentT0_ && segm->ist0Valid()) thisHit.timeCorr=segm->t0();
204  else thisHit.timeCorr=0.;
205  thisHit.timeCorr += theTimeOffset_;
206 
207  // signal propagation along the wire correction for unmached theta or phi segment hits
208  if (doWireCorr_ && !bothProjections && tsos.first.isValid()) {
209  const DTLayer* layer = theDTGeom->layer(hiti->wireId());
210  float propgL = layer->toLocal( tsos.first.globalPosition() ).y();
211  float wirePropCorr = propgL/24.4*0.00543; // signal propagation speed along the wire
212  if (thisHit.isLeft) wirePropCorr=-wirePropCorr;
213  thisHit.posInLayer += wirePropCorr;
214  const DTSuperLayer *sl = layer->superLayer();
215  float tofCorr = sl->surface().position().mag()-tsos.first.globalPosition().mag();
216  tofCorr = (tofCorr/29.979)*0.00543;
217  if (thisHit.isLeft) tofCorr=-tofCorr;
218  thisHit.posInLayer += tofCorr;
219  } else {
220  // non straight-line path correction
221  float slCorr = (dist_straight-dist)/29.979*0.00543;
222  if (thisHit.isLeft) slCorr=-slCorr;
223  thisHit.posInLayer += slCorr;
224  }
225 
226  if (debug) std::cout << " dist: " << dist << " t0: " << thisHit.posInLayer << std::endl;
227 
228  tms.push_back(thisHit);
229  }
230  } // phi = (0,1)
231  } // rechit
232 
233  bool modified = false;
234  std::vector <double> dstnc, dsegm, dtraj, hitWeightVertex, hitWeightInvbeta, left;
235 
236  // Now loop over the measurements, calculate 1/beta and cut away outliers
237  do {
238 
239  modified = false;
240  dstnc.clear();
241  dsegm.clear();
242  dtraj.clear();
243  hitWeightVertex.clear();
244  hitWeightInvbeta.clear();
245  left.clear();
246 
247  std::vector <int> hit_idx;
248  totalWeightInvbeta=0;
249  totalWeightVertex=0;
250 
251  // Rebuild segments
252  for (int sta=1;sta<5;sta++)
253  for (int phi=0;phi<2;phi++) {
254  std::vector <TimeMeasurement> seg;
255  std::vector <int> seg_idx;
256  int tmpos=0;
257  for (std::vector<TimeMeasurement>::iterator tm=tms.begin(); tm!=tms.end(); ++tm) {
258  if ((tm->station==sta) && (tm->isPhi==phi)) {
259  seg.push_back(*tm);
260  seg_idx.push_back(tmpos);
261  }
262  tmpos++;
263  }
264 
265  unsigned int segsize = seg.size();
266  if (segsize<theHitsMin_) continue;
267 
268  double a=0, b=0;
269  std::vector <double> hitxl,hitxr,hityl,hityr;
270 
271  for (std::vector<TimeMeasurement>::iterator tm=seg.begin(); tm!=seg.end(); ++tm) {
272 
273  DetId id = tm->driftCell;
274  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
275  DTChamberId chamberId(id.rawId());
276  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
277 
278  double celly=dtcham->toLocal(dtcell->position()).z();
279 
280  if (tm->isLeft) {
281  hitxl.push_back(celly);
282  hityl.push_back(tm->posInLayer);
283  } else {
284  hitxr.push_back(celly);
285  hityr.push_back(tm->posInLayer);
286  }
287  }
288 
289  if (!fitT0(a,b,hitxl,hityl,hitxr,hityr)) {
290  if (debug)
291  std::cout << " t0 = zero, Left hits: " << hitxl.size() << " Right hits: " << hitxr.size() << std::endl;
292  continue;
293  }
294 
295  // a segment must have at least one left and one right hit
296  if ((!hitxl.size()) || (!hityl.size())) continue;
297 
298  int segidx=0;
299  for (std::vector<TimeMeasurement>::const_iterator tm=seg.begin(); tm!=seg.end(); ++tm) {
300 
301  DetId id = tm->driftCell;
302  const GeomDet* dtcell = theTrackingGeometry->idToDet(id);
303  DTChamberId chamberId(id.rawId());
304  const GeomDet* dtcham = theTrackingGeometry->idToDet(chamberId);
305 
306  double layerZ = dtcham->toLocal(dtcell->position()).z();
307  double segmLocalPos = b+layerZ*a;
308  double hitLocalPos = tm->posInLayer;
309  int hitSide = -tm->isLeft*2+1;
310  double t0_segm = (-(hitSide*segmLocalPos)+(hitSide*hitLocalPos))/0.00543+tm->timeCorr;
311 
312  if (debug) std::cout << " Segm hit. dstnc: " << tm->distIP << " t0: " << t0_segm << std::endl;
313 
314  dstnc.push_back(tm->distIP);
315  dsegm.push_back(t0_segm);
316  left.push_back(hitSide);
317  hitWeightInvbeta.push_back(((double)seg.size()-2.)*tm->distIP*tm->distIP/((double)seg.size()*30.*30.*theError_*theError_));
318  hitWeightVertex.push_back(((double)seg.size()-2.)/((double)seg.size()*theError_*theError_));
319  hit_idx.push_back(seg_idx.at(segidx));
320  segidx++;
321  totalWeightInvbeta+=((double)seg.size()-2.)*tm->distIP*tm->distIP/((double)seg.size()*30.*30.*theError_*theError_);
322  totalWeightVertex+=((double)seg.size()-2.)/((double)seg.size()*theError_*theError_);
323  }
324  }
325 
326  if (totalWeightInvbeta==0) break;
327 
328  // calculate the value and error of 1/beta from the complete set of 1D hits
329  if (debug)
330  std::cout << " Points for global fit: " << dstnc.size() << std::endl;
331 
332  // inverse beta - weighted average of the contributions from individual hits
333  invbeta=0;
334  for (unsigned int i=0;i<dstnc.size();i++)
335  invbeta+=(1.+dsegm.at(i)/dstnc.at(i)*30.)*hitWeightInvbeta.at(i)/totalWeightInvbeta;
336 
337  double chimax=0.;
338  std::vector<TimeMeasurement>::iterator tmmax;
339 
340  // the dispersion of inverse beta
341  double diff;
342  for (unsigned int i=0;i<dstnc.size();i++) {
343  diff=(1.+dsegm.at(i)/dstnc.at(i)*30.)-invbeta;
344  diff=diff*diff*hitWeightInvbeta.at(i);
345  invbetaerr+=diff;
346  if (diff>chimax) {
347  tmmax=tms.begin()+hit_idx.at(i);
348  chimax=diff;
349  }
350  }
351 
352  invbetaerr=sqrt(invbetaerr/totalWeightInvbeta);
353 
354  // cut away the outliers
355  if (chimax>thePruneCut_) {
356  tms.erase(tmmax);
357  modified=true;
358  }
359 
360  if (debug)
361  std::cout << " Measured 1/beta: " << invbeta << " +/- " << invbetaerr << std::endl;
362 
363  } while (modified);
364 
365  for (unsigned int i=0;i<dstnc.size();i++) {
366  tmSequence.dstnc.push_back(dstnc.at(i));
367  tmSequence.local_t0.push_back(dsegm.at(i));
368  tmSequence.weightInvbeta.push_back(hitWeightInvbeta.at(i));
369  tmSequence.weightVertex.push_back(hitWeightVertex.at(i));
370  }
371 
372  tmSequence.totalWeightInvbeta=totalWeightInvbeta;
373  tmSequence.totalWeightVertex=totalWeightVertex;
374 
375 }
void update(const edm::EventSetup &setup)
update the services each event
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
edm::ESHandle< MagneticField > magneticField() const
get the magnetic field
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
float float float z
int TrackCharge
Definition: TrackCharge.h:4
virtual const GeomDet * idToDet(DetId) const
T mag() const
Definition: PV3DBase.h:67
const Surface::PositionType & position() const
The position (origin of the R.F.)
Definition: GeomDet.h:46
T sqrt(T t)
Definition: SSEVec.h:48
std::vector< double > weightInvbeta
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
unsigned int theHitsMin_
Definition: DetId.h:18
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:86
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
tuple cout
Definition: gather_cfg.py:121
DetId geographicalId() const
double t0() const
Get the segment t0 (if recomputed, 0 is returned otherwise)
std::vector< double > dstnc
if(conf.exists("allCellsPositionCalc"))
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
double 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 
)
private

Definition at line 378 of file DTTimingExtractor.cc.

References delta, i, alignCSCRings::s, and contentValuesCheck::ss.

Referenced by fillTiming().

378  {
379 
380  double sx=0,sy=0,sxy=0,sxx=0,ssx=0,ssy=0,s=0,ss=0;
381 
382  for (unsigned int i=0; i<xl.size(); i++) {
383  sx+=xl[i];
384  sy+=yl[i];
385  sxy+=xl[i]*yl[i];
386  sxx+=xl[i]*xl[i];
387  s++;
388  ssx+=xl[i];
389  ssy+=yl[i];
390  ss++;
391  }
392 
393  for (unsigned int i=0; i<xr.size(); i++) {
394  sx+=xr[i];
395  sy+=yr[i];
396  sxy+=xr[i]*yr[i];
397  sxx+=xr[i]*xr[i];
398  s++;
399  ssx-=xr[i];
400  ssy-=yr[i];
401  ss--;
402  }
403 
404  double delta = ss*ss*sxx+s*sx*sx+s*ssx*ssx-s*s*sxx-2*ss*sx*ssx;
405 
406  double t0_corr=0.;
407 
408  if (delta) {
409  a=(ssy*s*ssx+sxy*ss*ss+sy*sx*s-sy*ss*ssx-ssy*sx*ss-sxy*s*s)/delta;
410  b=(ssx*sy*ssx+sxx*ssy*ss+sx*sxy*s-sxx*sy*s-ssx*sxy*ss-sx*ssy*ssx)/delta;
411  t0_corr=(ssx*s*sxy+sxx*ss*sy+sx*sx*ssy-sxx*s*ssy-sx*ss*sxy-ssx*sx*sy)/delta;
412  }
413 
414  // convert drift distance to time
415  t0_corr/=-0.00543;
416 
417  return t0_corr;
418 }
dbl * delta
Definition: mlp_gen.cc:36
int i
Definition: DBlmapReader.cc:9
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121

Member Data Documentation

bool DTTimingExtractor::debug
private

Definition at line 90 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::doWireCorr_
private

Definition at line 87 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::dropTheta_
private

Definition at line 88 of file DTTimingExtractor.h.

Referenced by fillTiming().

edm::InputTag DTTimingExtractor::DTSegmentTags_
private

Definition at line 81 of file DTTimingExtractor.h.

bool DTTimingExtractor::requireBothProjections_
private

Definition at line 89 of file DTTimingExtractor.h.

Referenced by fillTiming().

double DTTimingExtractor::theError_
private

Definition at line 85 of file DTTimingExtractor.h.

Referenced by fillTiming().

unsigned int DTTimingExtractor::theHitsMin_
private

Definition at line 82 of file DTTimingExtractor.h.

Referenced by fillTiming().

MuonSegmentMatcher* DTTimingExtractor::theMatcher
private

Definition at line 94 of file DTTimingExtractor.h.

Referenced by DTTimingExtractor(), fillTiming(), and ~DTTimingExtractor().

double DTTimingExtractor::thePruneCut_
private

Definition at line 83 of file DTTimingExtractor.h.

Referenced by fillTiming().

MuonServiceProxy* DTTimingExtractor::theService
private

Definition at line 92 of file DTTimingExtractor.h.

Referenced by DTTimingExtractor(), fillTiming(), and ~DTTimingExtractor().

double DTTimingExtractor::theTimeOffset_
private

Definition at line 84 of file DTTimingExtractor.h.

Referenced by fillTiming().

bool DTTimingExtractor::useSegmentT0_
private

Definition at line 86 of file DTTimingExtractor.h.

Referenced by fillTiming().