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