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