CMS 3D CMS Logo

HitEff.cc
Go to the documentation of this file.
1 // Package: CalibTracker/SiStripHitEfficiency
3 // Class: HitEff
4 // Original Author: Keith Ulmer--University of Colorado
5 // keith.ulmer@colorado.edu
6 //
8 
9 // system include files
10 #include <memory>
11 #include <string>
12 #include <iostream>
13 
14 // user includes
61 
62 // ROOT includes
63 #include "TMath.h"
64 #include "TH1F.h"
65 
66 // custom made printout
67 #define LOGPRINT edm::LogPrint("SiStripHitEfficiency:HitEff")
68 
69 //
70 // constructors and destructor
71 //
72 
73 using namespace std;
75  : scalerToken_(consumes<LumiScalersCollection>(conf.getParameter<edm::InputTag>("lumiScalers"))),
76  metaDataToken_(consumes<OnlineLuminosityRecord>(conf.getParameter<edm::InputTag>("metadata"))),
77  commonModeToken_(mayConsume<edm::DetSetVector<SiStripRawDigi> >(conf.getParameter<edm::InputTag>("commonMode"))),
78  siStripClusterInfo_(consumesCollector()),
79  combinatorialTracks_token_(
80  consumes<reco::TrackCollection>(conf.getParameter<edm::InputTag>("combinatorialTracks"))),
81  trajectories_token_(consumes<std::vector<Trajectory> >(conf.getParameter<edm::InputTag>("trajectories"))),
82  trajTrackAsso_token_(consumes<TrajTrackAssociationCollection>(conf.getParameter<edm::InputTag>("trajectories"))),
83  clusters_token_(
84  consumes<edmNew::DetSetVector<SiStripCluster> >(conf.getParameter<edm::InputTag>("siStripClusters"))),
85  digisCol_token_(consumes(conf.getParameter<edm::InputTag>("siStripDigis"))),
86  digisVec_token_(consumes(conf.getParameter<edm::InputTag>("siStripDigis"))),
87  trackerEvent_token_(consumes<MeasurementTrackerEvent>(conf.getParameter<edm::InputTag>("trackerEvent"))),
88  topoToken_(esConsumes()),
89  geomToken_(esConsumes()),
90  cpeToken_(esConsumes(edm::ESInputTag("", "StripCPEfromTrackAngle"))),
91  siStripQualityToken_(esConsumes()),
92  magFieldToken_(esConsumes()),
93  measurementTkToken_(esConsumes()),
94  chi2MeasurementEstimatorToken_(esConsumes(edm::ESInputTag("", "Chi2"))),
95  propagatorToken_(esConsumes(edm::ESInputTag("", "PropagatorWithMaterial"))),
96  conf_(conf) {
97  usesResource(TFileService::kSharedResource);
98  compSettings = conf_.getUntrackedParameter<int>("CompressionSettings", -1);
99  layers = conf_.getParameter<int>("Layer");
100  DEBUG = conf_.getParameter<bool>("Debug");
101  addLumi_ = conf_.getUntrackedParameter<bool>("addLumi", false);
102  addCommonMode_ = conf_.getUntrackedParameter<bool>("addCommonMode", false);
103  cutOnTracks_ = conf_.getUntrackedParameter<bool>("cutOnTracks", false);
104  trackMultiplicityCut_ = conf.getUntrackedParameter<unsigned int>("trackMultiplicity", 100);
105  useFirstMeas_ = conf_.getUntrackedParameter<bool>("useFirstMeas", false);
106  useLastMeas_ = conf_.getUntrackedParameter<bool>("useLastMeas", false);
108  conf_.getUntrackedParameter<bool>("useAllHitsFromTracksWithMissingHits", false);
109  doMissingHitsRecovery_ = conf_.getUntrackedParameter<bool>("doMissingHitsRecovery", false);
110 
111  hitRecoveryCounters.resize(k_END_OF_LAYERS, 0);
112  hitTotalCounters.resize(k_END_OF_LAYERS, 0);
113 }
114 
117  if (compSettings > 0) {
118  edm::LogInfo("SiStripHitEfficiency:HitEff") << "the compressions settings are:" << compSettings << std::endl;
119  fs->file().SetCompressionSettings(compSettings);
120  }
121 
122  traj = fs->make<TTree>("traj", "tree of trajectory positions");
123 #ifdef ExtendedCALIBTree
124  traj->Branch("timeDT", &timeDT, "timeDT/F");
125  traj->Branch("timeDTErr", &timeDTErr, "timeDTErr/F");
126  traj->Branch("timeDTDOF", &timeDTDOF, "timeDTDOF/I");
127  traj->Branch("timeECAL", &timeECAL, "timeECAL/F");
128  traj->Branch("dedx", &dedx, "dedx/F");
129  traj->Branch("dedxNOM", &dedxNOM, "dedxNOM/I");
130  traj->Branch("nLostHits", &nLostHits, "nLostHits/I");
131  traj->Branch("chi2", &chi2, "chi2/F");
132  traj->Branch("p", &p, "p/F");
133 #endif
134  traj->Branch("TrajGlbX", &TrajGlbX, "TrajGlbX/F");
135  traj->Branch("TrajGlbY", &TrajGlbY, "TrajGlbY/F");
136  traj->Branch("TrajGlbZ", &TrajGlbZ, "TrajGlbZ/F");
137  traj->Branch("TrajLocX", &TrajLocX, "TrajLocX/F");
138  traj->Branch("TrajLocY", &TrajLocY, "TrajLocY/F");
139  traj->Branch("TrajLocAngleX", &TrajLocAngleX, "TrajLocAngleX/F");
140  traj->Branch("TrajLocAngleY", &TrajLocAngleY, "TrajLocAngleY/F");
141  traj->Branch("TrajLocErrX", &TrajLocErrX, "TrajLocErrX/F");
142  traj->Branch("TrajLocErrY", &TrajLocErrY, "TrajLocErrY/F");
143  traj->Branch("ClusterLocX", &ClusterLocX, "ClusterLocX/F");
144  traj->Branch("ClusterLocY", &ClusterLocY, "ClusterLocY/F");
145  traj->Branch("ClusterLocErrX", &ClusterLocErrX, "ClusterLocErrX/F");
146  traj->Branch("ClusterLocErrY", &ClusterLocErrY, "ClusterLocErrY/F");
147  traj->Branch("ClusterStoN", &ClusterStoN, "ClusterStoN/F");
148  traj->Branch("ResX", &ResX, "ResX/F");
149  traj->Branch("ResXSig", &ResXSig, "ResXSig/F");
150  traj->Branch("ModIsBad", &ModIsBad, "ModIsBad/i");
151  traj->Branch("SiStripQualBad", &SiStripQualBad, "SiStripQualBad/i");
152  traj->Branch("withinAcceptance", &withinAcceptance, "withinAcceptance/O");
153  traj->Branch("nHits", &nHits, "nHits/I");
154  traj->Branch("pT", &pT, "pT/F");
155  traj->Branch("highPurity", &highPurity, "highPurity/O");
156  traj->Branch("trajHitValid", &trajHitValid, "trajHitValid/i");
157  traj->Branch("Id", &Id, "Id/i");
158  traj->Branch("run", &run, "run/i");
159  traj->Branch("event", &event, "event/i");
160  traj->Branch("layer", &whatlayer, "layer/i");
161  traj->Branch("tquality", &tquality, "tquality/I");
162  traj->Branch("bunchx", &bunchx, "bunchx/I");
163  if (addLumi_) {
164  traj->Branch("instLumi", &instLumi, "instLumi/F");
165  traj->Branch("PU", &PU, "PU/F");
166  }
167  if (addCommonMode_)
168  traj->Branch("commonMode", &commonMode, "commonMode/F");
169 
170  events = 0;
171  EventTrackCKF = 0;
172 
173  totalNbHits = 0;
174  missHitPerLayer.resize(k_END_OF_LAYERS, 0);
175 }
176 
177 void HitEff::analyze(const edm::Event& e, const edm::EventSetup& es) {
178  //Retrieve tracker topology from geometry
179  const TrackerTopology* tTopo = &es.getData(topoToken_);
181 
182  // bool DEBUG = false;
183 
184  LogDebug("SiStripHitEfficiency:HitEff") << "beginning analyze from HitEff" << endl;
185 
186  using namespace edm;
187  using namespace reco;
188  // Step A: Get Inputs
189 
190  int run_nr = e.id().run();
191  int ev_nr = e.id().event();
192  int bunch_nr = e.bunchCrossing();
193 
194  // Luminosity informations
197 
198  instLumi = 0;
199  PU = 0;
200  if (addLumi_) {
201  if (lumiScalers.isValid() && !lumiScalers->empty()) {
202  if (lumiScalers->begin() != lumiScalers->end()) {
203  instLumi = lumiScalers->begin()->instantLumi();
204  PU = lumiScalers->begin()->pileup();
205  }
206  } else if (metaData.isValid()) {
207  instLumi = metaData->instLumi();
208  PU = metaData->avgPileUp();
209  } else {
210  edm::LogWarning("SiStripHitEfficiencyWorker") << "could not find a source for the Luminosity and PU";
211  }
212  }
213 
214  // CM
216  if (addCommonMode_)
217  e.getByToken(commonModeToken_, commonModeDigis);
218 
219  //CombinatoriaTrack
220  edm::Handle<reco::TrackCollection> trackCollectionCKF;
221  //edm::InputTag TkTagCKF = conf_.getParameter<edm::InputTag>("combinatorialTracks");
222  e.getByToken(combinatorialTracks_token_, trackCollectionCKF);
223 
224  edm::Handle<std::vector<Trajectory> > TrajectoryCollectionCKF;
225  //edm::InputTag TkTrajCKF = conf_.getParameter<edm::InputTag>("trajectories");
226  e.getByToken(trajectories_token_, TrajectoryCollectionCKF);
227 
228  edm::Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
229  e.getByToken(trajTrackAsso_token_, trajTrackAssociationHandle);
230 
231  // Clusters
232  // get the SiStripClusters from the event
234  //e.getByLabel("siStripClusters", theClusters);
235  e.getByToken(clusters_token_, theClusters);
236 
237  //get tracker geometry
239  const TrackerGeometry* tkgeom = &(*tracker);
240 
241  //get Cluster Parameter Estimator
242  //std::string cpe = conf_.getParameter<std::string>("StripCPE");
244  const StripClusterParameterEstimator& stripcpe(*parameterestimator);
245 
246  // get the SiStripQuality records
248 
249  const MagneticField* magField_ = &es.getData(magFieldToken_);
250 
251  // get the list of module IDs with FED-detected errors
252  // - In Aug-2023, the data format was changed from DetIdCollection to DetIdVector.
253  // - To provide some level of backward-compatibility,
254  // the plugin checks for both types giving preference to the new format.
255  // - If only the old format is available, the collection is
256  // converted to the new format, then used downstream.
257  auto const& fedErrorIdsCol_h = e.getHandle(digisCol_token_);
258  auto const& fedErrorIdsVec_h = e.getHandle(digisVec_token_);
259  if (not fedErrorIdsCol_h.isValid() and not fedErrorIdsVec_h.isValid()) {
260  throw cms::Exception("InvalidProductSiStripDetIdsWithFEDErrors")
261  << "no valid product for SiStrip DetIds with FED errors (see parameter \"siStripDigis\"), "
262  "neither for new format (DetIdVector) nor old format (DetIdCollection)";
263  }
264  auto const& fedErrorIds = fedErrorIdsVec_h.isValid() ? *fedErrorIdsVec_h : fedErrorIdsCol_h->as_vector();
265 
266  edm::ESHandle<MeasurementTracker> measurementTrackerHandle = es.getHandle(measurementTkToken_);
267 
269  //e.getByLabel("MeasurementTrackerEvent", measurementTrackerEvent);
271 
273  const Propagator* thePropagator = &es.getData(propagatorToken_);
274 
275  events++;
276 
277  // *************** SiStripCluster Collection
278  const edmNew::DetSetVector<SiStripCluster>& input = *theClusters;
279 
280  //go through clusters to write out global position of good clusters for the layer understudy for comparison
281  // Loop through clusters just to print out locations
282  // Commented out to avoid discussion, should really be deleted.
283  /*
284  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
285  // DSViter is a vector of SiStripClusters located on a single module
286  unsigned int ClusterId = DSViter->id();
287  DetId ClusterDetId(ClusterId);
288  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
289 
290  edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
291  edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end();
292  for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
293  //iter is a single SiStripCluster
294  StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
295 
296  const Surface* surface;
297  surface = &(tracker->idToDet(ClusterDetId)->surface());
298  LocalPoint lp = parameters.first;
299  GlobalPoint gp = surface->toGlobal(lp);
300  unsigned int layer = ::checkLayer(ClusterId, tTopo);
301  if(DEBUG) LOGPRINT << "Found hit in cluster collection layer = " << layer << " with id = " << ClusterId << " local X position = " << lp.x() << " +- " << sqrt(parameters.second.xx()) << " matched/stereo/rphi = " << ((ClusterId & 0x3)==0) << "/" << ((ClusterId & 0x3)==1) << "/" << ((ClusterId & 0x3)==2) << endl;
302  }
303  }
304  */
305 
306  // Tracking
307  const reco::TrackCollection* tracksCKF = trackCollectionCKF.product();
308  LogDebug("SiStripHitEfficiency:HitEff") << "number ckf tracks found = " << tracksCKF->size() << endl;
309  //if (tracksCKF->size() == 1 ){
310  if (!tracksCKF->empty()) {
311  if (cutOnTracks_ && (tracksCKF->size() >= trackMultiplicityCut_))
312  return;
313  if (cutOnTracks_)
314  LogDebug("SiStripHitEfficiency:HitEff")
315  << "starting checking good event with < " << trackMultiplicityCut_ << " tracks" << endl;
316 
317  EventTrackCKF++;
318 
319 #ifdef ExtendedCALIBTree
320  //get dEdx info if available
321  edm::Handle<ValueMap<DeDxData> > dEdxUncalibHandle;
322  if (e.getByLabel("dedxMedianCTF", dEdxUncalibHandle)) {
323  const ValueMap<DeDxData> dEdxTrackUncalib = *dEdxUncalibHandle.product();
324 
325  reco::TrackRef itTrack = reco::TrackRef(trackCollectionCKF, 0);
326  dedx = dEdxTrackUncalib[itTrack].dEdx();
327  dedxNOM = dEdxTrackUncalib[itTrack].numberOfMeasurements();
328  } else {
329  dedx = -999.0;
330  dedxNOM = -999;
331  }
332 
333  //get muon and ecal timing info if available
335  if (e.getByLabel("muonsWitht0Correction", muH)) {
336  const MuonCollection& muonsT0 = *muH.product();
337  if (!muonsT0.empty()) {
338  MuonTime mt0 = muonsT0[0].time();
339  timeDT = mt0.timeAtIpInOut;
340  timeDTErr = mt0.timeAtIpInOutErr;
341  timeDTDOF = mt0.nDof;
342 
343  bool hasCaloEnergyInfo = muonsT0[0].isEnergyValid();
344  if (hasCaloEnergyInfo)
345  timeECAL = muonsT0[0].calEnergy().ecal_time;
346  }
347  } else {
348  timeDT = -999.0;
349  timeDTErr = -999.0;
350  timeDTDOF = -999;
351  timeECAL = -999.0;
352  }
353 
354 #endif
355  // actually should do a loop over all the tracks in the event here
356 
357  // Looping over traj-track associations to be able to get traj & track informations
358  for (TrajTrackAssociationCollection::const_iterator it = trajTrackAssociationHandle->begin();
359  it != trajTrackAssociationHandle->end();
360  it++) {
361  edm::Ref<std::vector<Trajectory> > itraj = it->key;
362  reco::TrackRef itrack = it->val;
363 
364  // for each track, fill some variables such as number of hits and momentum
365  nHits = itraj->foundHits();
366 #ifdef ExtendedCALIBTree
367  nLostHits = itraj->lostHits();
368  chi2 = (itraj->chiSquared() / itraj->ndof());
369  p = itraj->lastMeasurement().updatedState().globalMomentum().mag();
370 #endif
371  pT = sqrt((itraj->lastMeasurement().updatedState().globalMomentum().x() *
372  itraj->lastMeasurement().updatedState().globalMomentum().x()) +
373  (itraj->lastMeasurement().updatedState().globalMomentum().y() *
374  itraj->lastMeasurement().updatedState().globalMomentum().y()));
375 
376  // track quality
378 
379  std::vector<TrajectoryMeasurement> TMeas = itraj->measurements();
380  totalNbHits += int(TMeas.size());
381  vector<TrajectoryMeasurement>::iterator itm;
382  double xloc = 0.;
383  double yloc = 0.;
384  double xErr = 0.;
385  double yErr = 0.;
386  double angleX = -999.;
387  double angleY = -999.;
388  double xglob, yglob, zglob;
389 
390  // Check whether the trajectory has some missing hits
391  bool hasMissingHits = false;
392  int previous_layer = 999;
393  vector<unsigned int> missedLayers;
394  for (const auto& itm : TMeas) {
395  auto theHit = itm.recHit();
396  unsigned int iidd = theHit->geographicalId().rawId();
397  int layer = ::checkLayer(iidd, tTopo);
398  int missedLayer = (layer + 1);
399  int previousMissedLayer = (layer + 2);
400  int diffPreviousLayer = (layer - previous_layer);
402  //Layers from TIB + TOB
403  if (diffPreviousLayer == -2 && missedLayer > k_LayersStart && missedLayer < k_LayersAtTOBEnd) {
404  missHitPerLayer[missedLayer] += 1;
405  hasMissingHits = true;
406  }
407  //Layers from TID
408  else if (diffPreviousLayer == -2 && (missedLayer > k_LayersAtTOBEnd + 1 && missedLayer <= k_LayersAtTIDEnd)) {
409  missHitPerLayer[missedLayer] += 1;
410  hasMissingHits = true;
411  }
412  //Layers from TEC
413  else if (diffPreviousLayer == -2 && missedLayer > k_LayersAtTIDEnd && missedLayer <= k_LayersAtTECEnd) {
414  missHitPerLayer[missedLayer] += 1;
415  hasMissingHits = true;
416  }
417 
418  //##### TID Layer 11 in transition TID -> TIB : layer is in TIB, previous layer = 12
419  if ((layer > k_LayersStart && layer <= k_LayersAtTIBEnd) && (previous_layer == 12)) {
420  missHitPerLayer[11] += 1;
421  hasMissingHits = true;
422  }
423 
424  //##### TEC Layer 14 in transition TEC -> TOB : layer is in TOB, previous layer = 15
425  if ((layer > k_LayersAtTIBEnd && layer <= k_LayersAtTOBEnd) && (previous_layer == 15)) {
426  missHitPerLayer[14] += 1;
427  hasMissingHits = true;
428  }
429 
430  //####### Consecutive missing hits case #######
431 
432  //##### Layers from TIB + TOB
433  if (diffPreviousLayer == -3 && missedLayer > k_LayersStart && missedLayer < k_LayersAtTOBEnd &&
434  previousMissedLayer > k_LayersStart && previousMissedLayer < k_LayersAtTOBEnd) {
435  missHitPerLayer[missedLayer] += 1;
436  missHitPerLayer[previousMissedLayer] += 1;
437  hasMissingHits = true;
438  }
439 
440  //##### Layers from TEC
441  else if (diffPreviousLayer == -3 && missedLayer > k_LayersAtTIDEnd && missedLayer <= k_LayersAtTECEnd &&
442  previousMissedLayer > k_LayersAtTIDEnd && previousMissedLayer <= k_LayersAtTECEnd) {
443  missHitPerLayer[missedLayer] += 1;
444  missHitPerLayer[previousMissedLayer] += 1;
445  hasMissingHits = true;
446  }
447  }
448 
449  if (theHit->getType() == TrackingRecHit::Type::missing)
450  hasMissingHits = true;
451 
452  if (hasMissingHits)
453  missedLayers.push_back(layer);
454  previous_layer = layer;
455  }
456 
457  // Loop on each measurement and take it into consideration
458  //--------------------------------------------------------
459  unsigned int prev_TKlayers = 0;
460  for (itm = TMeas.begin(); itm != TMeas.end(); itm++) {
461  auto theInHit = (*itm).recHit();
462 
463  LogDebug("SiStripHitEfficiency:HitEff") << "theInHit is valid = " << theInHit->isValid() << endl;
464 
465  unsigned int iidd = theInHit->geographicalId().rawId();
466  bool foundConsMissingHits = false;
467  unsigned int TKlayers = ::checkLayer(iidd, tTopo);
468  LogDebug("SiStripHitEfficiency:HitEff") << "TKlayer from trajectory: " << TKlayers << " from module = " << iidd
469  << " matched/stereo/rphi = " << ((iidd & 0x3) == 0) << "/"
470  << ((iidd & 0x3) == 1) << "/" << ((iidd & 0x3) == 2) << endl;
471 
472  // Test first and last points of the trajectory
473  // the list of measurements starts from outer layers !!! This could change -> should add a check
474  bool isFirstMeas = (itm == (TMeas.end() - 1));
475  bool isLastMeas = (itm == (TMeas.begin()));
476 
477  if (!useFirstMeas_ && isFirstMeas)
478  continue;
479  if (!useLastMeas_ && isLastMeas)
480  continue;
481 
482  // In case of missing hit in the track, check whether to use the other hits or not.
483  if (hasMissingHits && theInHit->getType() != TrackingRecHit::Type::missing &&
485  continue;
486 
487  // If Trajectory measurement from TOB 6 or TEC 9, skip it because it's always valid they are filled later
488  if (TKlayers == 10 || TKlayers == 22) {
489  LogDebug("SiStripHitEfficiency:HitEff") << "skipping original TM for TOB 6 or TEC 9" << endl;
490  continue;
491  }
492 
493  // Make vector of TrajectoryAtInvalidHits to hold the trajectories
494  std::vector<TrajectoryAtInvalidHit> TMs;
495 
496  // Make AnalyticalPropagator to use in TAVH constructor
498 
499  // for double sided layers check both sensors--if no hit was found on either sensor surface,
500  // the trajectory measurements only have one invalid hit entry on the matched surface
501  // so get the TrajectoryAtInvalidHit for both surfaces and include them in the study
502  if (::isDoubleSided(iidd, tTopo) && ((iidd & 0x3) == 0)) {
503  // do hit eff check twice--once for each sensor
504  //add a TM for each surface
505  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 1));
506  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 2));
507  } else if (::isDoubleSided(iidd, tTopo) && (!::check2DPartner(iidd, TMeas))) {
508  // if only one hit was found the trajectory measurement is on that sensor surface, and the other surface from
509  // the matched layer should be added to the study as well
510  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 1));
511  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 2));
512  LogDebug("SiStripHitEfficiency:HitEff") << " found a hit with a missing partner";
513  } else {
514  //only add one TM for the single surface and the other will be added in the next iteration
515  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator));
516  }
517  bool missingHitAdded = false;
518 
519  vector<TrajectoryMeasurement> tmpTmeas, prev_tmpTmeas;
520  unsigned int misLayer = TKlayers + 1;
521  unsigned int previousMisLayer = TKlayers + 2;
522  //Use bool doMissingHitsRecovery to add possible missing hits based on actual/previous hit
524  if (int(TKlayers) - int(prev_TKlayers) == -2) {
525  const DetLayer* detlayer = itm->layer();
526  const LayerMeasurements layerMeasurements{*measurementTrackerHandle, *measurementTrackerEvent};
527  const TrajectoryStateOnSurface tsos = itm->updatedState();
528  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, *thePropagator, *estimator);
529 
530  if (misLayer > k_LayersAtTIDEnd && misLayer < k_LayersAtTECEnd) { //TEC
531  std::vector<ForwardDetLayer const*> negTECLayers =
532  measurementTrackerHandle->geometricSearchTracker()->negTecLayers();
533  std::vector<ForwardDetLayer const*> posTECLayers =
534  measurementTrackerHandle->geometricSearchTracker()->posTecLayers();
535  const DetLayer* tecLayerneg = negTECLayers[misLayer - k_LayersAtTIDEnd - 1];
536  const DetLayer* tecLayerpos = posTECLayers[misLayer - k_LayersAtTIDEnd - 1];
537  if (tTopo->tecSide(iidd) == 1) {
538  tmpTmeas = layerMeasurements.measurements(*tecLayerneg, tsos, *thePropagator, *estimator);
539  } else if (tTopo->tecSide(iidd) == 2) {
540  tmpTmeas = layerMeasurements.measurements(*tecLayerpos, tsos, *thePropagator, *estimator);
541  }
542  }
543 
544  else if (misLayer == (k_LayersAtTIDEnd - 1) ||
545  misLayer == k_LayersAtTIDEnd) { // This is for TID layers 12 and 13
546 
547  std::vector<ForwardDetLayer const*> negTIDLayers =
548  measurementTrackerHandle->geometricSearchTracker()->negTidLayers();
549  std::vector<ForwardDetLayer const*> posTIDLayers =
550  measurementTrackerHandle->geometricSearchTracker()->posTidLayers();
551  const DetLayer* tidLayerneg = negTIDLayers[misLayer - k_LayersAtTOBEnd - 1];
552  const DetLayer* tidLayerpos = posTIDLayers[misLayer - k_LayersAtTOBEnd - 1];
553 
554  if (tTopo->tidSide(iidd) == 1) {
555  tmpTmeas = layerMeasurements.measurements(*tidLayerneg, tsos, *thePropagator, *estimator);
556  } else if (tTopo->tidSide(iidd) == 2) {
557  tmpTmeas = layerMeasurements.measurements(*tidLayerpos, tsos, *thePropagator, *estimator);
558  }
559  }
560 
561  if (misLayer > k_LayersStart && misLayer < k_LayersAtTOBEnd) { // Barrel
562 
563  std::vector<BarrelDetLayer const*> barrelTIBLayers =
564  measurementTrackerHandle->geometricSearchTracker()->tibLayers();
565  std::vector<BarrelDetLayer const*> barrelTOBLayers =
566  measurementTrackerHandle->geometricSearchTracker()->tobLayers();
567 
568  if (misLayer > k_LayersStart && misLayer <= k_LayersAtTIBEnd) {
569  const DetLayer* tibLayer = barrelTIBLayers[misLayer - k_LayersStart - 1];
570  tmpTmeas = layerMeasurements.measurements(*tibLayer, tsos, *thePropagator, *estimator);
571  } else if (misLayer > k_LayersAtTIBEnd && misLayer < k_LayersAtTOBEnd) {
572  const DetLayer* tobLayer = barrelTOBLayers[misLayer - k_LayersAtTIBEnd - 1];
573  tmpTmeas = layerMeasurements.measurements(*tobLayer, tsos, *thePropagator, *estimator);
574  }
575  }
576  }
577  if ((int(TKlayers) > k_LayersStart && int(TKlayers) <= k_LayersAtTIBEnd) && int(prev_TKlayers) == 12) {
578  const DetLayer* detlayer = itm->layer();
579  const LayerMeasurements layerMeasurements{*measurementTrackerHandle, *measurementTrackerEvent};
580  const TrajectoryStateOnSurface tsos = itm->updatedState();
581  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, *thePropagator, *estimator);
582  std::vector<ForwardDetLayer const*> negTIDLayers =
583  measurementTrackerHandle->geometricSearchTracker()->negTidLayers();
584  std::vector<ForwardDetLayer const*> posTIDLayers =
585  measurementTrackerHandle->geometricSearchTracker()->posTidLayers();
586 
587  const DetLayer* tidLayerneg = negTIDLayers[k_LayersStart];
588  const DetLayer* tidLayerpos = posTIDLayers[k_LayersStart];
589  if (tTopo->tidSide(iidd) == 1) {
590  tmpTmeas = layerMeasurements.measurements(*tidLayerneg, tsos, *thePropagator, *estimator);
591  } else if (tTopo->tidSide(iidd) == 2) {
592  tmpTmeas = layerMeasurements.measurements(*tidLayerpos, tsos, *thePropagator, *estimator);
593  }
594  }
595 
596  if ((int(TKlayers) > k_LayersAtTIBEnd && int(TKlayers) <= k_LayersAtTOBEnd) && int(prev_TKlayers) == 15) {
597  const DetLayer* detlayer = itm->layer();
598  const LayerMeasurements layerMeasurements{*measurementTrackerHandle, *measurementTrackerEvent};
599  const TrajectoryStateOnSurface tsos = itm->updatedState();
600  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, *thePropagator, *estimator);
601 
602  std::vector<ForwardDetLayer const*> negTECLayers =
603  measurementTrackerHandle->geometricSearchTracker()->negTecLayers();
604  std::vector<ForwardDetLayer const*> posTECLayers =
605  measurementTrackerHandle->geometricSearchTracker()->posTecLayers();
606 
607  const DetLayer* tecLayerneg = negTECLayers[k_LayersStart];
608  const DetLayer* tecLayerpos = posTECLayers[k_LayersStart];
609  if (tTopo->tecSide(iidd) == 1) {
610  tmpTmeas = layerMeasurements.measurements(*tecLayerneg, tsos, *thePropagator, *estimator);
611  } else if (tTopo->tecSide(iidd) == 2) {
612  tmpTmeas = layerMeasurements.measurements(*tecLayerpos, tsos, *thePropagator, *estimator);
613  }
614  }
615 
616  //Test for two consecutive missing hits
617  if (int(TKlayers) - int(prev_TKlayers) == -3) {
618  foundConsMissingHits = true;
619  const DetLayer* detlayer = itm->layer();
620  const LayerMeasurements layerMeasurements{*measurementTrackerHandle, *measurementTrackerEvent};
621  const TrajectoryStateOnSurface tsos = itm->updatedState();
622  std::vector<DetLayer::DetWithState> compatDets = detlayer->compatibleDets(tsos, *thePropagator, *estimator);
623 
624  if (misLayer > k_LayersStart && misLayer <= k_LayersAtTOBEnd && previousMisLayer > k_LayersStart &&
625  previousMisLayer <= k_LayersAtTOBEnd) { //Barrel case
626  std::vector<BarrelDetLayer const*> barrelTIBLayers =
627  measurementTrackerHandle->geometricSearchTracker()->tibLayers();
628  std::vector<BarrelDetLayer const*> barrelTOBLayers =
629  measurementTrackerHandle->geometricSearchTracker()->tobLayers();
630  if (misLayer > k_LayersStart && misLayer < k_LayersAtTIBEnd) {
631  const DetLayer* tibLayer = barrelTIBLayers[misLayer - k_LayersStart - 1];
632  const DetLayer* prevTibLayer = barrelTIBLayers[previousMisLayer - k_LayersStart - 1];
633 
634  tmpTmeas = layerMeasurements.measurements(*tibLayer, tsos, *thePropagator, *estimator);
635  prev_tmpTmeas = layerMeasurements.measurements(*prevTibLayer, tsos, *thePropagator, *estimator);
636  } else if (misLayer > k_LayersAtTIBEnd && misLayer < k_LayersAtTOBEnd) {
637  const DetLayer* tobLayer = barrelTOBLayers[misLayer - k_LayersAtTIBEnd - 1];
638  const DetLayer* prevTobLayer = barrelTOBLayers[previousMisLayer - k_LayersAtTIBEnd - 1];
639  tmpTmeas = layerMeasurements.measurements(*tobLayer, tsos, *thePropagator, *estimator);
640  prev_tmpTmeas = layerMeasurements.measurements(*prevTobLayer, tsos, *thePropagator, *estimator);
641  }
642  } else if (misLayer > k_LayersAtTIDEnd && misLayer < k_LayersAtTECEnd &&
643  previousMisLayer > k_LayersAtTIDEnd && previousMisLayer < k_LayersAtTECEnd) { //TEC
644  std::vector<ForwardDetLayer const*> negTECLayers =
645  measurementTrackerHandle->geometricSearchTracker()->negTecLayers();
646  std::vector<ForwardDetLayer const*> posTECLayers =
647  measurementTrackerHandle->geometricSearchTracker()->posTecLayers();
648 
649  const DetLayer* tecLayerneg = negTECLayers[misLayer - k_LayersAtTIDEnd - 1];
650  const DetLayer* prevTecLayerneg = negTECLayers[previousMisLayer - k_LayersAtTIDEnd - 1];
651 
652  const DetLayer* tecLayerpos = posTECLayers[misLayer - k_LayersAtTIDEnd - 1];
653  const DetLayer* prevTecLayerpos = posTECLayers[previousMisLayer - k_LayersAtTIDEnd - 1];
654 
655  if (tTopo->tecSide(iidd) == 1) {
656  tmpTmeas = layerMeasurements.measurements(*tecLayerneg, tsos, *thePropagator, *estimator);
657  prev_tmpTmeas = layerMeasurements.measurements(*prevTecLayerneg, tsos, *thePropagator, *estimator);
658  } else if (tTopo->tecSide(iidd) == 2) {
659  tmpTmeas = layerMeasurements.measurements(*tecLayerpos, tsos, *thePropagator, *estimator);
660  prev_tmpTmeas = layerMeasurements.measurements(*prevTecLayerpos, tsos, *thePropagator, *estimator);
661  }
662  }
663  }
664  if (!tmpTmeas.empty() && !foundConsMissingHits) {
665  TrajectoryMeasurement TM_tmp(tmpTmeas.back());
666  unsigned int iidd_tmp = TM_tmp.recHit()->geographicalId().rawId();
667  if (iidd_tmp != 0) {
668  LogDebug("SiStripHitEfficiency:HitEff") << " hit actually being added to TM vector";
669  if ((!useAllHitsFromTracksWithMissingHits_ || (!useFirstMeas_ && isFirstMeas)))
670  TMs.clear();
671  if (::isDoubleSided(iidd_tmp, tTopo)) {
672  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp, tTopo, tkgeom, propagator, 1));
673  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp, tTopo, tkgeom, propagator, 2));
674  } else
675  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp, tTopo, tkgeom, propagator));
676  missingHitAdded = true;
677  hitRecoveryCounters[misLayer] += 1;
678  }
679  }
680 
681  if (!tmpTmeas.empty() && !prev_tmpTmeas.empty() &&
682  foundConsMissingHits) { //Found two consecutive missing hits
683  TrajectoryMeasurement TM_tmp1(tmpTmeas.back());
684  TrajectoryMeasurement TM_tmp2(prev_tmpTmeas.back());
685  //Inner and outer hits module IDs
686  unsigned int modIdInner = TM_tmp1.recHit()->geographicalId().rawId();
687  unsigned int modIdOuter = TM_tmp2.recHit()->geographicalId().rawId();
688  bool innerModInactive = false, outerModInactive = false;
689  for (const auto& tm : tmpTmeas) { //Check if inner module is inactive
690  unsigned int tmModId = tm.recHit()->geographicalId().rawId();
691  if (tmModId == modIdInner && tm.recHit()->getType() == 2) {
692  innerModInactive = true;
693  break;
694  }
695  }
696  for (const auto& tm : prev_tmpTmeas) { //Check if outer module is inactive
697  unsigned int tmModId = tm.recHit()->geographicalId().rawId();
698  if (tmModId == modIdOuter && tm.recHit()->getType() == 2) {
699  outerModInactive = true;
700  break; //Found the inactive module
701  }
702  }
703 
704  if (outerModInactive) { //If outer missing hit is in inactive module, recover the inner one
705  if (modIdInner != 0) {
706  LogDebug("SiStripHitEfficiency:HitEff") << " hit actually being added to TM vector";
707  if ((!useAllHitsFromTracksWithMissingHits_ || (!useFirstMeas_ && isFirstMeas)))
708  TMs.clear();
709  if (::isDoubleSided(modIdInner, tTopo)) {
710  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp1, tTopo, tkgeom, propagator, 1));
711  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp1, tTopo, tkgeom, propagator, 2));
712  } else
713  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp1, tTopo, tkgeom, propagator));
714  missingHitAdded = true;
715  hitRecoveryCounters[misLayer] += 1;
716  }
717  }
718  if (innerModInactive) { //If inner missing hit is in inactive module, recover the outer one
719  if (modIdOuter != 0) {
720  LogDebug("SiStripHitEfficiency:HitEff") << " hit actually being added to TM vector";
721  if ((!useAllHitsFromTracksWithMissingHits_ || (!useFirstMeas_ && isFirstMeas)))
722  TMs.clear();
723  if (::isDoubleSided(modIdOuter, tTopo)) {
724  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp2, tTopo, tkgeom, propagator, 1));
725  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp2, tTopo, tkgeom, propagator, 2));
726  } else
727  TMs.push_back(TrajectoryAtInvalidHit(TM_tmp2, tTopo, tkgeom, propagator));
728  missingHitAdded = true;
729  hitRecoveryCounters[previousMisLayer] += 1;
730  }
731  }
732  }
733  }
734 
735  prev_TKlayers = TKlayers;
736  if (!useFirstMeas_ && isFirstMeas && !missingHitAdded)
737  continue;
738  if (!useLastMeas_ && isLastMeas)
739  continue;
740  bool hitsWithBias = false;
741  for (auto ilayer : missedLayers) {
742  if (ilayer < TKlayers)
743  hitsWithBias = true;
744  }
745  if (hasMissingHits && theInHit->getType() != TrackingRecHit::Type::missing && !missingHitAdded &&
746  hitsWithBias && !useAllHitsFromTracksWithMissingHits_) {
747  continue;
748  }
750  //Now check for tracks at TOB6 and TEC9
751 
752  // to make sure we only propagate on the last TOB5 hit check the next entry isn't also in TOB5
753  // to avoid bias, make sure the TOB5 hit is valid (an invalid hit on TOB5 could only exist with a valid hit on TOB6)
754 
755  bool isValid = theInHit->isValid();
756  bool isLast = (itm == (TMeas.end() - 1));
757  bool isLastTOB5 = true;
758  if (!isLast) {
759  if (::checkLayer((++itm)->recHit()->geographicalId().rawId(), tTopo) == 9)
760  isLastTOB5 = false;
761  else
762  isLastTOB5 = true;
763  --itm;
764  }
765 
766  if (TKlayers == 9 && isValid && isLastTOB5) {
767  // if ( TKlayers==9 && itm==TMeas.rbegin()) {
768  // if ( TKlayers==9 && (itm==TMeas.back()) ) { // to check for only the last entry in the trajectory for propagation
769  std::vector<BarrelDetLayer const*> barrelTOBLayers =
770  measurementTrackerHandle->geometricSearchTracker()->tobLayers();
771  const DetLayer* tob6 = barrelTOBLayers[barrelTOBLayers.size() - 1];
772  const LayerMeasurements layerMeasurements{*measurementTrackerHandle, *measurementTrackerEvent};
773  const TrajectoryStateOnSurface tsosTOB5 = itm->updatedState();
774  auto tmp = layerMeasurements.measurements(*tob6, tsosTOB5, *thePropagator, *estimator);
775 
776  if (!tmp.empty()) {
777  LogDebug("SiStripHitEfficiency:HitEff") << "size of TM from propagation = " << tmp.size() << endl;
778 
779  // take the last of the TMs, which is always an invalid hit
780  // if no detId is available, ie detId==0, then no compatible layer was crossed
781  // otherwise, use that TM for the efficiency measurement
782  TrajectoryMeasurement tob6TM(tmp.back());
783  const auto& tob6Hit = tob6TM.recHit();
784 
785  if (tob6Hit->geographicalId().rawId() != 0) {
786  LogDebug("SiStripHitEfficiency:HitEff") << "tob6 hit actually being added to TM vector" << endl;
787  TMs.push_back(TrajectoryAtInvalidHit(tob6TM, tTopo, tkgeom, propagator));
788  }
789  }
790  }
791 
792  bool isLastTEC8 = true;
793  if (!isLast) {
794  if (::checkLayer((++itm)->recHit()->geographicalId().rawId(), tTopo) == 21)
795  isLastTEC8 = false;
796  else
797  isLastTEC8 = true;
798  --itm;
799  }
800 
801  if (TKlayers == 21 && isValid && isLastTEC8) {
802  std::vector<const ForwardDetLayer*> posTecLayers =
803  measurementTrackerHandle->geometricSearchTracker()->posTecLayers();
804  const DetLayer* tec9pos = posTecLayers[posTecLayers.size() - 1];
805  std::vector<const ForwardDetLayer*> negTecLayers =
806  measurementTrackerHandle->geometricSearchTracker()->negTecLayers();
807  const DetLayer* tec9neg = negTecLayers[negTecLayers.size() - 1];
808  const LayerMeasurements layerMeasurements{*measurementTrackerHandle, *measurementTrackerEvent};
809  const TrajectoryStateOnSurface tsosTEC9 = itm->updatedState();
810 
811  // check if track on positive or negative z
812  if (!(iidd == StripSubdetector::TEC))
813  LogDebug("SiStripHitEfficiency:HitEff") << "there is a problem with TEC 9 extrapolation" << endl;
814 
815  //LOGPRINT << " tec9 id = " << iidd << " and side = " << tTopo->tecSide(iidd) << endl;
816  vector<TrajectoryMeasurement> tmp;
817  if (tTopo->tecSide(iidd) == 1) {
818  tmp = layerMeasurements.measurements(*tec9neg, tsosTEC9, *thePropagator, *estimator);
819  //LOGPRINT << "on negative side" << endl;
820  }
821  if (tTopo->tecSide(iidd) == 2) {
822  tmp = layerMeasurements.measurements(*tec9pos, tsosTEC9, *thePropagator, *estimator);
823  //LOGPRINT << "on positive side" << endl;
824  }
825 
826  if (!tmp.empty()) {
827  // take the last of the TMs, which is always an invalid hit
828  // if no detId is available, ie detId==0, then no compatible layer was crossed
829  // otherwise, use that TM for the efficiency measurement
830  TrajectoryMeasurement tec9TM(tmp.back());
831  const auto& tec9Hit = tec9TM.recHit();
832 
833  unsigned int tec9id = tec9Hit->geographicalId().rawId();
834  LogDebug("SiStripHitEfficiency:HitEff")
835  << "tec9id = " << tec9id << " is Double sided = " << ::isDoubleSided(tec9id, tTopo)
836  << " and 0x3 = " << (tec9id & 0x3) << endl;
837 
838  if (tec9Hit->geographicalId().rawId() != 0) {
839  LogDebug("SiStripHitEfficiency:HitEff") << "tec9 hit actually being added to TM vector" << endl;
840  // in tec the hit can be single or doubled sided. whenever the invalid hit at the end of vector of TMs is
841  // double sided it is always on the matched surface, so we need to split it into the true sensor surfaces
842  if (::isDoubleSided(tec9id, tTopo)) {
843  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator, 1));
844  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator, 2));
845  } else
846  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator));
847  }
848  } //else LOGPRINT << "tec9 tmp empty" << endl;
849  }
850  hitTotalCounters[TKlayers] += 1;
851 
853 
854  // Modules Constraints
855 
856  for (std::vector<TrajectoryAtInvalidHit>::const_iterator TM = TMs.begin(); TM != TMs.end(); ++TM) {
857  // --> Get trajectory from combinatedState
858  iidd = TM->monodet_id();
859  LogDebug("SiStripHitEfficiency:HitEff") << "setting iidd = " << iidd << " before checking efficiency and ";
860 
861  xloc = TM->localX();
862  yloc = TM->localY();
863 
864  angleX = atan(TM->localDxDz());
865  angleY = atan(TM->localDyDz());
866 
867  TrajLocErrX = 0.0;
868  TrajLocErrY = 0.0;
869 
870  xglob = TM->globalX();
871  yglob = TM->globalY();
872  zglob = TM->globalZ();
873  xErr = TM->localErrorX();
874  yErr = TM->localErrorY();
875 
876  TrajGlbX = 0.0;
877  TrajGlbY = 0.0;
878  TrajGlbZ = 0.0;
879  withinAcceptance = TM->withinAcceptance();
880 
881  trajHitValid = TM->validHit();
882  int TrajStrip = -1;
883 
884  // reget layer from iidd here, to account for TOB 6 and TEC 9 TKlayers being off
885  TKlayers = ::checkLayer(iidd, tTopo);
886 
887  if ((layers == TKlayers) || (layers == 0)) { // Look at the layer not used to reconstruct the track
888  whatlayer = TKlayers;
889  LogDebug("SiStripHitEfficiency:HitEff") << "Looking at layer under study" << endl;
890  ModIsBad = 2;
891  Id = 0;
892  SiStripQualBad = 0;
893  run = 0;
894  event = 0;
895  TrajLocX = 0.0;
896  TrajLocY = 0.0;
897  TrajLocAngleX = -999.0;
898  TrajLocAngleY = -999.0;
899  ResX = 0.0;
900  ResXSig = 0.0;
901  ClusterLocX = 0.0;
902  ClusterLocY = 0.0;
903  ClusterLocErrX = 0.0;
904  ClusterLocErrY = 0.0;
905  ClusterStoN = 0.0;
906  bunchx = 0;
907  commonMode = -100;
908 
909  // RPhi RecHit Efficiency
910 
911  if (!input.empty()) {
912  LogDebug("SiStripHitEfficiency:HitEff") << "Checking clusters with size = " << input.size() << endl;
913  int nClusters = 0;
914  std::vector<std::vector<float> >
915  VCluster_info; //fill with X residual, X residual pull, local X, sig(X), local Y, sig(Y), StoN
916  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end();
917  DSViter++) {
918  // DSViter is a vector of SiStripClusters located on a single module
919  //if (DEBUG) LOGPRINT << "the ID from the DSViter = " << DSViter->id() << endl;
920  unsigned int ClusterId = DSViter->id();
921  if (ClusterId == iidd) {
922  LogDebug("SiStripHitEfficiency:HitEff")
923  << "found (ClusterId == iidd) with ClusterId = " << ClusterId << " and iidd = " << iidd << endl;
924  DetId ClusterDetId(ClusterId);
925  const StripGeomDetUnit* stripdet = (const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
926  const StripTopology& Topo = stripdet->specificTopology();
927 
928  float hbedge = 0.0;
929  float htedge = 0.0;
930  float hapoth = 0.0;
931  float uylfac = 0.0;
932  float uxlden = 0.0;
933  if (TKlayers >= 11) {
934  const BoundPlane& plane = stripdet->surface();
935  const TrapezoidalPlaneBounds* trapezoidalBounds(
936  dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
937  std::array<const float, 4> const& parameterTrap =
938  (*trapezoidalBounds).parameters(); // el bueno aqui
939  hbedge = parameterTrap[0];
940  htedge = parameterTrap[1];
941  hapoth = parameterTrap[3];
942  uylfac = (htedge - hbedge) / (htedge + hbedge) / hapoth;
943  uxlden = 1 + yloc * uylfac;
944  }
945 
946  // Need to know position of trajectory in strip number for selecting the right APV later
947  if (TrajStrip == -1) {
948  int nstrips = Topo.nstrips();
949  float pitch = stripdet->surface().bounds().width() / nstrips;
950  TrajStrip = xloc / pitch + nstrips / 2.0;
951  // Need additionnal corrections for endcap
952  if (TKlayers >= 11) {
953  float TrajLocXMid = xloc / (1 + (htedge - hbedge) * yloc / (htedge + hbedge) /
954  hapoth); // radialy extrapolated x loc position at middle
955  TrajStrip = TrajLocXMid / pitch + nstrips / 2.0;
956  }
957  //LOGPRINT<<" Layer "<<TKlayers<<" TrajStrip: "<<nstrips<<" "<<pitch<<" "<<TrajStrip<<endl;
958  }
959 
960  for (edmNew::DetSet<SiStripCluster>::const_iterator iter = DSViter->begin(); iter != DSViter->end();
961  ++iter) {
962  //iter is a single SiStripCluster
964  float res = (parameters.first.x() - xloc);
965  float sigma = ::checkConsistency(parameters, xloc, xErr);
966  // The consistency is probably more accurately measured with the Chi2MeasurementEstimator. To use it
967  // you need a TransientTrackingRecHit instead of the cluster
968  //theEstimator= new Chi2MeasurementEstimator(30);
969  //const Chi2MeasurementEstimator *theEstimator(100);
970  //theEstimator->estimate(TM->tsos(), TransientTrackingRecHit);
971 
972  if (TKlayers >= 11) {
973  res = parameters.first.x() - xloc / uxlden; // radialy extrapolated x loc position at middle
974  sigma = abs(res) /
975  sqrt(parameters.second.xx() + xErr * xErr / uxlden / uxlden +
976  yErr * yErr * xloc * xloc * uylfac * uylfac / uxlden / uxlden / uxlden / uxlden);
977  }
978 
979  siStripClusterInfo_.setCluster(*iter, ClusterId);
980  // signal to noise from SiStripClusterInfo not working in 225. I'll fix this after the interface
981  // redesign in 300 -ku
982  //float cluster_info[7] = {res, sigma, parameters.first.x(), sqrt(parameters.second.xx()), parameters.first.y(), sqrt(parameters.second.yy()), signal_to_noise};
983  std::vector<float> cluster_info;
984  cluster_info.push_back(res);
985  cluster_info.push_back(sigma);
986  cluster_info.push_back(parameters.first.x());
987  cluster_info.push_back(sqrt(parameters.second.xx()));
988  cluster_info.push_back(parameters.first.y());
989  cluster_info.push_back(sqrt(parameters.second.yy()));
990  cluster_info.push_back(siStripClusterInfo_.signalOverNoise());
991  VCluster_info.push_back(cluster_info);
992  nClusters++;
993  LogDebug("SiStripHitEfficiency:HitEff") << "Have ID match. residual = " << VCluster_info.back()[0]
994  << " res sigma = " << VCluster_info.back()[1] << endl;
995  LogDebug("SiStripHitEfficiency:HitEff")
996  << "trajectory measurement compatability estimate = " << (*itm).estimate() << endl;
997  LogDebug("SiStripHitEfficiency:HitEff")
998  << "hit position = " << parameters.first.x() << " hit error = " << sqrt(parameters.second.xx())
999  << " trajectory position = " << xloc << " traj error = " << xErr << endl;
1000  }
1001  }
1002  }
1003  float FinalResSig = 1000.0;
1004  float FinalCluster[7] = {1000.0, 1000.0, 0.0, 0.0, 0.0, 0.0, 0.0};
1005  if (nClusters > 0) {
1006  LogDebug("SiStripHitEfficiency:HitEff") << "found clusters > 0" << endl;
1007  if (nClusters > 1) {
1008  //get the smallest one
1009  vector<vector<float> >::iterator ires;
1010  for (ires = VCluster_info.begin(); ires != VCluster_info.end(); ires++) {
1011  if (abs((*ires)[1]) < abs(FinalResSig)) {
1012  FinalResSig = (*ires)[1];
1013  for (unsigned int i = 0; i < ires->size(); i++) {
1014  LogDebug("SiStripHitEfficiency:HitEff")
1015  << "filling final cluster. i = " << i << " before fill FinalCluster[i]=" << FinalCluster[i]
1016  << " and (*ires)[i] =" << (*ires)[i] << endl;
1017  FinalCluster[i] = (*ires)[i];
1018  LogDebug("SiStripHitEfficiency:HitEff")
1019  << "filling final cluster. i = " << i << " after fill FinalCluster[i]=" << FinalCluster[i]
1020  << " and (*ires)[i] =" << (*ires)[i] << endl;
1021  }
1022  }
1023  LogDebug("SiStripHitEfficiency:HitEff")
1024  << "iresidual = " << (*ires)[0] << " isigma = " << (*ires)[1]
1025  << " and FinalRes = " << FinalCluster[0] << endl;
1026  }
1027  } else {
1028  FinalResSig = VCluster_info.at(0)[1];
1029  for (unsigned int i = 0; i < VCluster_info.at(0).size(); i++) {
1030  FinalCluster[i] = VCluster_info.at(0)[i];
1031  }
1032  }
1033  VCluster_info.clear();
1034  }
1035 
1036  LogDebug("SiStripHitEfficiency:HitEff")
1037  << "Final residual in X = " << FinalCluster[0] << "+-" << (FinalCluster[0] / FinalResSig) << endl;
1038  LogDebug("SiStripHitEfficiency:HitEff") << "Checking location of trajectory: abs(yloc) = " << abs(yloc)
1039  << " abs(xloc) = " << abs(xloc) << endl;
1040  LogDebug("SiStripHitEfficiency:HitEff")
1041  << "Checking location of cluster hit: yloc = " << FinalCluster[4] << "+-" << FinalCluster[5]
1042  << " xloc = " << FinalCluster[2] << "+-" << FinalCluster[3] << endl;
1043  LogDebug("SiStripHitEfficiency:HitEff") << "Final cluster signal to noise = " << FinalCluster[6] << endl;
1044 
1045  float exclusionWidth = 0.4;
1046  float TOBexclusion = 0.0;
1047  float TECexRing5 = -0.89;
1048  float TECexRing6 = -0.56;
1049  float TECexRing7 = 0.60;
1050  //Added by Chris Edelmaier to do TEC bonding exclusion
1051  int subdetector = ((iidd >> 25) & 0x7);
1052  int ringnumber = ((iidd >> 5) & 0x7);
1053 
1054  //New TOB and TEC bonding region exclusion zone
1055  if ((TKlayers >= 5 && TKlayers < 11) ||
1056  ((subdetector == 6) && ((ringnumber >= 5) && (ringnumber <= 7)))) {
1057  //There are only 2 cases that we need to exclude for
1058  float highzone = 0.0;
1059  float lowzone = 0.0;
1060  float higherr = yloc + 5.0 * yErr;
1061  float lowerr = yloc - 5.0 * yErr;
1062  if (TKlayers >= 5 && TKlayers < 11) {
1063  //TOB zone
1064  highzone = TOBexclusion + exclusionWidth;
1065  lowzone = TOBexclusion - exclusionWidth;
1066  } else if (ringnumber == 5) {
1067  //TEC ring 5
1068  highzone = TECexRing5 + exclusionWidth;
1069  lowzone = TECexRing5 - exclusionWidth;
1070  } else if (ringnumber == 6) {
1071  //TEC ring 6
1072  highzone = TECexRing6 + exclusionWidth;
1073  lowzone = TECexRing6 - exclusionWidth;
1074  } else if (ringnumber == 7) {
1075  //TEC ring 7
1076  highzone = TECexRing7 + exclusionWidth;
1077  lowzone = TECexRing7 - exclusionWidth;
1078  }
1079  //Now that we have our exclusion region, we just have to properly identify it
1080  if ((highzone <= higherr) && (highzone >= lowerr))
1081  withinAcceptance = false;
1082  if ((lowzone >= lowerr) && (lowzone <= higherr))
1083  withinAcceptance = false;
1084  if ((higherr <= highzone) && (higherr >= lowzone))
1085  withinAcceptance = false;
1086  if ((lowerr >= lowzone) && (lowerr <= highzone))
1087  withinAcceptance = false;
1088  }
1089 
1090  // fill ntuple varibles
1091  //get global position from module id number iidd
1092  TrajGlbX = xglob;
1093  TrajGlbY = yglob;
1094  TrajGlbZ = zglob;
1095 
1096  TrajLocErrX = xErr;
1097  TrajLocErrY = yErr;
1098 
1099  Id = iidd;
1100  run = run_nr;
1101  event = ev_nr;
1102  bunchx = bunch_nr;
1103  //if ( SiStripQuality_->IsModuleBad(iidd) ) {
1104  if (SiStripQuality_->getBadApvs(iidd) != 0) {
1105  SiStripQualBad = 1;
1106  LogDebug("SiStripHitEfficiency:HitEff") << "strip is bad from SiStripQuality" << endl;
1107  } else {
1108  SiStripQualBad = 0;
1109  LogDebug("SiStripHitEfficiency:HitEff") << "strip is good from SiStripQuality" << endl;
1110  }
1111 
1112  //check for FED-detected errors and include those in SiStripQualBad
1113  for (unsigned int ii = 0; ii < fedErrorIds.size(); ii++) {
1114  if (iidd == fedErrorIds[ii].rawId())
1115  SiStripQualBad = 1;
1116  }
1117 
1118  TrajLocX = xloc;
1119  TrajLocY = yloc;
1120  TrajLocAngleX = angleX;
1121  TrajLocAngleY = angleY;
1122  ResX = FinalCluster[0];
1123  ResXSig = FinalResSig;
1124  if (FinalResSig != FinalCluster[1])
1125  LogDebug("SiStripHitEfficiency:HitEff")
1126  << "Problem with best cluster selection because FinalResSig = " << FinalResSig
1127  << " and FinalCluster[1] = " << FinalCluster[1] << endl;
1128  ClusterLocX = FinalCluster[2];
1129  ClusterLocY = FinalCluster[4];
1130  ClusterLocErrX = FinalCluster[3];
1131  ClusterLocErrY = FinalCluster[5];
1132  ClusterStoN = FinalCluster[6];
1133 
1134  // CM of APV crossed by traj
1135  if (addCommonMode_)
1136  if (commonModeDigis.isValid() && TrajStrip >= 0 && TrajStrip <= 768) {
1137  edm::DetSetVector<SiStripRawDigi>::const_iterator digiframe = commonModeDigis->find(iidd);
1138  if (digiframe != commonModeDigis->end())
1139  if ((unsigned)TrajStrip / 128 < digiframe->data.size())
1140  commonMode = digiframe->data.at(TrajStrip / 128).adc();
1141  }
1142 
1143  LogDebug("SiStripHitEfficiency:HitEff") << "before check good" << endl;
1144 
1145  if (FinalResSig < 999.0) { //could make requirement on track/hit consistency, but for
1146  //now take anything with a hit on the module
1147  LogDebug("SiStripHitEfficiency:HitEff")
1148  << "hit being counted as good " << FinalCluster[0] << " FinalRecHit " << iidd << " TKlayers "
1149  << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd
1150  << " matched/stereo/rphi = " << ((iidd & 0x3) == 0) << "/" << ((iidd & 0x3) == 1) << "/"
1151  << ((iidd & 0x3) == 2) << endl;
1152  ModIsBad = 0;
1153  traj->Fill();
1154  } else {
1155  LogDebug("SiStripHitEfficiency:HitEff")
1156  << "hit being counted as bad ######### Invalid RPhi FinalResX " << FinalCluster[0]
1157  << " FinalRecHit " << iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc
1158  << " module " << iidd << " matched/stereo/rphi = " << ((iidd & 0x3) == 0) << "/"
1159  << ((iidd & 0x3) == 1) << "/" << ((iidd & 0x3) == 2) << endl;
1160  ModIsBad = 1;
1161  traj->Fill();
1162 
1163  LogDebug("SiStripHitEfficiency:HitEff") << " RPhi Error " << sqrt(xErr * xErr + yErr * yErr)
1164  << " ErrorX " << xErr << " yErr " << yErr << endl;
1165  }
1166  LogDebug("SiStripHitEfficiency:HitEff") << "after good location check" << endl;
1167  }
1168  LogDebug("SiStripHitEfficiency:HitEff") << "after list of clusters" << endl;
1169  }
1170  LogDebug("SiStripHitEfficiency:HitEff") << "After layers=TKLayers if" << endl;
1171  }
1172  LogDebug("SiStripHitEfficiency:HitEff") << "After looping over TrajAtValidHit list" << endl;
1173  }
1174  LogDebug("SiStripHitEfficiency:HitEff") << "end TMeasurement loop" << endl;
1175  }
1176  LogDebug("SiStripHitEfficiency:HitEff") << "end of trajectories loop" << endl;
1177  }
1178 }
1179 
1181  traj->GetDirectory()->cd();
1182  traj->Write();
1183 
1184  LogDebug("SiStripHitEfficiency:HitEff") << " Events Analysed " << events << endl;
1185  LogDebug("SiStripHitEfficiency:HitEff") << " Number Of Tracked events " << EventTrackCKF << endl;
1186 
1187  if (doMissingHitsRecovery_) {
1188  float totTIB = 0.0;
1189  float totTOB = 0.0;
1190  float totTID = 0.0;
1191  float totTEC = 0.0;
1192 
1193  float totTIBrepro = 0.0;
1194  float totTOBrepro = 0.0;
1195  float totTIDrepro = 0.0;
1196  float totTECrepro = 0.0;
1197 
1198  edm::LogInfo("SiStripHitEfficiency:HitEff") << "Within TIB :";
1199  for (int i = 0; i <= k_LayersAtTIBEnd; i++) {
1200  edm::LogInfo("SiStripHitEfficiency:HitEff")
1201  << "Layer " << i << " has : " << missHitPerLayer[i] << "/" << totalNbHits << " = "
1202  << (missHitPerLayer[i] * 1.0 / totalNbHits) * 100 << " % of missing hit";
1203  totTIB += missHitPerLayer[i];
1204  edm::LogInfo("SiStripHitEfficiency:HitEff")
1205  << "Removing recovered hits : layer " << i << " has : " << missHitPerLayer[i] - hitRecoveryCounters[i] << "/"
1206  << totalNbHits << " = " << ((missHitPerLayer[i] - hitRecoveryCounters[i]) * 1.0 / totalNbHits) * 100
1207  << " % of missing hit";
1208  totTIBrepro += (missHitPerLayer[i] - hitRecoveryCounters[i]);
1209  }
1210  edm::LogInfo("SiStripHitEfficiency:HitEff")
1211  << "TOTAL % of missing hits within TIB :" << (totTIB * 1.0 / totalNbHits) * 100 << "%";
1212  edm::LogInfo("SiStripHitEfficiency:HitEff")
1213  << "AFTER repropagation :" << (totTIBrepro * 1.0 / totalNbHits) * 100 << "%";
1214 
1215  edm::LogInfo("SiStripHitEfficiency:HitEff") << "Within TOB :";
1216  for (int i = k_LayersAtTIBEnd + 1; i <= k_LayersAtTOBEnd; i++) {
1217  edm::LogInfo("SiStripHitEfficiency:HitEff")
1218  << "Layer " << i << " has : " << missHitPerLayer[i] << "/" << totalNbHits << " = "
1219  << (missHitPerLayer[i] * 1.0 / totalNbHits) * 100 << " % of missing hit";
1220  totTOB += missHitPerLayer[i];
1221  edm::LogInfo("SiStripHitEfficiency:HitEff")
1222  << "Removing recovered hits : layer " << i << " has : " << missHitPerLayer[i] - hitRecoveryCounters[i] << "/"
1223  << totalNbHits << " = " << ((missHitPerLayer[i] - hitRecoveryCounters[i]) * 1.0 / totalNbHits) * 100
1224  << " % of missing hit";
1225  totTOBrepro += (missHitPerLayer[i] - hitRecoveryCounters[i]);
1226  }
1227  edm::LogInfo("SiStripHitEfficiency:HitEff")
1228  << "TOTAL % of missing hits within TOB :" << (totTOB * 1.0 / totalNbHits) * 100 << "%";
1229  edm::LogInfo("SiStripHitEfficiency:HitEff")
1230  << "AFTER repropagation :" << (totTOBrepro * 1.0 / totalNbHits) * 100 << "%";
1231 
1232  edm::LogInfo("SiStripHitEfficiency:HitEff") << "Within TID :";
1233  for (int i = k_LayersAtTOBEnd + 1; i <= k_LayersAtTIDEnd; i++) {
1234  edm::LogInfo("SiStripHitEfficiency:HitEff")
1235  << "Layer " << i << " has : " << missHitPerLayer[i] << "/" << totalNbHits << " = "
1236  << (missHitPerLayer[i] * 1.0 / totalNbHits) * 100 << " % of missing hit";
1237  totTID += missHitPerLayer[i];
1238  edm::LogInfo("SiStripHitEfficiency:HitEff")
1239  << "Removing recovered hits : layer " << i << " has : " << missHitPerLayer[i] - hitRecoveryCounters[i] << "/"
1240  << totalNbHits << " = " << ((missHitPerLayer[i] - hitRecoveryCounters[i]) * 1.0 / totalNbHits) * 100
1241  << " % of missing hit";
1242  totTIDrepro += (missHitPerLayer[i] - hitRecoveryCounters[i]);
1243  }
1244  edm::LogInfo("SiStripHitEfficiency:HitEff")
1245  << "TOTAL % of missing hits within TID :" << (totTID * 1.0 / totalNbHits) * 100 << "%";
1246  edm::LogInfo("SiStripHitEfficiency:HitEff")
1247  << "AFTER repropagation :" << (totTIDrepro * 1.0 / totalNbHits) * 100 << "%";
1248 
1249  edm::LogInfo("SiStripHitEfficiency:HitEff") << "Within TEC :";
1250  for (int i = k_LayersAtTIDEnd + 1; i < k_END_OF_LAYERS; i++) {
1251  edm::LogInfo("SiStripHitEfficiency:HitEff")
1252  << "Layer " << i << " has : " << missHitPerLayer[i] << "/" << totalNbHits << " = "
1253  << (missHitPerLayer[i] * 1.0 / totalNbHits) * 100 << " % of missing hit";
1254  totTEC += missHitPerLayer[i];
1255  edm::LogInfo("SiStripHitEfficiency:HitEff")
1256  << "Removing recovered hits : layer " << i << " has : " << missHitPerLayer[i] - hitRecoveryCounters[i] << "/"
1257  << totalNbHits << " = " << ((missHitPerLayer[i] - hitRecoveryCounters[i]) * 1.0 / totalNbHits) * 100
1258  << " % of missing hit";
1259  totTECrepro += (missHitPerLayer[i] - hitRecoveryCounters[i]);
1260  }
1261  edm::LogInfo("SiStripHitEfficiency:HitEff")
1262  << "TOTAL % of missing hits within TEC :" << (totTEC * 1.0 / totalNbHits) * 100 << "%";
1263  edm::LogInfo("SiStripHitEfficiency:HitEff")
1264  << "AFTER repropagation :" << (totTECrepro * 1.0 / totalNbHits) * 100 << "%";
1265 
1266  edm::LogInfo("SiStripHitEfficiency:HitEff") << " Hit recovery summary:";
1267 
1268  for (int ilayer = 0; ilayer < k_END_OF_LAYERS; ilayer++) {
1269  edm::LogInfo("SiStripHitEfficiency:HitEff")
1270  << " layer " << ilayer << ": " << hitRecoveryCounters[ilayer] << " / " << hitTotalCounters[ilayer];
1271  }
1272  }
1273 }
1274 
1275 //define this as a plug-in
std::pair< LocalPoint, LocalError > LocalValues
static const std::string kSharedResource
Definition: TFileService.h:76
static constexpr auto TEC
virtual int nstrips() const =0
float TrajGlbZ
Definition: HitEff.h:123
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
Definition: HitEff.h:94
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
float ClusterLocErrX
Definition: HitEff.h:126
std::vector< unsigned int > hitRecoveryCounters
Definition: HitEff.h:109
float ClusterLocErrY
Definition: HitEff.h:126
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
const edm::EDGetTokenT< OnlineLuminosityRecord > metaDataToken_
Definition: HitEff.h:66
void endJob() override
Definition: HitEff.cc:1180
void setCluster(const SiStripCluster &cluster, int detId)
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
float TrajGlbY
Definition: HitEff.h:123
std::vector< unsigned int > hitTotalCounters
Definition: HitEff.h:110
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > topoToken_
Definition: HitEff.h:90
float TrajLocX
Definition: HitEff.h:124
float ResXSig
Definition: HitEff.h:127
virtual void localParameters(AClusters const &clusters, ALocalValues &retValues, const GeomDetUnit &gd, const LocalTrajectoryParameters &ltp) const
bool useAllHitsFromTracksWithMissingHits_
Definition: HitEff.h:77
Definition: HitEff.h:53
unsigned int tidSide(const DetId &id) const
TString subdetector
int compSettings
Definition: HitEff.h:104
bool doMissingHitsRecovery_
Definition: HitEff.h:78
float instLumi
Definition: HitEff.h:137
const edm::EDGetTokenT< std::vector< Trajectory > > trajectories_token_
Definition: HitEff.h:81
T const * product() const
Definition: Handle.h:70
Class to contain the online luminosity from soft FED 1022.
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
const edm::EDGetTokenT< DetIdVector > digisVec_token_
Definition: HitEff.h:85
HitEff(const edm::ParameterSet &conf)
Definition: HitEff.cc:74
SiStripClusterInfo siStripClusterInfo_
Definition: HitEff.h:69
float ClusterStoN
Definition: HitEff.h:126
data_type const * const_iterator
Definition: DetSetNew.h:31
float TrajLocErrX
Definition: HitEff.h:125
float instLumi() const
Return the luminosity for the current nibble.
float commonMode
Definition: HitEff.h:138
const edm::EDGetTokenT< DetIdCollection > digisCol_token_
Definition: HitEff.h:84
const edm::ESGetToken< SiStripQuality, SiStripQualityRcd > siStripQualityToken_
Definition: HitEff.h:93
Definition: Electron.h:6
missing
Definition: combine.py:5
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
bool DEBUG
Definition: HitEff.h:106
const edm::EDGetTokenT< MeasurementTrackerEvent > trackerEvent_token_
Definition: HitEff.h:86
unsigned int bunchx
Definition: HitEff.h:135
static std::string const input
Definition: EdmProvDump.cc:50
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
Definition: HitEff.h:91
const_iterator end() const
last iterator over the map (read only)
float signalOverNoise() const
T getUntrackedParameter(std::string const &, T const &) const
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: HitEff.cc:177
bool highPurity
Definition: HitEff.h:132
size_type size() const
Return the number of contained DetSets.
Definition: DetSetVector.h:250
const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clusters_token_
Definition: HitEff.h:83
unsigned int tecSide(const DetId &id) const
float PU
Definition: HitEff.h:137
int nHits
Definition: HitEff.h:133
unsigned int Id
Definition: HitEff.h:129
const edm::EDGetTokenT< reco::TrackCollection > combinatorialTracks_token_
Definition: HitEff.h:80
T sqrt(T t)
Definition: SSEVec.h:19
float timeAtIpInOutErr
Definition: MuonTime.h:14
const edm::ESGetToken< Propagator, TrackingComponentsRecord > propagatorToken_
Definition: HitEff.h:97
int nDof
number of muon stations used
Definition: MuonTime.h:9
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int trajHitValid
Definition: HitEff.h:135
unsigned int ModIsBad
Definition: HitEff.h:128
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float TrajGlbX
Definition: HitEff.h:123
const edm::ESGetToken< MeasurementTracker, CkfComponentsRecord > measurementTkToken_
Definition: HitEff.h:95
bool withinAcceptance
Definition: HitEff.h:131
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
unsigned int run
Definition: HitEff.h:135
float ClusterLocX
Definition: HitEff.h:126
ii
Definition: cuy.py:589
int events
Definition: HitEff.h:102
float TrajLocAngleX
Definition: HitEff.h:124
const edm::ESGetToken< Chi2MeasurementEstimatorBase, TrackingComponentsRecord > chi2MeasurementEstimatorToken_
Definition: HitEff.h:96
Log< level::Info, false > LogInfo
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
Definition: DetId.h:17
int totalNbHits
Definition: HitEff.h:121
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
bool useLastMeas_
Definition: HitEff.h:76
unsigned int whatlayer
Definition: HitEff.h:107
const edm::ESGetToken< StripClusterParameterEstimator, TkStripCPERecord > cpeToken_
Definition: HitEff.h:92
virtual const std::array< const float, 4 > parameters() const
unsigned int trackMultiplicityCut_
Definition: HitEff.h:74
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
unsigned int layers
Definition: HitEff.h:105
float pT
Definition: HitEff.h:134
float ResX
Definition: HitEff.h:127
void beginJob() override
Definition: HitEff.cc:115
bool addLumi_
Definition: HitEff.h:71
bool isValid() const
Definition: HandleBase.h:70
unsigned int SiStripQualBad
Definition: HitEff.h:130
bool addCommonMode_
Definition: HitEff.h:72
const_iterator begin() const
first iterator over the map (read only)
const edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackAsso_token_
Definition: HitEff.h:82
fixed size matrix
HLT enums.
float TrajLocErrY
Definition: HitEff.h:125
float TrajLocAngleY
Definition: HitEff.h:124
std::vector< LumiScalers > LumiScalersCollection
Definition: LumiScalers.h:144
bool cutOnTracks_
Definition: HitEff.h:73
bool useFirstMeas_
Definition: HitEff.h:75
Log< level::Warning, false > LogWarning
float ClusterLocY
Definition: HitEff.h:126
void initEvent(const edm::EventSetup &iSetup)
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
A Digi for the silicon strip detector, containing only adc information, and suitable for storing raw ...
tmp
align.sh
Definition: createJobs.py:716
const edm::EDGetTokenT< LumiScalersCollection > scalerToken_
Definition: HitEff.h:65
TTree * traj
Definition: HitEff.h:101
virtual float width() const =0
const edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > commonModeToken_
Definition: HitEff.h:67
edm::ParameterSet conf_
Definition: HitEff.h:99
float timeAtIpInOut
Definition: MuonTime.h:13
int EventTrackCKF
Definition: HitEff.h:102
float avgPileUp() const
Return the average pileup for th current nibble.
float TrajLocY
Definition: HitEff.h:124
int tquality
Definition: HitEff.h:136
std::vector< int > missHitPerLayer
Definition: HitEff.h:122
Definition: event.py:1
ConstRecHitPointer const & recHit() const
#define LogDebug(id)
short getBadApvs(uint32_t detid) const
const Bounds & bounds() const
Definition: Surface.h:87