CMS 3D CMS Logo

OverlapValidation.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Alignment/OfflineValidation
4 // Class: OverlapValidation
5 //
11 //
12 // Original Authors: Wolfgang Adam, Keith Ulmer
13 // Created: Thu Oct 11 14:53:32 CEST 2007
14 // $Id: OverlapValidation.cc,v 1.16 2009/11/04 19:40:27 kaulmer Exp $
15 //
16 //
17 
18 // system include files
19 #include <memory>
20 #include <vector>
21 #include <utility>
22 
23 // user include files
71 
72 // ROOT includes
73 #include "TFile.h"
74 #include "TTree.h"
75 
76 using namespace std;
77 //
78 // class decleration
79 //
80 
81 class OverlapValidation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
82 public:
83  explicit OverlapValidation(const edm::ParameterSet&);
84  ~OverlapValidation() override;
85  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
86 
87 private:
88  typedef vector<Trajectory> TrajectoryCollection;
89 
92  void analyze(const edm::Event&, const edm::EventSetup&) override;
93  void endJob() override;
94 
95  virtual void analyzeTrajectory(const Trajectory&,
96  const Propagator&,
98  const TrackerTopology* const tTopo);
99  int layerFromId(const DetId&, const TrackerTopology* const tTopo) const;
100 
101  // ----------member data ---------------------------
102 
106 
110  bool doSimHit_;
113 
115  int overlapCounts_[3];
116 
117  TTree* rootTree_;
120  const bool addExtraBranches_;
121  const int minHitsCut_;
122  const float chi2ProbCut_;
123 
126  unsigned short hitCounts_[2];
127  float chi2_[2];
128  unsigned int overlapIds_[2];
129  float predictedPositions_[3][2];
130  float predictedLocalParameters_[5][2];
131  float predictedLocalErrors_[5][2];
136  float hitPositions_[2];
137  float hitErrors_[2];
138  float hitPositionsY_[2];
139  float hitErrorsY_[2];
140  float simHitPositions_[2];
141  float simHitPositionsY_[2];
142  float clusterWidthX_[2];
143  float clusterWidthY_[2];
144  float clusterSize_[2];
146  int edge_[2];
147 
148  vector<bool> acceptLayer;
149  float momentum_;
150  uint run_, event_;
152 
153  //added by Heshy and Jared
154  float moduleX_[2];
155  float moduleY_[2];
156  float moduleZ_[2];
157  int subdetID;
161  //added by Jason
162  float localxdotglobalphi_[2];
163  float localxdotglobalr_[2];
164  float localxdotglobalz_[2];
165  float localxdotglobalx_[2];
166  float localxdotglobaly_[2];
167  float localydotglobalphi_[2];
168  float localydotglobalr_[2];
169  float localydotglobalz_[2];
170  float localydotglobalx_[2];
171  float localydotglobaly_[2];
172 };
173 
174 //
175 // constants, enums and typedefs
176 //
177 
178 //
179 // static data member definitions
180 //
181 
182 using std::vector;
183 //
184 // constructors and destructor
185 //
187  : geomToken_(esConsumes()),
188  magFieldToken_(esConsumes()),
189  topoToken_(esConsumes()),
190  config_(iConfig),
191  rootTree_(nullptr),
193  compressionSettings_(iConfig.getUntrackedParameter<int>("compressionSettings", -1)),
194  addExtraBranches_(false),
195  minHitsCut_(6),
196  chi2ProbCut_(0.001) {
197  usesResource(TFileService::kSharedResource);
199  //now do what ever initialization is needed
200  trajectoryTag_ = iConfig.getParameter<edm::InputTag>("trajectories");
202  doSimHit_ = iConfig.getParameter<bool>("associateStrip");
204 
205  overlapCounts_[0] = 0; // #trajectories
206  overlapCounts_[1] = 0; // #hits
207  overlapCounts_[2] = 0; // #overlap hits
208  acceptLayer.resize(7, false);
209  acceptLayer[PixelSubdetector::PixelBarrel] = iConfig.getParameter<bool>("usePXB");
210  acceptLayer[PixelSubdetector::PixelEndcap] = iConfig.getParameter<bool>("usePXF");
211  acceptLayer[StripSubdetector::TIB] = iConfig.getParameter<bool>("useTIB");
212  acceptLayer[StripSubdetector::TOB] = iConfig.getParameter<bool>("useTOB");
213  acceptLayer[StripSubdetector::TID] = iConfig.getParameter<bool>("useTID");
214  acceptLayer[StripSubdetector::TEC] = iConfig.getParameter<bool>("useTEC");
215  barrelOnly_ = iConfig.getParameter<bool>("barrelOnly");
216 
218  //
219  // root output
220  //
221  if (compressionSettings_ > 0) {
222  fs->file().SetCompressionSettings(compressionSettings_);
223  }
224 
225  rootTree_ = fs->make<TTree>("Overlaps", "Overlaps");
226  if (addExtraBranches_) {
227  rootTree_->Branch("hitCounts", hitCounts_, "found/s:lost/s");
228  rootTree_->Branch("chi2", chi2_, "chi2/F:ndf/F");
229  rootTree_->Branch("path", &overlapPath_, "path/F");
230  }
231  rootTree_->Branch("layer", &layer_, "layer/i");
232  rootTree_->Branch("detids", overlapIds_, "id[2]/i");
233  rootTree_->Branch("predPos", predictedPositions_, "gX[2]/F:gY[2]/F:gZ[2]/F");
234  rootTree_->Branch("predPar", predictedLocalParameters_, "predQP[2]/F:predDX[2]/F:predDY[2]/F:predX[2]/F:predY[2]/F");
235  rootTree_->Branch("predErr", predictedLocalErrors_, "predEQP[2]/F:predEDX[2]/F:predEDY[2]/F:predEX[2]/F:predEY[2]/F");
236  rootTree_->Branch("predEDeltaX", &predictedDeltaXError_, "sigDeltaX/F");
237  rootTree_->Branch("predEDeltaY", &predictedDeltaYError_, "sigDeltaY/F");
238  rootTree_->Branch("relSignX", &relativeXSign_, "relSignX/B");
239  rootTree_->Branch("relSignY", &relativeYSign_, "relSignY/B");
240  rootTree_->Branch("hitX", hitPositions_, "hitX[2]/F");
241  rootTree_->Branch("hitEX", hitErrors_, "hitEX[2]/F");
242  rootTree_->Branch("hitY", hitPositionsY_, "hitY[2]/F");
243  rootTree_->Branch("hitEY", hitErrorsY_, "hitEY[2]/F");
244  if (addExtraBranches_) {
245  rootTree_->Branch("simX", simHitPositions_, "simX[2]/F");
246  rootTree_->Branch("simY", simHitPositionsY_, "simY[2]/F");
247  rootTree_->Branch("clusterSize", clusterSize_, "clusterSize[2]/F");
248  rootTree_->Branch("clusterWidthX", clusterWidthX_, "clusterWidthX[2]/F");
249  rootTree_->Branch("clusterWidthY", clusterWidthY_, "clusterWidthY[2]/F");
250  rootTree_->Branch("clusterCharge", clusterCharge_, "clusterCharge[2]/i");
251  rootTree_->Branch("edge", edge_, "edge[2]/I");
252  }
253  rootTree_->Branch("momentum", &momentum_, "momentum/F");
254  rootTree_->Branch("run", &run_, "run/i");
255  rootTree_->Branch("event", &event_, "event/i");
256  rootTree_->Branch("subdetID", &subdetID, "subdetID/I");
257  rootTree_->Branch("moduleX", moduleX_, "moduleX[2]/F");
258  rootTree_->Branch("moduleY", moduleY_, "moduleY[2]/F");
259  rootTree_->Branch("moduleZ", moduleZ_, "moduleZ[2]/F");
260  rootTree_->Branch("localxdotglobalphi", localxdotglobalphi_, "localxdotglobalphi[2]/F");
261  rootTree_->Branch("localxdotglobalr", localxdotglobalr_, "localxdotglobalr[2]/F");
262  rootTree_->Branch("localxdotglobalz", localxdotglobalz_, "localxdotglobalz[2]/F");
263  rootTree_->Branch("localxdotglobalx", localxdotglobalx_, "localxdotglobalx[2]/F");
264  rootTree_->Branch("localxdotglobaly", localxdotglobaly_, "localxdotglobaly[2]/F");
265  rootTree_->Branch("localydotglobalphi", localydotglobalphi_, "localydotglobalphi[2]/F");
266  rootTree_->Branch("localydotglobalr", localydotglobalr_, "localydotglobalr[2]/F");
267  rootTree_->Branch("localydotglobalz", localydotglobalz_, "localydotglobalz[2]/F");
268  rootTree_->Branch("localydotglobalx", localydotglobalx_, "localydotglobalx[2]/F");
269  rootTree_->Branch("localydotglobaly", localydotglobaly_, "localydotglobaly[2]/F");
270 }
271 
272 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
275  desc.setComment("Overlap Validation analyzer plugin.");
277  //desc.add<bool>("associateStrip",false); // NB this comes from the fillPSetDescription
278  desc.addUntracked<int>("compressionSettings", -1);
279  desc.add<edm::InputTag>("trajectories", edm::InputTag("FinalTrackRefitter"));
280  desc.add<bool>("barrelOnly", false);
281  desc.add<bool>("usePXB", true);
282  desc.add<bool>("usePXF", true);
283  desc.add<bool>("useTIB", true);
284  desc.add<bool>("useTOB", true);
285  desc.add<bool>("useTID", true);
286  desc.add<bool>("useTEC", true);
287  descriptions.addWithDefaultLabel(desc);
288 }
289 
291  edm::LogWarning w("Overlaps");
292  // do anything here that needs to be done at desctruction time
293  // (e.g. close files, deallocate resources etc.)
294  w << "Counters =";
295  w << " Number of tracks: " << overlapCounts_[0];
296  w << " Number of valid hits: " << overlapCounts_[1];
297  w << " Number of overlaps: " << overlapCounts_[2];
298 }
299 
300 //
301 // member functions
302 //
303 
304 // ------------ method called to for each event ------------
306  using namespace edm;
307  //
308  // mag field & search tracker
309  //
310  magField_ = &iSetup.getData(magFieldToken_);
311  //
312  // propagator
313  //
315  //
316  // geometry
317  //
319  //
320  // make associator for SimHits
321  //
323  if (doSimHit_) {
325  associator = new TrackerHitAssociator(iEvent, hitassociatorconfig);
326  } else {
327  associator = nullptr;
328  }
329 
330  //if(doSimHit_) associator = new TrackerHitAssociator(iEvent, config_); else associator = 0;
331 
332  //
333  // trajectories (from refit)
334  //
335  //typedef vector<Trajectory> TrajectoryCollection;
336  edm::Handle<TrajectoryCollection> trajectoryCollectionHandle;
337  iEvent.getByToken(trajectoryToken_, trajectoryCollectionHandle);
338  const TrajectoryCollection* const trajectoryCollection = trajectoryCollectionHandle.product();
339 
340  //
341  // loop over trajectories from refit
342  const TrackerTopology* const tTopo = &iSetup.getData(topoToken_);
343  for (const auto& trajectory : *trajectoryCollection)
344  analyzeTrajectory(trajectory, propagator, *associator, tTopo);
345 
346  run_ = iEvent.id().run();
347  event_ = iEvent.id().event();
348 }
349 
351  const Propagator& propagator,
353  const TrackerTopology* const tTopo) {
354  typedef std::pair<const TrajectoryMeasurement*, const TrajectoryMeasurement*> Overlap;
355  typedef vector<Overlap> OverlapContainer;
356  ++overlapCounts_[0];
357 
358  OverlapContainer overlapHits;
359 
360  // quality cuts on trajectory
361  // min. # hits / matched hits
362 
363  if (trajectory.foundHits() < minHitsCut_)
364  return;
365  if (ChiSquaredProbability((double)(trajectory.chiSquared()), (double)(trajectory.ndof(false))) < chi2ProbCut_)
366  return;
367  //
368  // loop over measurements in the trajectory and calculate residuals
369  //
370 
371  vector<TrajectoryMeasurement> measurements(trajectory.measurements());
372  for (vector<TrajectoryMeasurement>::const_iterator itm = measurements.begin(); itm != measurements.end(); ++itm) {
373  //
374  // skip "invalid" (i.e. missing) hits
375  //
376  ConstRecHitPointer hit = itm->recHit();
377  DetId id = hit->geographicalId();
378  int layer(layerFromId(id, tTopo));
379  int subDet = id.subdetId();
380 
381  if (!hit->isValid()) {
382  edm::LogVerbatim("OverlapValidation") << "Invalid";
383  continue;
384  }
385  if (barrelOnly_ && (subDet == StripSubdetector::TID || subDet == StripSubdetector::TEC))
386  return;
387 
388  //edm::LogVerbatim("OverlapValidation") << "Check " <<subDet << ", layer = " << layer<<" stereo: "<< ((subDet > 2)?(SiStripDetId(id).stereo()):2);
389  //cout << "Check SubID " <<subDet << ", layer = " << layer<<" stereo: "<< ((subDet > 2)?(SiStripDetId(id).stereo()):2) << endl;
390 
391  //
392  // check for overlap: same subdet-id && layer number for
393  // two consecutive hits
394  //
395  ++overlapCounts_[1];
396  if ((layer != -1) && (acceptLayer[subDet])) {
397  for (vector<TrajectoryMeasurement>::const_iterator itmCompare = itm - 1;
398  itmCompare >= measurements.begin() && itmCompare > itm - 4;
399  --itmCompare) {
400  DetId compareId = itmCompare->recHit()->geographicalId();
401 
402  if (subDet != compareId.subdetId() || layer != layerFromId(compareId, tTopo))
403  break;
404  if (!itmCompare->recHit()->isValid())
405  continue;
406  if ((subDet == PixelSubdetector::PixelBarrel || subDet == PixelSubdetector::PixelEndcap) ||
407  (SiStripDetId(id).stereo() == SiStripDetId(compareId).stereo())) {
408  overlapHits.push_back(std::make_pair(&(*itmCompare), &(*itm)));
409  //edm::LogVerbatim("OverlapValidation") << "adding pair "<< ((subDet > 2)?(SiStripDetId(id).stereo()) : 2)
410  // << " from layer = " << layer;
411  //cout << "adding pair "<< ((subDet > 2)?(SiStripDetId(id).stereo()) : 2) << " from subDet = " << subDet << " and layer = " << layer;
412  //cout << " \t"<<run_<< "\t"<<event_<<"\t";
413  //cout << min(id.rawId(),compareId.rawId())<<"\t"<<max(id.rawId(),compareId.rawId())<<endl;
414  if (SiStripDetId(id).glued() == id.rawId())
415  edm::LogInfo("Overlaps") << "BAD GLUED: Have glued layer with id = " << id.rawId()
416  << " and glued id = " << SiStripDetId(id).glued()
417  << " and stereo = " << SiStripDetId(id).stereo() << endl;
418  if (SiStripDetId(compareId).glued() == compareId.rawId())
419  edm::LogInfo("Overlaps") << "BAD GLUED: Have glued layer with id = " << compareId.rawId()
420  << " and glued id = " << SiStripDetId(compareId).glued()
421  << " and stereo = " << SiStripDetId(compareId).stereo() << endl;
422  break;
423  }
424  }
425  }
426  }
427 
428  //
429  // Loop over all overlap pairs.
430  //
431  hitCounts_[0] = trajectory.foundHits();
432  hitCounts_[1] = trajectory.lostHits();
433  chi2_[0] = trajectory.chiSquared();
434  chi2_[1] = trajectory.ndof(false);
435 
436  for (const auto& overlapHit : overlapHits) {
437  //
438  // create reference state @ module 1 (no info from overlap hits)
439  //
440  ++overlapCounts_[2];
441  // backward predicted state at module 1
442  TrajectoryStateOnSurface bwdPred1 = overlapHit.first->backwardPredictedState();
443  if (!bwdPred1.isValid())
444  continue;
445  //cout << "momentum from backward predicted state = " << bwdPred1.globalMomentum().mag() << endl;
446  // forward predicted state at module 2
447  TrajectoryStateOnSurface fwdPred2 = overlapHit.second->forwardPredictedState();
448  //cout << "momentum from forward predicted state = " << fwdPred2.globalMomentum().mag() << endl;
449  if (!fwdPred2.isValid())
450  continue;
451  // extrapolate fwdPred2 to module 1
452  TrajectoryStateOnSurface fwdPred2At1 = propagator.propagate(fwdPred2, bwdPred1.surface());
453  if (!fwdPred2At1.isValid())
454  continue;
455  // combine fwdPred2At1 with bwdPred1 (ref. state, best estimate without hits 1 and 2)
456  TrajectoryStateOnSurface comb1 = combiner_.combine(bwdPred1, fwdPred2At1);
457  if (!comb1.isValid())
458  continue;
459  //
460  // propagation of reference parameters to module 2
461  //
462  std::pair<TrajectoryStateOnSurface, double> tsosWithS = propagator.propagateWithPath(comb1, fwdPred2.surface());
463  TrajectoryStateOnSurface comb1At2 = tsosWithS.first;
464  if (!comb1At2.isValid())
465  continue;
466  //distance of propagation from one surface to the next==could cut here
467  overlapPath_ = tsosWithS.second;
468  if (abs(overlapPath_) > 15)
469  continue; //cut to remove hit pairs > 15 cm apart
470  // global position on module 1
472  predictedPositions_[0][0] = position.x();
473  predictedPositions_[1][0] = position.y();
474  predictedPositions_[2][0] = position.z();
475  momentum_ = comb1.globalMomentum().mag();
476  //cout << "momentum from combination = " << momentum_ << endl;
477  //cout << "magnetic field from TSOS = " << comb1.magneticField()->inTesla(position).mag() << endl;
478  // local parameters and errors on module 1
479  AlgebraicVector5 pars = comb1.localParameters().vector();
480  AlgebraicSymMatrix55 errs = comb1.localError().matrix();
481  for (int i = 0; i < 5; ++i) {
482  predictedLocalParameters_[i][0] = pars[i];
483  predictedLocalErrors_[i][0] = sqrt(errs(i, i));
484  }
485  // global position on module 2
486  position = comb1At2.globalPosition();
487  predictedPositions_[0][1] = position.x();
488  predictedPositions_[1][1] = position.y();
489  predictedPositions_[2][1] = position.z();
490  // local parameters and errors on module 2
491  pars = comb1At2.localParameters().vector();
492  errs = comb1At2.localError().matrix();
493  for (int i = 0; i < 5; ++i) {
494  predictedLocalParameters_[i][1] = pars[i];
495  predictedLocalErrors_[i][1] = sqrt(errs(i, i));
496  }
497 
498  //print out local errors in X to check
499  //cout << "Predicted local error in X at 1 = " << predictedLocalErrors_[3][0] << " and predicted local error in X at 2 is = " << predictedLocalErrors_[3][1] << endl;
500  //cout << "Predicted local error in Y at 1 = " << predictedLocalErrors_[4][0] << " and predicted local error in Y at 2 is = " << predictedLocalErrors_[4][1] << endl;
501 
502  //
503  // jacobians (local-to-global@1,global 1-2,global-to-local@2)
504  //
505  JacobianLocalToCurvilinear jacLocToCurv(comb1.surface(), comb1.localParameters(), *magField_);
506  AnalyticalCurvilinearJacobian jacCurvToCurv(
507  comb1.globalParameters(), comb1At2.globalPosition(), comb1At2.globalMomentum(), tsosWithS.second);
508  JacobianCurvilinearToLocal jacCurvToLoc(comb1At2.surface(), comb1At2.localParameters(), *magField_);
509  // combined jacobian local-1-to-local-2
510  AlgebraicMatrix55 jacobian = jacLocToCurv.jacobian() * jacCurvToCurv.jacobian() * jacCurvToLoc.jacobian();
511  // covariance on module 1
512  AlgebraicSymMatrix55 covComb1 = comb1.localError().matrix();
513  // variance and correlations for predicted local_x on modules 1 and 2
514  double c00 = covComb1(3, 3);
515  double c10(0.);
516  double c11(0.);
517  for (int i = 1; i < 5; ++i) {
518  c10 += jacobian(3, i) * covComb1(i, 3);
519  for (int j = 1; j < 5; ++j)
520  c11 += jacobian(3, i) * covComb1(i, j) * jacobian(3, j);
521  }
522  // choose relative sign in order to minimize error on difference
523  double diff = c00 - 2 * fabs(c10) + c11;
524  diff = diff > 0 ? sqrt(diff) : -sqrt(-diff);
526  relativeXSign_ = c10 > 0 ? -1 : 1;
527  //
528  // now find variance and correlations for predicted local_y
529  double c00Y = covComb1(4, 4);
530  double c10Y(0.);
531  double c11Y(0.);
532  for (int i = 1; i < 5; ++i) {
533  c10Y += jacobian(4, i) * covComb1(i, 4);
534  for (int j = 1; j < 5; ++j)
535  c11Y += jacobian(4, i) * covComb1(i, j) * jacobian(4, j);
536  }
537  double diffY = c00Y - 2 * fabs(c10Y) + c11Y;
538  diffY = diffY > 0 ? sqrt(diffY) : -sqrt(-diffY);
539  predictedDeltaYError_ = diffY;
540  relativeYSign_ = c10Y > 0 ? -1 : 1;
541 
542  // information on modules and hits
543  overlapIds_[0] = overlapHit.first->recHit()->geographicalId().rawId();
544  overlapIds_[1] = overlapHit.second->recHit()->geographicalId().rawId();
545 
546  //added by Heshy and Jared
547  moduleX_[0] = overlapHit.first->recHit()->det()->surface().position().x();
548  moduleX_[1] = overlapHit.second->recHit()->det()->surface().position().x();
549  moduleY_[0] = overlapHit.first->recHit()->det()->surface().position().y();
550  moduleY_[1] = overlapHit.second->recHit()->det()->surface().position().y();
551  moduleZ_[0] = overlapHit.first->recHit()->det()->surface().position().z();
552  moduleZ_[1] = overlapHit.second->recHit()->det()->surface().position().z();
553  subdetID = overlapHit.first->recHit()->geographicalId().subdetId();
554  localxdotglobalphi_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(onezero).phi() -
555  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).phi();
556  localxdotglobalphi_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(onezero).phi() -
557  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).phi();
558  //added by Jason
559  localxdotglobalr_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(onezero).perp() -
560  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).perp();
561  localxdotglobalr_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(onezero).perp() -
562  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).perp();
563  localxdotglobalz_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(onezero).z() -
564  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).z();
565  localxdotglobalz_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(onezero).z() -
566  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).z();
567  localxdotglobalx_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(onezero).x() -
568  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).x();
569  localxdotglobalx_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(onezero).x() -
570  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).x();
571  localxdotglobaly_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(onezero).y() -
572  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).y();
573  localxdotglobaly_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(onezero).y() -
574  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).y();
575  localydotglobalr_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(zeroone).perp() -
576  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).perp();
577  localydotglobalr_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(zeroone).perp() -
578  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).perp();
579  localydotglobalz_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(zeroone).z() -
580  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).z();
581  localydotglobalz_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(zeroone).z() -
582  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).z();
583  localydotglobalx_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(zeroone).x() -
584  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).x();
585  localydotglobalx_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(zeroone).x() -
586  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).x();
587  localydotglobaly_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(zeroone).y() -
588  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).y();
589  localydotglobaly_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(zeroone).y() -
590  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).y();
591  localydotglobalphi_[0] = overlapHit.first->recHit()->det()->surface().toGlobal(zeroone).phi() -
592  overlapHit.first->recHit()->det()->surface().toGlobal(zerozero).phi();
593  localydotglobalphi_[1] = overlapHit.second->recHit()->det()->surface().toGlobal(zeroone).phi() -
594  overlapHit.second->recHit()->det()->surface().toGlobal(zerozero).phi();
595 
596  if (overlapHit.first->recHit()->geographicalId().subdetId() == StripSubdetector::TIB)
597  layer_ = layerFromId(overlapHit.first->recHit()->geographicalId().rawId(), tTopo);
598  else if (overlapHit.first->recHit()->geographicalId().subdetId() == StripSubdetector::TOB)
599  layer_ = layerFromId(overlapHit.first->recHit()->geographicalId().rawId(), tTopo) + 4;
600  else if (overlapHit.first->recHit()->geographicalId().subdetId() == StripSubdetector::TID)
601  layer_ = layerFromId(overlapHit.first->recHit()->geographicalId().rawId(), tTopo) + 10;
602  else if (overlapHit.first->recHit()->geographicalId().subdetId() == StripSubdetector::TEC)
603  layer_ = layerFromId(overlapHit.first->recHit()->geographicalId().rawId(), tTopo) + 13;
604  else if (overlapHit.first->recHit()->geographicalId().subdetId() == PixelSubdetector::PixelBarrel)
605  layer_ = layerFromId(overlapHit.first->recHit()->geographicalId().rawId(), tTopo) + 20;
606  else if (overlapHit.first->recHit()->geographicalId().subdetId() == PixelSubdetector::PixelEndcap)
607  layer_ = layerFromId(overlapHit.first->recHit()->geographicalId().rawId(), tTopo) + 30;
608  else
609  layer_ = 99;
610 
611  if (overlapIds_[0] == SiStripDetId(overlapHit.first->recHit()->geographicalId()).glued())
612  edm::LogWarning("Overlaps") << "BAD GLUED: First Id = " << overlapIds_[0] << " has glued = "
613  << SiStripDetId(overlapHit.first->recHit()->geographicalId()).glued()
614  << " and stereo = "
615  << SiStripDetId(overlapHit.first->recHit()->geographicalId()).stereo() << endl;
616  if (overlapIds_[1] == SiStripDetId(overlapHit.second->recHit()->geographicalId()).glued())
617  edm::LogWarning("Overlaps") << "BAD GLUED: Second Id = " << overlapIds_[1] << " has glued = "
618  << SiStripDetId(overlapHit.second->recHit()->geographicalId()).glued()
619  << " and stereo = "
620  << SiStripDetId(overlapHit.second->recHit()->geographicalId()).stereo() << endl;
621 
622  const TransientTrackingRecHit::ConstRecHitPointer firstRecHit = overlapHit.first->recHit();
623  const TransientTrackingRecHit::ConstRecHitPointer secondRecHit = overlapHit.second->recHit();
624  hitPositions_[0] = firstRecHit->localPosition().x();
625  hitErrors_[0] = sqrt(firstRecHit->localPositionError().xx());
626  hitPositions_[1] = secondRecHit->localPosition().x();
627  hitErrors_[1] = sqrt(secondRecHit->localPositionError().xx());
628 
629  hitPositionsY_[0] = firstRecHit->localPosition().y();
630  hitErrorsY_[0] = sqrt(firstRecHit->localPositionError().yy());
631  hitPositionsY_[1] = secondRecHit->localPosition().y();
632  hitErrorsY_[1] = sqrt(secondRecHit->localPositionError().yy());
633 
634  //cout << "printing local X hit position and error for the overlap hits. Hit 1 = " << hitPositions_[0] << "+-" << hitErrors_[0] << " and hit 2 is " << hitPositions_[1] << "+-" << hitErrors_[1] << endl;
635 
636  DetId id1 = overlapHit.first->recHit()->geographicalId();
637  DetId id2 = overlapHit.second->recHit()->geographicalId();
638  int layer1 = layerFromId(id1, tTopo);
639  int subDet1 = id1.subdetId();
640  int layer2 = layerFromId(id2, tTopo);
641  int subDet2 = id2.subdetId();
642  if (abs(hitPositions_[0]) > 5)
643  edm::LogInfo("Overlaps") << "BAD: Bad hit position: Id = " << id1.rawId()
644  << " stereo = " << SiStripDetId(id1).stereo()
645  << " glued = " << SiStripDetId(id1).glued() << " from subdet = " << subDet1
646  << " and layer = " << layer1 << endl;
647  if (abs(hitPositions_[1]) > 5)
648  edm::LogInfo("Overlaps") << "BAD: Bad hit position: Id = " << id2.rawId()
649  << " stereo = " << SiStripDetId(id2).stereo()
650  << " glued = " << SiStripDetId(id2).glued() << " from subdet = " << subDet2
651  << " and layer = " << layer2 << endl;
652 
653  // get track momentum
654  momentum_ = comb1.globalMomentum().mag();
655 
656  // get cluster size
657  if (!(subDet1 == PixelSubdetector::PixelBarrel || subDet1 == PixelSubdetector::PixelEndcap)) { //strip
658  const TransientTrackingRecHit::ConstRecHitPointer thit1 = overlapHit.first->recHit();
659  const SiStripRecHit1D* hit1 = dynamic_cast<const SiStripRecHit1D*>((*thit1).hit());
660  if (hit1) {
661  //check cluster width
662  const SiStripRecHit1D::ClusterRef& cluster1 = hit1->cluster();
663  clusterSize_[0] = (cluster1->amplitudes()).size();
664  clusterWidthX_[0] = (cluster1->amplitudes()).size();
665  clusterWidthY_[0] = -1;
666 
667  //check if cluster at edge of sensor
668  uint16_t firstStrip = cluster1->firstStrip();
669  uint16_t lastStrip = firstStrip + (cluster1->amplitudes()).size() - 1;
670  unsigned short Nstrips;
671  Nstrips = detInfo_.getNumberOfApvsAndStripLength(id1).first * 128;
672  bool atEdge = false;
673  if (firstStrip == 0 || lastStrip == (Nstrips - 1))
674  atEdge = true;
675  if (atEdge)
676  edge_[0] = 1;
677  else
678  edge_[0] = -1;
679 
680  // get cluster total charge
681  const auto& stripCharges = cluster1->amplitudes();
682  uint16_t charge = 0;
683  for (uint i = 0; i < stripCharges.size(); i++) {
684  charge += stripCharges[i];
685  }
686  clusterCharge_[0] = charge;
687  }
688 
689  const TransientTrackingRecHit::ConstRecHitPointer thit2 = overlapHit.second->recHit();
690  const SiStripRecHit1D* hit2 = dynamic_cast<const SiStripRecHit1D*>((*thit2).hit());
691  if (hit2) {
692  const SiStripRecHit1D::ClusterRef& cluster2 = hit2->cluster();
693  clusterSize_[1] = (cluster2->amplitudes()).size();
694  clusterWidthX_[1] = (cluster2->amplitudes()).size();
695  clusterWidthY_[1] = -1;
696 
697  uint16_t firstStrip = cluster2->firstStrip();
698  uint16_t lastStrip = firstStrip + (cluster2->amplitudes()).size() - 1;
699  unsigned short Nstrips;
700  Nstrips = detInfo_.getNumberOfApvsAndStripLength(id2).first * 128;
701  bool atEdge = false;
702  if (firstStrip == 0 || lastStrip == (Nstrips - 1))
703  atEdge = true;
704  if (atEdge)
705  edge_[1] = 1;
706  else
707  edge_[1] = -1;
708 
709  // get cluster total charge
710  const auto& stripCharges = cluster2->amplitudes();
711  uint16_t charge = 0;
712  for (uint i = 0; i < stripCharges.size(); i++) {
713  charge += stripCharges[i];
714  }
715  clusterCharge_[1] = charge;
716  }
717  //cout << "strip cluster size2 = " << clusterWidthX_[0] << " and size 2 = " << clusterWidthX_[1] << endl;
718  }
719 
720  if (subDet2 == PixelSubdetector::PixelBarrel || subDet2 == PixelSubdetector::PixelEndcap) { //pixel
721 
722  const TransientTrackingRecHit::ConstRecHitPointer thit1 = overlapHit.first->recHit();
723  const SiPixelRecHit* recHitPix1 = dynamic_cast<const SiPixelRecHit*>((*thit1).hit());
724  if (recHitPix1) {
725  // check for cluster size and width
726  SiPixelRecHit::ClusterRef const& cluster1 = recHitPix1->cluster();
727 
728  clusterSize_[0] = cluster1->size();
729  clusterWidthX_[0] = cluster1->sizeX();
730  clusterWidthY_[0] = cluster1->sizeY();
731 
732  // check for cluster at edge
733  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>((*trackerGeometry_).idToDet(id1));
734  const PixelTopology* topol = (&(theGeomDet->specificTopology()));
735 
736  int minPixelRow = cluster1->minPixelRow(); //x
737  int maxPixelRow = cluster1->maxPixelRow();
738  int minPixelCol = cluster1->minPixelCol(); //y
739  int maxPixelCol = cluster1->maxPixelCol();
740 
741  bool edgeHitX = (topol->isItEdgePixelInX(minPixelRow)) || (topol->isItEdgePixelInX(maxPixelRow));
742  bool edgeHitY = (topol->isItEdgePixelInY(minPixelCol)) || (topol->isItEdgePixelInY(maxPixelCol));
743  if (edgeHitX || edgeHitY)
744  edge_[0] = 1;
745  else
746  edge_[0] = -1;
747 
748  clusterCharge_[0] = (uint)cluster1->charge();
749 
750  } else {
751  edm::LogWarning("Overlaps") << "didn't find pixel cluster" << endl;
752  continue;
753  }
754 
755  const TransientTrackingRecHit::ConstRecHitPointer thit2 = overlapHit.second->recHit();
756  const SiPixelRecHit* recHitPix2 = dynamic_cast<const SiPixelRecHit*>((*thit2).hit());
757  if (recHitPix2) {
758  SiPixelRecHit::ClusterRef const& cluster2 = recHitPix2->cluster();
759 
760  clusterSize_[1] = cluster2->size();
761  clusterWidthX_[1] = cluster2->sizeX();
762  clusterWidthY_[1] = cluster2->sizeY();
763  //cout << "second pixel cluster is " << clusterSize_[1] << " pixels with x width = " << clusterWidthX_[1] << " and y width = " << clusterWidthY_[1] << endl;
764 
765  const PixelGeomDetUnit* theGeomDet = dynamic_cast<const PixelGeomDetUnit*>((*trackerGeometry_).idToDet(id2));
766  const PixelTopology* topol = (&(theGeomDet->specificTopology()));
767 
768  int minPixelRow = cluster2->minPixelRow(); //x
769  int maxPixelRow = cluster2->maxPixelRow();
770  int minPixelCol = cluster2->minPixelCol(); //y
771  int maxPixelCol = cluster2->maxPixelCol();
772 
773  bool edgeHitX = (topol->isItEdgePixelInX(minPixelRow)) || (topol->isItEdgePixelInX(maxPixelRow));
774  bool edgeHitY = (topol->isItEdgePixelInY(minPixelCol)) || (topol->isItEdgePixelInY(maxPixelCol));
775  if (edgeHitX || edgeHitY)
776  edge_[1] = 1;
777  else
778  edge_[1] = -1;
779 
780  clusterCharge_[1] = (uint)cluster2->charge();
781 
782  } else {
783  edm::LogWarning("Overlaps") << "didn't find pixel cluster" << endl;
784  continue;
785  }
786  }
787 
788  //also check for edge pixels
789 
790  //try writing out the SimHit info (for MC only)
791  if (doSimHit_) {
792  std::vector<PSimHit> psimHits1;
793  std::vector<PSimHit> psimHits2;
794  //calculate layer
795  DetId id = overlapHit.first->recHit()->geographicalId();
796  int layer(-1);
797  layer = layerFromId(id, tTopo);
798  int subDet = id.subdetId();
799  edm::LogVerbatim("OverlapValidation") << "Subdet = " << subDet << " ; layer = " << layer;
800 
801  psimHits1 = associator.associateHit(*(firstRecHit->hit()));
802  edm::LogVerbatim("OverlapValidation") << "single hit ";
803  edm::LogVerbatim("OverlapValidation") << "length of psimHits1: " << psimHits1.size();
804  if (!psimHits1.empty()) {
805  float closest_dist = 99999.9;
806  std::vector<PSimHit>::const_iterator closest_simhit = psimHits1.begin();
807  for (std::vector<PSimHit>::const_iterator m = psimHits1.begin(); m < psimHits1.end(); m++) {
808  //find closest simHit to the recHit
809  float simX = (*m).localPosition().x();
810  float dist = fabs(simX - (overlapHit.first->recHit()->localPosition().x()));
811  edm::LogVerbatim("OverlapValidation")
812  << "simHit1 simX = " << simX << " hitX = " << overlapHit.first->recHit()->localPosition().x()
813  << " distX = " << dist << " layer = " << layer;
814  if (dist < closest_dist) {
815  //cout << "found newest closest dist for simhit1" << endl;
816  closest_dist = dist;
817  closest_simhit = m;
818  }
819  }
820  //if glued layer, convert sim hit position to matchedhit surface
821  //layer index from 1-4 for TIB, 1-6 for TOB
822  // Are the sim hits on the glued layers or are they split???
823  if (subDet > 2 && !SiStripDetId(id).glued()) {
824  const GluedGeomDet* gluedDet =
825  (const GluedGeomDet*)(*trackerGeometry_).idToDet((*firstRecHit).hit()->geographicalId());
826  const StripGeomDetUnit* stripDet = (StripGeomDetUnit*)gluedDet->monoDet();
827  GlobalPoint gp = stripDet->surface().toGlobal((*closest_simhit).localPosition());
828  LocalPoint lp = gluedDet->surface().toLocal(gp);
829  LocalVector localdirection = (*closest_simhit).localDirection();
830  GlobalVector globaldirection = stripDet->surface().toGlobal(localdirection);
831  LocalVector direction = gluedDet->surface().toLocal(globaldirection);
832  float scale = -lp.z() / direction.z();
833  LocalPoint projectedPos = lp + scale * direction;
834  simHitPositions_[0] = projectedPos.x();
835  edm::LogVerbatim("OverlapValidation") << "simhit position from matched layer = " << simHitPositions_[0];
836  simHitPositionsY_[0] = projectedPos.y();
837  } else {
838  simHitPositions_[0] = (*closest_simhit).localPosition().x();
839  simHitPositionsY_[0] = (*closest_simhit).localPosition().y();
840  edm::LogVerbatim("OverlapValidation") << "simhit position from non-matched layer = " << simHitPositions_[0];
841  }
842  edm::LogVerbatim("OverlapValidation") << "hit position = " << hitPositions_[0];
843  } else {
844  simHitPositions_[0] = -99.;
845  simHitPositionsY_[0] = -99.;
846  //cout << " filling simHitX: " << -99 << endl;
847  }
848 
849  psimHits2 = associator.associateHit(*(secondRecHit->hit()));
850  if (!psimHits2.empty()) {
851  float closest_dist = 99999.9;
852  std::vector<PSimHit>::const_iterator closest_simhit = psimHits2.begin();
853  for (std::vector<PSimHit>::const_iterator m = psimHits2.begin(); m < psimHits2.end(); m++) {
854  float simX = (*m).localPosition().x();
855  float dist = fabs(simX - (overlapHit.second->recHit()->localPosition().x()));
856  if (dist < closest_dist) {
857  closest_dist = dist;
858  closest_simhit = m;
859  }
860  }
861  //if glued layer, convert sim hit position to matchedhit surface
862  // if no sim hits on matched layers then this section can be removed
863  if (subDet > 2 && !SiStripDetId(id).glued()) {
864  const GluedGeomDet* gluedDet =
865  (const GluedGeomDet*)(*trackerGeometry_).idToDet((*secondRecHit).hit()->geographicalId());
866  const StripGeomDetUnit* stripDet = (StripGeomDetUnit*)gluedDet->monoDet();
867  GlobalPoint gp = stripDet->surface().toGlobal((*closest_simhit).localPosition());
868  LocalPoint lp = gluedDet->surface().toLocal(gp);
869  LocalVector localdirection = (*closest_simhit).localDirection();
870  GlobalVector globaldirection = stripDet->surface().toGlobal(localdirection);
871  LocalVector direction = gluedDet->surface().toLocal(globaldirection);
872  float scale = -lp.z() / direction.z();
873  LocalPoint projectedPos = lp + scale * direction;
874  simHitPositions_[1] = projectedPos.x();
875  simHitPositionsY_[1] = projectedPos.y();
876  } else {
877  simHitPositions_[1] = (*closest_simhit).localPosition().x();
878  simHitPositionsY_[1] = (*closest_simhit).localPosition().y();
879  }
880  } else {
881  simHitPositions_[1] = -99.;
882  simHitPositionsY_[1] = -99.;
883  }
884  }
885  rootTree_->Fill();
886  }
887 }
888 
889 int OverlapValidation::layerFromId(const DetId& id, const TrackerTopology* const tTopo) const {
890  TrackerAlignableId aliid;
891  std::pair<int, int> subdetandlayer = aliid.typeAndLayerFromDetId(id, tTopo);
892  int layer = subdetandlayer.second;
893 
894  return layer;
895 }
896 
898  if (rootTree_) {
899  rootTree_->GetDirectory()->cd();
900  rootTree_->Write();
901  delete rootTree_;
902  }
903 }
904 
905 //define this as a plug-in
ClusterRef cluster() const
Definition: SiPixelRecHit.h:47
edm::InputTag trajectoryTag_
size
Write out results.
static const std::string kSharedResource
Definition: TFileService.h:76
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Log< level::Info, true > LogVerbatim
const Point2DBase< float, LocalTag > zerozero
static constexpr auto TEC
int layerFromId(const DetId &, const TrackerTopology *const tTopo) const
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
TrajectoryStateCombiner combiner_
edm::EDGetTokenT< TrajectoryCollection > trajectoryToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int overlapIds_[2]
~OverlapValidation() override
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static void fillPSetDescription(edm::ParameterSetDescription &descriptions)
void endJob() override
const LocalTrajectoryError & localError() const
vector< bool > acceptLayer
std::string fullPath() const
Definition: FileInPath.cc:161
const bool addExtraBranches_
T w() const
T z() const
Definition: PV3DBase.h:61
float predictedPositions_[3][2]
const Point2DBase< float, LocalTag > onezero
T const * product() const
Definition: Handle.h:70
const GlobalTrajectoryParameters & globalParameters() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepStd< double, 5, 5 > > AlgebraicMatrix55
uint32_t glued() const
Definition: SiStripDetId.h:163
float chiSquared() const
Definition: Trajectory.h:241
int lostHits() const
Definition: Trajectory.h:217
const LocalTrajectoryParameters & localParameters() const
const SurfaceType & surface() const
int foundHits() const
Definition: Trajectory.h:206
LocalPoint toLocal(const GlobalPoint &gp) const
const Point2DBase< float, LocalTag > zeroone
unsigned short hitCounts_[2]
example_stream void analyze(const edm::Event &, const edm::EventSetup &) override
void analyze(const edm::Event &, const edm::EventSetup &) override
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
DataContainer const & measurements() const
Definition: Trajectory.h:178
const MagneticField * magField_
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
virtual bool isItEdgePixelInX(int ixbin) const =0
int iEvent
Definition: GenABIO.cc:224
int ndof(bool bon=true) const
Definition: Trajectory.cc:97
virtual void analyzeTrajectory(const Trajectory &, const Propagator &, TrackerHitAssociator &, const TrackerTopology *const tTopo)
GlobalPoint globalPosition() const
SiStripDetInfo detInfo_
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
vector< Trajectory > TrajectoryCollection
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
T sqrt(T t)
Definition: SSEVec.h:23
AlgebraicVector5 vector() const
uint32_t stereo() const
Definition: SiStripDetId.h:168
T mag() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ROOT::Math::SVector< double, 5 > AlgebraicVector5
float ChiSquaredProbability(double chiSquared, double nrDOF)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
const TrackerGeometry * trackerGeometry_
SiStripDetInfo read(std::string filePath)
float predictedLocalErrors_[5][2]
static constexpr auto TOB
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
const GeomDetUnit * monoDet() const
Definition: GluedGeomDet.h:19
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
Log< level::Info, false > LogInfo
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:18
OverlapValidation(const edm::ParameterSet &)
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
const AlgebraicMatrix55 & jacobian() const
Definition: DetId.h:17
float predictedLocalParameters_[5][2]
static const char clusterCharge_[]
static constexpr auto TIB
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
edm::FileInPath FileInPath_
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
const int compressionSettings_
GlobalVector globalMomentum() const
edm::ParameterSet config_
virtual bool isItEdgePixelInY(int iybin) const =0
HLT enums.
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
static int position[264][3]
Definition: ReadPGInfo.cc:289
const AlgebraicSymMatrix55 & matrix() const
std::pair< int, int > typeAndLayerFromDetId(const DetId &detId, const TrackerTopology *tTopo) const
static constexpr char const *const kDefaultFile
Log< level::Warning, false > LogWarning
ClusterRef cluster() const
static constexpr auto TID
TSOS combine(const TSOS &pTsos1, const TSOS &pTsos2) const
Our base class.
Definition: SiPixelRecHit.h:23