CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HitEff Class Reference

#include <HitEff.h>

Inheritance diagram for HitEff:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

bool check2DPartner (unsigned int iidd, const std::vector< TrajectoryMeasurement > &traj)
 
double checkConsistency (const StripClusterParameterEstimator::LocalValues &parameters, double xx, double xerr)
 
unsigned int checkLayer (unsigned int iidd, const TrackerTopology *tTopo)
 
 HitEff (const edm::ParameterSet &conf)
 
bool isDoubleSided (unsigned int iidd, const TrackerTopology *tTopo) const
 
 ~HitEff () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c) override
 
void beginJob () override
 
void endJob () override
 

Private Attributes

bool addCommonMode_
 
bool addLumi_
 
unsigned int bunchx
 
float ClusterLocErrX
 
float ClusterLocErrY
 
float ClusterLocX
 
float ClusterLocY
 
const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clusters_token_
 
float ClusterStoN
 
const edm::EDGetTokenT< reco::TrackCollectioncombinatorialTracks_token_
 
float commonMode
 
const edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > commonModeToken_
 
int compSettings
 
edm::ParameterSet conf_
 
bool cutOnTracks_
 
bool DEBUG
 
const edm::EDGetTokenT< DetIdCollectiondigis_token_
 
unsigned int event
 
int events
 
int EventTrackCKF
 
bool highPurity
 
unsigned int Id
 
float instLumi
 
unsigned int layers
 
unsigned int ModIsBad
 
int nHits
 
float pT
 
float PU
 
float ResX
 
float ResXSig
 
unsigned int run
 
const edm::EDGetTokenT< LumiScalersCollectionscalerToken_
 
unsigned int SiStripQualBad
 
int tquality
 
const edm::EDGetTokenT< MeasurementTrackerEventtrackerEvent_token_
 
unsigned int trackMultiplicityCut_
 
TTree * traj
 
const edm::EDGetTokenT< std::vector< Trajectory > > trajectories_token_
 
float TrajGlbX
 
float TrajGlbY
 
float TrajGlbZ
 
unsigned int trajHitValid
 
float TrajLocAngleX
 
float TrajLocAngleY
 
float TrajLocErrX
 
float TrajLocErrY
 
float TrajLocX
 
float TrajLocY
 
const edm::EDGetTokenT< TrajTrackAssociationCollectiontrajTrackAsso_token_
 
bool useAllHitsFromTracksWithMissingHits_
 
bool useFirstMeas_
 
bool useLastMeas_
 
unsigned int whatlayer
 
bool withinAcceptance
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 48 of file HitEff.h.

Constructor & Destructor Documentation

HitEff::HitEff ( const edm::ParameterSet conf)
explicit

Definition at line 74 of file HitEff.cc.

References addCommonMode_, addLumi_, compSettings, conf_, cutOnTracks_, DEBUG, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), layers, trackMultiplicityCut_, useAllHitsFromTracksWithMissingHits_, useFirstMeas_, and useLastMeas_.

74  :
75  scalerToken_( consumes< LumiScalersCollection >(conf.getParameter<edm::InputTag>("lumiScalers")) ),
77  combinatorialTracks_token_( consumes< reco::TrackCollection >(conf.getParameter<edm::InputTag>("combinatorialTracks")) ),
78  trajectories_token_( consumes< std::vector<Trajectory> >(conf.getParameter<edm::InputTag>("trajectories")) ),
79  trajTrackAsso_token_( consumes< TrajTrackAssociationCollection >(conf.getParameter<edm::InputTag>("trajectories")) ),
81  digis_token_( consumes< DetIdCollection >(conf.getParameter<edm::InputTag>("siStripDigis")) ),
82  trackerEvent_token_( consumes< MeasurementTrackerEvent>(conf.getParameter<edm::InputTag>("trackerEvent")) ),
83  conf_(conf)
84 {
85  compSettings=conf_.getUntrackedParameter<int>("CompressionSettings",-1);
86  layers =conf_.getParameter<int>("Layer");
87  DEBUG = conf_.getParameter<bool>("Debug");
88  addLumi_ = conf_.getUntrackedParameter<bool>("addLumi", false);
89  addCommonMode_ = conf_.getUntrackedParameter<bool>("addCommonMode", false);
90  cutOnTracks_ = conf_.getUntrackedParameter<bool>("cutOnTracks", false);
91  trackMultiplicityCut_ = conf.getUntrackedParameter<unsigned int>("trackMultiplicity",100);
92  useFirstMeas_ = conf_.getUntrackedParameter<bool>("useFirstMeas", false);
93  useLastMeas_ = conf_.getUntrackedParameter<bool>("useLastMeas", false);
94  useAllHitsFromTracksWithMissingHits_ = conf_.getUntrackedParameter<bool>("useAllHitsFromTracksWithMissingHits", false);
95 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
bool useAllHitsFromTracksWithMissingHits_
Definition: HitEff.h:73
int compSettings
Definition: HitEff.h:87
bool DEBUG
Definition: HitEff.h:89
const edm::EDGetTokenT< MeasurementTrackerEvent > trackerEvent_token_
Definition: HitEff.h:80
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EDGetTokenT< DetIdCollection > digis_token_
Definition: HitEff.h:79
const edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackAsso_token_
Definition: HitEff.h:77
bool useLastMeas_
Definition: HitEff.h:72
unsigned int trackMultiplicityCut_
Definition: HitEff.h:70
unsigned int layers
Definition: HitEff.h:88
const edm::EDGetTokenT< reco::TrackCollection > combinatorialTracks_token_
Definition: HitEff.h:75
bool addLumi_
Definition: HitEff.h:67
bool addCommonMode_
Definition: HitEff.h:68
bool cutOnTracks_
Definition: HitEff.h:69
bool useFirstMeas_
Definition: HitEff.h:71
const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clusters_token_
Definition: HitEff.h:78
const edm::EDGetTokenT< LumiScalersCollection > scalerToken_
Definition: HitEff.h:64
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
const edm::EDGetTokenT< std::vector< Trajectory > > trajectories_token_
Definition: HitEff.h:76
const edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > commonModeToken_
Definition: HitEff.h:65
edm::ParameterSet conf_
Definition: HitEff.h:82
HitEff::~HitEff ( )
override

Definition at line 98 of file HitEff.cc.

98 { }

Member Function Documentation

void HitEff::analyze ( const edm::Event e,
const edm::EventSetup c 
)
overrideprivate

Definition at line 161 of file HitEff.cc.

References funct::abs(), addCommonMode_, addLumi_, AnalyticalPropagator_cfi::AnalyticalPropagator, anyDirection, edm::AssociationMap< Tag >::begin(), edmNew::DetSetVector< T >::begin(), Surface::bounds(), edm::EventBase::bunchCrossing(), bunchx, check2DPartner(), checkConsistency(), checkLayer(), vertices_cff::chi2, ClusterLocErrX, ClusterLocErrY, ClusterLocX, ClusterLocY, clusters_token_, ClusterStoN, combinatorialTracks_token_, commonMode, commonModeToken_, edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > >::const_iterator, cutOnTracks_, digis_token_, edmNew::DetSetVector< T >::empty(), edm::AssociationMap< Tag >::end(), edmNew::DetSetVector< T >::end(), edm::EventID::event(), events, EventTrackCKF, edm::EventSetup::get(), SiStripQuality::getBadApvs(), edm::Event::getByLabel(), edm::Event::getByToken(), highPurity, muons_cff::highPurity, mps_fire::i, edm::EventBase::id(), Id, TrackerGeometry::idToDetUnit(), cuy::ii, input, instLumi, ires, isDoubleSided(), edm::HandleBase::isValid(), edm::Ref< C, T, F >::key(), layers, StripClusterParameterEstimator::localParameters(), LogDebug, LayerMeasurements::measurements(), ecalDrivenElectronSeedsParameters_cff::measurementTrackerEvent, ModIsBad, reco::MuonTime::nDof, nHits, StripTopology::nstrips(), AlCaHLTBitMon_ParallelJobs::p, TrapezoidalPlaneBounds::parameters(), edm::Handle< T >::product(), edm::ESHandle< T >::product(), PhotonConversionTrajectorySeedProducerFromQuadruplets_cfi::propagator, pT, PU, rpcPointValidation_cfi::recHit, TrajectoryMeasurement::recHit(), ResX, ResXSig, edm::EventID::run(), run, scalerToken_, SiStripQualBad, edm::EDCollection< T >::size(), edm::DetSetVector< T >::size(), edmNew::DetSetVector< T >::size(), StripGeomDetUnit::specificTopology(), mathSSE::sqrt(), subdetector, GeomDet::surface(), StripSubdetector::TEC, TrackerTopology::tecSide(), reco::MuonTime::timeAtIpInOut, reco::MuonTime::timeAtIpInOutErr, tmp, trackingTruthProducer_cfi::tracker, trackerEvent_token_, trackMultiplicityCut_, traj, trajectories_token_, TrajGlbX, TrajGlbY, TrajGlbZ, trajHitValid, TrajLocAngleX, TrajLocAngleY, TrajLocErrX, TrajLocErrY, TrajLocX, TrajLocY, trajTrackAsso_token_, useAllHitsFromTracksWithMissingHits_, useFirstMeas_, useLastMeas_, whatlayer, Bounds::width(), and withinAcceptance.

161  {
162  //Retrieve tracker topology from geometry
163  edm::ESHandle<TrackerTopology> tTopoHandle;
164  es.get<TrackerTopologyRcd>().get(tTopoHandle);
165  const TrackerTopology* const tTopo = tTopoHandle.product();
166 
167  // bool DEBUG = false;
168 
169  LogDebug("SiStripHitEfficiency:HitEff") << "beginning analyze from HitEff" << endl;
170 
171  using namespace edm;
172  using namespace reco;
173  // Step A: Get Inputs
174 
175  int run_nr = e.id().run();
176  int ev_nr = e.id().event();
177  int bunch_nr = e.bunchCrossing();
178 
179  // Luminosity informations
181  instLumi=0; PU=0;
182  if(addLumi_) {
183  e.getByToken(scalerToken_, lumiScalers);
184  if (lumiScalers->begin() != lumiScalers->end()) {
185  instLumi = lumiScalers->begin()->instantLumi();
186  PU = lumiScalers->begin()->pileup();
187  }
188  }
189 
190  // CM
192  if(addCommonMode_) e.getByToken(commonModeToken_, commonModeDigis);
193 
194  //CombinatoriaTrack
195  edm::Handle<reco::TrackCollection> trackCollectionCKF;
196  //edm::InputTag TkTagCKF = conf_.getParameter<edm::InputTag>("combinatorialTracks");
197  e.getByToken(combinatorialTracks_token_,trackCollectionCKF);
198 
199  edm::Handle<std::vector<Trajectory> > TrajectoryCollectionCKF;
200  //edm::InputTag TkTrajCKF = conf_.getParameter<edm::InputTag>("trajectories");
201  e.getByToken(trajectories_token_,TrajectoryCollectionCKF);
202 
203  edm::Handle<TrajTrackAssociationCollection> trajTrackAssociationHandle;
204  e.getByToken(trajTrackAsso_token_, trajTrackAssociationHandle);
205 
206  // Clusters
207  // get the SiStripClusters from the event
209  //e.getByLabel("siStripClusters", theClusters);
210  e.getByToken(clusters_token_, theClusters);
211 
212  //get tracker geometry
214  es.get<TrackerDigiGeometryRecord>().get(tracker);
215  const TrackerGeometry * tkgeom=&(* tracker);
216 
217  //get Cluster Parameter Estimator
218  //std::string cpe = conf_.getParameter<std::string>("StripCPE");
220  es.get<TkStripCPERecord>().get("StripCPEfromTrackAngle", parameterestimator);
221  const StripClusterParameterEstimator &stripcpe(*parameterestimator);
222 
223  // get the SiStripQuality records
224  edm::ESHandle<SiStripQuality> SiStripQuality_;
225 //LQ commenting the try/catch that causes problem in 74X calibTree production
226 // try {
227 // es.get<SiStripQualityRcd>().get("forCluster",SiStripQuality_);
228 // }
229 // catch (...) {
230  es.get<SiStripQualityRcd>().get(SiStripQuality_);
231 // }
232 
233  edm::ESHandle<MagneticField> magFieldHandle;
234  es.get<IdealMagneticFieldRecord>().get(magFieldHandle);
235  const MagneticField* magField_ = magFieldHandle.product();
236 
237  // get the list of module IDs with FED-detected errors
238  edm::Handle< DetIdCollection > fedErrorIds;
239  //e.getByLabel("siStripDigis", fedErrorIds );
240  e.getByToken(digis_token_, fedErrorIds );
241 
242  ESHandle<MeasurementTracker> measurementTrackerHandle;
243  es.get<CkfComponentsRecord>().get(measurementTrackerHandle);
244 
246  //e.getByLabel("MeasurementTrackerEvent", measurementTrackerEvent);
247  e.getByToken(trackerEvent_token_,measurementTrackerEvent);
248 
250  es.get<TrackingComponentsRecord>().get("Chi2",est);
251 
253  es.get<TrackingComponentsRecord>().get("PropagatorWithMaterial",prop);
254  const Propagator* thePropagator = prop.product();
255 
256  events++;
257 
258  // *************** SiStripCluster Collection
259  const edmNew::DetSetVector<SiStripCluster>& input = *theClusters;
260 
261  //go through clusters to write out global position of good clusters for the layer understudy for comparison
262  // Loop through clusters just to print out locations
263  // Commented out to avoid discussion, should really be deleted.
264  /*
265  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
266  // DSViter is a vector of SiStripClusters located on a single module
267  unsigned int ClusterId = DSViter->id();
268  DetId ClusterDetId(ClusterId);
269  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
270 
271  edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
272  edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end();
273  for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
274  //iter is a single SiStripCluster
275  StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
276 
277  const Surface* surface;
278  surface = &(tracker->idToDet(ClusterDetId)->surface());
279  LocalPoint lp = parameters.first;
280  GlobalPoint gp = surface->toGlobal(lp);
281  unsigned int layer = checkLayer(ClusterId, tTopo);
282  if(DEBUG) cout << "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;
283  }
284  }
285  */
286 
287  // Tracking
288  const reco::TrackCollection *tracksCKF=trackCollectionCKF.product();
289  LogDebug("SiStripHitEfficiency:HitEff") << "number ckf tracks found = " << tracksCKF->size() << endl;
290  //if (tracksCKF->size() == 1 ){
291  if (!tracksCKF->empty()) {
292  if( cutOnTracks_ && (tracksCKF->size() >= trackMultiplicityCut_) ) return;
293  if( cutOnTracks_ ) LogDebug("SiStripHitEfficiency:HitEff") << "starting checking good event with < "<< trackMultiplicityCut_ <<" tracks" << endl;
294 
295  EventTrackCKF++;
296 
297 #ifdef ExtendedCALIBTree
298  //get dEdx info if available
299  Handle<ValueMap<DeDxData> > dEdxUncalibHandle;
300  if (e.getByLabel("dedxMedianCTF", dEdxUncalibHandle)) {
301  const ValueMap<DeDxData> dEdxTrackUncalib = *dEdxUncalibHandle.product();
302 
303  reco::TrackRef itTrack = reco::TrackRef( trackCollectionCKF, 0 );
304  dedx = dEdxTrackUncalib[itTrack].dEdx();
305  dedxNOM = dEdxTrackUncalib[itTrack].numberOfMeasurements();
306  } else {
307  dedx = -999.0; dedxNOM = -999;
308  }
309 
310 
311 
312  //get muon and ecal timing info if available
314  if(e.getByLabel("muonsWitht0Correction",muH)){
315  const MuonCollection & muonsT0 = *muH.product();
316  if(!muonsT0.empty()) {
317  MuonTime mt0 = muonsT0[0].time();
318  timeDT = mt0.timeAtIpInOut;
319  timeDTErr = mt0.timeAtIpInOutErr;
320  timeDTDOF = mt0.nDof;
321 
322  bool hasCaloEnergyInfo = muonsT0[0].isEnergyValid();
323  if (hasCaloEnergyInfo) timeECAL = muonsT0[0].calEnergy().ecal_time;
324  }
325  } else {
326  timeDT = -999.0; timeDTErr = -999.0; timeDTDOF = -999; timeECAL = -999.0;
327  }
328 
329 #endif
330  // actually should do a loop over all the tracks in the event here
331 
332  // Looping over traj-track associations to be able to get traj & track informations
333  for(TrajTrackAssociationCollection::const_iterator it = trajTrackAssociationHandle->begin();
334  it!=trajTrackAssociationHandle->end();
335  it++) {
336 
338  reco::TrackRef itrack = it->val;
339 
340  // for each track, fill some variables such as number of hits and momentum
341  nHits = itraj->foundHits();
342 #ifdef ExtendedCALIBTree
343  nLostHits = itraj->lostHits();
344  chi2 = (itraj->chiSquared()/itraj->ndof());
345  p = itraj->lastMeasurement().updatedState().globalMomentum().mag();
346 #endif
347  pT = sqrt( ( itraj->lastMeasurement().updatedState().globalMomentum().x() *
348  itraj->lastMeasurement().updatedState().globalMomentum().x()) +
349  ( itraj->lastMeasurement().updatedState().globalMomentum().y() *
350  itraj->lastMeasurement().updatedState().globalMomentum().y()) );
351 
352  // track quality
354 
355 
356 
357  std::vector<TrajectoryMeasurement> TMeas=itraj->measurements();
358  vector<TrajectoryMeasurement>::iterator itm;
359  double xloc = 0.;
360  double yloc = 0.;
361  double xErr = 0.;
362  double yErr = 0.;
363  double angleX = -999.;
364  double angleY = -999.;
365  double xglob,yglob,zglob;
366 
367 
368  // Check whether the trajectory has some missing hits
369  bool hasMissingHits=false;
370  for (itm=TMeas.begin();itm!=TMeas.end();itm++){
371  auto theHit = (*itm).recHit();
372  if(theHit->getType()==TrackingRecHit::Type::missing) hasMissingHits=true;
373  }
374 
375 
376  // Loop on each measurement and take it into consideration
377  //--------------------------------------------------------
378 
379  for (itm=TMeas.begin();itm!=TMeas.end();itm++){
380  auto theInHit = (*itm).recHit();
381 
382  LogDebug("SiStripHitEfficiency:HitEff") << "theInHit is valid = " << theInHit->isValid() << endl;
383 
384  unsigned int iidd = theInHit->geographicalId().rawId();
385 
386  unsigned int TKlayers = checkLayer(iidd, tTopo);
387  LogDebug("SiStripHitEfficiency:HitEff") << "TKlayer from trajectory: " << TKlayers << " from module = " << iidd << " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
388 
389 
390  // Test first and last points of the trajectory
391  // the list of measurements starts from outer layers !!! This could change -> should add a check
392  bool isFirstMeas = (itm==(TMeas.end()-1));
393  bool isLastMeas = (itm==(TMeas.begin()));
394 
395  if(!useFirstMeas_ && isFirstMeas) continue;
396  if(!useLastMeas_ && isLastMeas) continue;
397 
398  // In case of missing hit in the track, check whether to use the other hits or not.
399  if(hasMissingHits && theInHit->getType()!=TrackingRecHit::Type::missing && !useAllHitsFromTracksWithMissingHits_) continue;
400 
401 
402 
403  // If Trajectory measurement from TOB 6 or TEC 9, skip it because it's always valid they are filled later
404  if ( TKlayers == 10 || TKlayers == 22 ) {
405  LogDebug("SiStripHitEfficiency:HitEff") << "skipping original TM for TOB 6 or TEC 9" << endl;
406  continue;
407  }
408 
409  // Make vector of TrajectoryAtInvalidHits to hold the trajectories
410  std::vector<TrajectoryAtInvalidHit> TMs;
411 
412  // Make AnalyticalPropagator to use in TAVH constructor
414 
415  // for double sided layers check both sensors--if no hit was found on either sensor surface,
416  // the trajectory measurements only have one invalid hit entry on the matched surface
417  // so get the TrajectoryAtInvalidHit for both surfaces and include them in the study
418  if (isDoubleSided(iidd, tTopo) && ((iidd & 0x3)==0) ) {
419  // do hit eff check twice--once for each sensor
420  //add a TM for each surface
421  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 1));
422  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 2));
423  } else if ( isDoubleSided(iidd, tTopo) && (!check2DPartner(iidd, TMeas)) ) {
424  // if only one hit was found the trajectory measurement is on that sensor surface, and the other surface from
425  // the matched layer should be added to the study as well
426  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 1));
427  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 2));
428  LogDebug("SiStripHitEfficiency:HitEff") << " found a hit with a missing partner" << endl;
429  } else {
430  //only add one TM for the single surface and the other will be added in the next iteration
431  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator));
432  }
433 
435  //Now check for tracks at TOB6 and TEC9
436 
437  // to make sure we only propagate on the last TOB5 hit check the next entry isn't also in TOB5
438  // to avoid bias, make sure the TOB5 hit is valid (an invalid hit on TOB5 could only exist with a valid hit on TOB6)
439 
440  bool isValid = theInHit->isValid();
441  bool isLast = (itm==(TMeas.end()-1));
442  bool isLastTOB5 = true;
443  if (!isLast) {
444  if ( checkLayer((++itm)->recHit()->geographicalId().rawId(), tTopo) == 9 ) isLastTOB5 = false;
445  else isLastTOB5 = true;
446  --itm;
447  }
448 
449  if ( TKlayers==9 && isValid && isLastTOB5 ) {
450  // if ( TKlayers==9 && itm==TMeas.rbegin()) {
451  // if ( TKlayers==9 && (itm==TMeas.back()) ) { // to check for only the last entry in the trajectory for propagation
452  std::vector< BarrelDetLayer const*> barrelTOBLayers = measurementTrackerHandle->geometricSearchTracker()->tobLayers() ;
453  const DetLayer* tob6 = barrelTOBLayers[barrelTOBLayers.size()-1];
454  const MeasurementEstimator* estimator = est.product();
455  const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(*measurementTrackerHandle, *measurementTrackerEvent);
456  const TrajectoryStateOnSurface tsosTOB5 = itm->updatedState();
457  vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements(*tob6, tsosTOB5, *thePropagator, *estimator);
458 
459  if ( !tmp.empty()) {
460  LogDebug("SiStripHitEfficiency:HitEff") << "size of TM from propagation = " << tmp.size() << endl;
461 
462  // take the last of the TMs, which is always an invalid hit
463  // if no detId is available, ie detId==0, then no compatible layer was crossed
464  // otherwise, use that TM for the efficiency measurement
465  TrajectoryMeasurement tob6TM(tmp.back());
466  const auto& tob6Hit = tob6TM.recHit();
467 
468  if (tob6Hit->geographicalId().rawId()!=0) {
469  LogDebug("SiStripHitEfficiency:HitEff") << "tob6 hit actually being added to TM vector" << endl;
470  TMs.push_back(TrajectoryAtInvalidHit(tob6TM, tTopo, tkgeom, propagator));
471  }
472  }
473  }
474 
475  bool isLastTEC8 = true;
476  if (!isLast) {
477  if ( checkLayer((++itm)->recHit()->geographicalId().rawId(), tTopo) == 21 ) isLastTEC8 = false;
478  else isLastTEC8 = true;
479  --itm;
480  }
481 
482  if ( TKlayers==21 && isValid && isLastTEC8 ) {
483 
484  std::vector< const ForwardDetLayer*> posTecLayers = measurementTrackerHandle->geometricSearchTracker()->posTecLayers() ;
485  const DetLayer* tec9pos = posTecLayers[posTecLayers.size()-1];
486  std::vector< const ForwardDetLayer*> negTecLayers = measurementTrackerHandle->geometricSearchTracker()->negTecLayers() ;
487  const DetLayer* tec9neg = negTecLayers[negTecLayers.size()-1];
488 
489  const MeasurementEstimator* estimator = est.product();
490  const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(*measurementTrackerHandle, *measurementTrackerEvent);
491  const TrajectoryStateOnSurface tsosTEC9 = itm->updatedState();
492 
493  // check if track on positive or negative z
494  if (!(iidd == StripSubdetector::TEC)) LogDebug("SiStripHitEfficiency:HitEff") << "there is a problem with TEC 9 extrapolation" << endl;
495 
496  //cout << " tec9 id = " << iidd << " and side = " << tTopo->tecSide(iidd) << endl;
497  vector<TrajectoryMeasurement> tmp;
498  if ( tTopo->tecSide(iidd) == 1 ) {
499  tmp = theLayerMeasurements->measurements(*tec9neg, tsosTEC9, *thePropagator, *estimator);
500  //cout << "on negative side" << endl;
501  }
502  if ( tTopo->tecSide(iidd) == 2 ) {
503  tmp = theLayerMeasurements->measurements(*tec9pos, tsosTEC9, *thePropagator, *estimator);
504  //cout << "on positive side" << endl;
505  }
506 
507  if ( !tmp.empty()) {
508  // take the last of the TMs, which is always an invalid hit
509  // if no detId is available, ie detId==0, then no compatible layer was crossed
510  // otherwise, use that TM for the efficiency measurement
511  TrajectoryMeasurement tec9TM(tmp.back());
512  const auto& tec9Hit = tec9TM.recHit();
513 
514  unsigned int tec9id = tec9Hit->geographicalId().rawId();
515  LogDebug("SiStripHitEfficiency:HitEff") << "tec9id = " << tec9id << " is Double sided = " << isDoubleSided(tec9id, tTopo) << " and 0x3 = " << (tec9id & 0x3) << endl;
516 
517  if (tec9Hit->geographicalId().rawId()!=0) {
518  LogDebug("SiStripHitEfficiency:HitEff") << "tec9 hit actually being added to TM vector" << endl;
519  // in tec the hit can be single or doubled sided. whenever the invalid hit at the end of vector of TMs is
520  // double sided it is always on the matched surface, so we need to split it into the true sensor surfaces
521  if (isDoubleSided(tec9id, tTopo)) {
522  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator, 1));
523  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator, 2));
524  } else
525  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator));
526  }
527  } //else cout << "tec9 tmp empty" << endl;
528  }
529 
531 
532  // Modules Constraints
533 
534  for(std::vector<TrajectoryAtInvalidHit>::const_iterator TM=TMs.begin();TM!=TMs.end();++TM) {
535 
536  // --> Get trajectory from combinatedState
537  iidd = TM->monodet_id();
538  LogDebug("SiStripHitEfficiency:HitEff") << "setting iidd = " << iidd << " before checking efficiency and ";
539 
540  xloc = TM->localX();
541  yloc = TM->localY();
542 
543  angleX = atan( TM->localDxDz() );
544  angleY = atan( TM->localDyDz() );
545 
546  TrajLocErrX = 0.0; TrajLocErrY = 0.0;
547 
548  xglob = TM->globalX();
549  yglob = TM->globalY();
550  zglob = TM->globalZ();
551  xErr = TM->localErrorX();
552  yErr = TM->localErrorY();
553 
554  TrajGlbX = 0.0; TrajGlbY = 0.0; TrajGlbZ = 0.0;
555  withinAcceptance = TM->withinAcceptance();
556 
557  trajHitValid = TM->validHit();
558  int TrajStrip = -1;
559 
560  // reget layer from iidd here, to account for TOB 6 and TEC 9 TKlayers being off
561  TKlayers = checkLayer(iidd, tTopo);
562 
563  if ( (layers == TKlayers) || (layers == 0) ) { // Look at the layer not used to reconstruct the track
564  whatlayer = TKlayers;
565  LogDebug("SiStripHitEfficiency:HitEff") << "Looking at layer under study" << endl;
566  ModIsBad = 2; Id = 0; SiStripQualBad = 0;
567  run = 0; event = 0; TrajLocX = 0.0; TrajLocY = 0.0;
568  TrajLocAngleX = -999.0; TrajLocAngleY = -999.0; ResX = 0.0; ResXSig = 0.0;
569  ClusterLocX = 0.0; ClusterLocY = 0.0; ClusterLocErrX = 0.0; ClusterLocErrY = 0.0; ClusterStoN = 0.0;
570  bunchx = 0; commonMode = -100;
571 
572  // RPhi RecHit Efficiency
573 
574  if (!input.empty() ) {
575  LogDebug("SiStripHitEfficiency:HitEff") << "Checking clusters with size = " << input.size() << endl;
576  int nClusters = 0;
577  std::vector< std::vector<float> > VCluster_info; //fill with X residual, X residual pull, local X, sig(X), local Y, sig(Y), StoN
578  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
579  // DSViter is a vector of SiStripClusters located on a single module
580  //if (DEBUG) cout << "the ID from the DSViter = " << DSViter->id() << endl;
581  unsigned int ClusterId = DSViter->id();
582  if (ClusterId == iidd) {
583  LogDebug("SiStripHitEfficiency:HitEff") << "found (ClusterId == iidd) with ClusterId = " << ClusterId << " and iidd = " << iidd << endl;
584  DetId ClusterDetId(ClusterId);
585  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
586  const StripTopology& Topo = stripdet->specificTopology();
587 
588  float hbedge = 0.0;
589  float htedge = 0.0;
590  float hapoth = 0.0;
591  float uylfac = 0.0;
592  float uxlden = 0.0;
593  if(TKlayers>=11) {
594  const BoundPlane& plane = stripdet->surface();
595  const TrapezoidalPlaneBounds* trapezoidalBounds( dynamic_cast<const TrapezoidalPlaneBounds*>(&(plane.bounds())));
596  std::array<const float, 4> const & parameterTrap = (*trapezoidalBounds).parameters(); // el bueno aqui
597  hbedge = parameterTrap[0];
598  htedge = parameterTrap[1];
599  hapoth = parameterTrap[3];
600  uylfac = (htedge-hbedge)/(htedge+hbedge)/hapoth;
601  uxlden = 1 + yloc*uylfac;
602  }
603 
604  // Need to know position of trajectory in strip number for selecting the right APV later
605  if( TrajStrip==-1 ) {
606  int nstrips = Topo.nstrips();
607  float pitch = stripdet->surface().bounds().width() / nstrips;
608  TrajStrip = xloc/pitch + nstrips/2.0;
609  // Need additionnal corrections for endcap
610  if(TKlayers>=11) {
611  float TrajLocXMid = xloc / (1 + (htedge-hbedge)*yloc/(htedge+hbedge)/hapoth) ; // radialy extrapolated x loc position at middle
612  TrajStrip = TrajLocXMid/pitch + nstrips/2.0;
613  }
614  //cout<<" Layer "<<TKlayers<<" TrajStrip: "<<nstrips<<" "<<pitch<<" "<<TrajStrip<<endl;
615  }
616 
617 
618  for(edmNew::DetSet<SiStripCluster>::const_iterator iter=DSViter->begin();iter!=DSViter->end();++iter) {
619  //iter is a single SiStripCluster
620  StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
621  float res = (parameters.first.x() - xloc);
622  float sigma = checkConsistency(parameters , xloc, xErr);
623  // The consistency is probably more accurately measured with the Chi2MeasurementEstimator. To use it
624  // you need a TransientTrackingRecHit instead of the cluster
625  //theEstimator= new Chi2MeasurementEstimator(30);
626  //const Chi2MeasurementEstimator *theEstimator(100);
627  //theEstimator->estimate(TM->tsos(), TransientTrackingRecHit);
628 
629  if(TKlayers>=11) {
630  res = parameters.first.x() - xloc/uxlden ; // radialy extrapolated x loc position at middle
631  sigma = abs(res) / sqrt(parameters.second.xx() + xErr*xErr/uxlden/uxlden + yErr*yErr*xloc*xloc*uylfac*uylfac/uxlden/uxlden/uxlden/uxlden);
632  }
633 
634  SiStripClusterInfo clusterInfo = SiStripClusterInfo(*iter, es, ClusterId);
635  // signal to noise from SiStripClusterInfo not working in 225. I'll fix this after the interface
636  // redesign in 300 -ku
637  //float cluster_info[7] = {res, sigma, parameters.first.x(), sqrt(parameters.second.xx()), parameters.first.y(), sqrt(parameters.second.yy()), signal_to_noise};
638  std::vector< float > cluster_info;
639  cluster_info.push_back(res);
640  cluster_info.push_back(sigma);
641  cluster_info.push_back(parameters.first.x());
642  cluster_info.push_back(sqrt(parameters.second.xx()));
643  cluster_info.push_back(parameters.first.y());
644  cluster_info.push_back(sqrt(parameters.second.yy()));
645  cluster_info.push_back( clusterInfo.signalOverNoise() );
646  //cluster_info.push_back( clusterInfo.getSignalOverNoise() );
647  VCluster_info.push_back(cluster_info);
648  nClusters++;
649  LogDebug("SiStripHitEfficiency:HitEff") << "Have ID match. residual = " << VCluster_info.back()[0] << " res sigma = " << VCluster_info.back()[1] << endl;
650  LogDebug("SiStripHitEfficiency:HitEff") << "trajectory measurement compatability estimate = " << (*itm).estimate() << endl;
651  LogDebug("SiStripHitEfficiency:HitEff") << "hit position = " << parameters.first.x() << " hit error = " << sqrt(parameters.second.xx()) << " trajectory position = " << xloc << " traj error = " << xErr << endl;
652  }
653  }
654  }
655  float FinalResSig = 1000.0;
656  float FinalCluster[7]= {1000.0, 1000.0, 0.0, 0.0, 0.0, 0.0, 0.0};
657  if (nClusters > 0) {
658  LogDebug("SiStripHitEfficiency:HitEff") << "found clusters > 0" << endl;
659  if (nClusters > 1) {
660  //get the smallest one
661  vector< vector<float> >::iterator ires;
662  for (ires=VCluster_info.begin(); ires!=VCluster_info.end(); ires++){
663  if ( abs((*ires)[1]) < abs(FinalResSig)) {
664  FinalResSig = (*ires)[1];
665  for (unsigned int i = 0; i<ires->size(); i++) {
666  LogDebug("SiStripHitEfficiency:HitEff") << "filling final cluster. i = " << i << " before fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
667  FinalCluster[i] = (*ires)[i];
668  LogDebug("SiStripHitEfficiency:HitEff") << "filling final cluster. i = " << i << " after fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
669  }
670  }
671  LogDebug("SiStripHitEfficiency:HitEff") << "iresidual = " << (*ires)[0] << " isigma = " << (*ires)[1] << " and FinalRes = " << FinalCluster[0] << endl;
672  }
673  }
674  else {
675  FinalResSig = VCluster_info.at(0)[1];
676  for (unsigned int i = 0; i<VCluster_info.at(0).size(); i++) {
677  FinalCluster[i] = VCluster_info.at(0)[i];
678  }
679  }
680  nClusters=0;
681  VCluster_info.clear();
682  }
683 
684  LogDebug("SiStripHitEfficiency:HitEff") << "Final residual in X = " << FinalCluster[0] << "+-" << (FinalCluster[0]/FinalResSig) << endl;
685  LogDebug("SiStripHitEfficiency:HitEff") << "Checking location of trajectory: abs(yloc) = " << abs(yloc) << " abs(xloc) = " << abs(xloc) << endl;
686  LogDebug("SiStripHitEfficiency:HitEff") << "Checking location of cluster hit: yloc = " << FinalCluster[4] << "+-" << FinalCluster[5] << " xloc = " << FinalCluster[2] << "+-" << FinalCluster[3] << endl;
687  LogDebug("SiStripHitEfficiency:HitEff") << "Final cluster signal to noise = " << FinalCluster[6] << endl;
688 
689  float exclusionWidth = 0.4;
690  float TOBexclusion = 0.0;
691  float TECexRing5 = -0.89;
692  float TECexRing6 = -0.56;
693  float TECexRing7 = 0.60;
694  //Added by Chris Edelmaier to do TEC bonding exclusion
695  int subdetector = ((iidd>>25) & 0x7);
696  int ringnumber = ((iidd>>5) & 0x7);
697 
698  //New TOB and TEC bonding region exclusion zone
699  if((TKlayers >= 5 && TKlayers < 11)||((subdetector == 6)&&( (ringnumber >= 5)&&(ringnumber <=7) ))) {
700  //There are only 2 cases that we need to exclude for
701  float highzone = 0.0;
702  float lowzone = 0.0;
703  float higherr = yloc + 5.0*yErr;
704  float lowerr = yloc - 5.0*yErr;
705  if(TKlayers >= 5 && TKlayers < 11) {
706  //TOB zone
707  highzone = TOBexclusion + exclusionWidth;
708  lowzone = TOBexclusion - exclusionWidth;
709  } else if (ringnumber == 5) {
710  //TEC ring 5
711  highzone = TECexRing5 + exclusionWidth;
712  lowzone = TECexRing5 - exclusionWidth;
713  } else if (ringnumber == 6) {
714  //TEC ring 6
715  highzone = TECexRing6 + exclusionWidth;
716  lowzone = TECexRing6 - exclusionWidth;
717  } else if (ringnumber == 7) {
718  //TEC ring 7
719  highzone = TECexRing7 + exclusionWidth;
720  lowzone = TECexRing7 - exclusionWidth;
721  }
722  //Now that we have our exclusion region, we just have to properly identify it
723  if((highzone <= higherr)&&(highzone >= lowerr)) withinAcceptance = false;
724  if((lowzone >= lowerr)&&(lowzone <= higherr)) withinAcceptance = false;
725  if((higherr <= highzone)&&(higherr >= lowzone)) withinAcceptance = false;
726  if((lowerr >= lowzone)&&(lowerr <= highzone)) withinAcceptance = false;
727  }
728 
729  // fill ntuple varibles
730  //get global position from module id number iidd
731  TrajGlbX = xglob;
732  TrajGlbY = yglob;
733  TrajGlbZ = zglob;
734 
735  TrajLocErrX = xErr;
736  TrajLocErrY = yErr;
737 
738  Id = iidd;
739  run = run_nr;
740  event = ev_nr;
741  bunchx = bunch_nr;
742  //if ( SiStripQuality_->IsModuleBad(iidd) ) {
743  if ( SiStripQuality_->getBadApvs(iidd)!=0 ) {
744  SiStripQualBad = 1;
745  LogDebug("SiStripHitEfficiency:HitEff") << "strip is bad from SiStripQuality" << endl;
746  } else {
747  SiStripQualBad = 0;
748  LogDebug("SiStripHitEfficiency:HitEff") << "strip is good from SiStripQuality" << endl;
749  }
750 
751  //check for FED-detected errors and include those in SiStripQualBad
752  for (unsigned int ii=0;ii< fedErrorIds->size();ii++) {
753  if (iidd == (*fedErrorIds)[ii].rawId() )
754  SiStripQualBad = 1;
755  }
756 
757  TrajLocX = xloc;
758  TrajLocY = yloc;
759  TrajLocAngleX = angleX;
760  TrajLocAngleY = angleY;
761  ResX = FinalCluster[0];
762  ResXSig = FinalResSig;
763  if (FinalResSig != FinalCluster[1]) LogDebug("SiStripHitEfficiency:HitEff") << "Problem with best cluster selection because FinalResSig = " << FinalResSig << " and FinalCluster[1] = " << FinalCluster[1] << endl;
764  ClusterLocX = FinalCluster[2];
765  ClusterLocY = FinalCluster[4];
766  ClusterLocErrX = FinalCluster[3];
767  ClusterLocErrY = FinalCluster[5];
768  ClusterStoN = FinalCluster[6];
769 
770  // CM of APV crossed by traj
771  if(addCommonMode_)
772  if(commonModeDigis.isValid() && TrajStrip>=0 && TrajStrip<=768) {
773  edm::DetSetVector<SiStripRawDigi>::const_iterator digiframe = commonModeDigis->find(iidd);
774  if(digiframe != commonModeDigis->end())
775  if( (unsigned) TrajStrip/128 < digiframe->data.size())
776  commonMode = digiframe->data.at(TrajStrip/128).adc();
777  }
778 
779  LogDebug("SiStripHitEfficiency:HitEff") << "before check good" << endl;
780 
781  if ( FinalResSig < 999.0) { //could make requirement on track/hit consistency, but for
782  //now take anything with a hit on the module
783  LogDebug("SiStripHitEfficiency:HitEff") << "hit being counted as good " << FinalCluster[0] << " FinalRecHit " <<
784  iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd <<
785  " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
786  ModIsBad = 0;
787  traj->Fill();
788  }
789  else {
790  LogDebug("SiStripHitEfficiency:HitEff") << "hit being counted as bad ######### Invalid RPhi FinalResX " << FinalCluster[0] << " FinalRecHit " <<
791  iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd <<
792  " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
793  ModIsBad = 1;
794  traj->Fill();
795 
796  LogDebug("SiStripHitEfficiency:HitEff") << " RPhi Error " << sqrt(xErr*xErr + yErr*yErr) << " ErrorX " << xErr << " yErr " << yErr << endl;
797  } LogDebug("SiStripHitEfficiency:HitEff") << "after good location check" << endl;
798  } LogDebug("SiStripHitEfficiency:HitEff") << "after list of clusters" << endl;
799  } LogDebug("SiStripHitEfficiency:HitEff") << "After layers=TKLayers if" << endl;
800  } LogDebug("SiStripHitEfficiency:HitEff") << "After looping over TrajAtValidHit list" << endl;
801  } LogDebug("SiStripHitEfficiency:HitEff") << "end TMeasurement loop" << endl;
802  } LogDebug("SiStripHitEfficiency:HitEff") << "end of trajectories loop" << endl;
803  }
804 }
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
size_type size() const
Definition: EDCollection.h:97
EventNumber_t event() const
Definition: EventID.h:41
std::vector< TrajectoryMeasurement > measurements(const DetLayer &layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
float TrajGlbZ
Definition: HitEff.h:102
bool check2DPartner(unsigned int iidd, const std::vector< TrajectoryMeasurement > &traj)
Definition: HitEff.cc:853
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
float ClusterLocErrX
Definition: HitEff.h:105
const_iterator end(bool update=false) const
float ClusterLocErrY
Definition: HitEff.h:105
ConstRecHitPointer const & recHit() const
float TrajGlbY
Definition: HitEff.h:102
float TrajLocX
Definition: HitEff.h:103
virtual const std::array< const float, 4 > parameters() const
const_iterator end() const
last iterator over the map (read only)
float ResXSig
Definition: HitEff.h:106
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
bool useAllHitsFromTracksWithMissingHits_
Definition: HitEff.h:73
TString subdetector
int ires[2]
float instLumi
Definition: HitEff.h:113
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:15
int bunchCrossing() const
Definition: EventBase.h:66
const Bounds & bounds() const
Definition: Surface.h:120
std::pair< LocalPoint, LocalError > LocalValues
float ClusterStoN
Definition: HitEff.h:105
data_type const * const_iterator
Definition: DetSetNew.h:30
float TrajLocErrX
Definition: HitEff.h:104
key_type key() const
Accessor for product key.
Definition: Ref.h:265
float commonMode
Definition: HitEff.h:114
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
Definition: Electron.h:6
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
const edm::EDGetTokenT< MeasurementTrackerEvent > trackerEvent_token_
Definition: HitEff.h:80
unsigned int bunchx
Definition: HitEff.h:111
static std::string const input
Definition: EdmProvDump.cc:44
virtual float width() const =0
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
bool isDoubleSided(unsigned int iidd, const TrackerTopology *tTopo) const
Definition: HitEff.cc:821
bool highPurity
Definition: HitEff.h:108
short getBadApvs(const uint32_t &detid) const
float PU
Definition: HitEff.h:113
int nHits
Definition: HitEff.h:109
unsigned int Id
Definition: HitEff.h:107
T sqrt(T t)
Definition: SSEVec.h:18
float timeAtIpInOutErr
Definition: MuonTime.h:15
const edm::EDGetTokenT< DetIdCollection > digis_token_
Definition: HitEff.h:79
int nDof
number of muon stations used
Definition: MuonTime.h:10
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int trajHitValid
Definition: HitEff.h:111
unsigned int ModIsBad
Definition: HitEff.h:107
float TrajGlbX
Definition: HitEff.h:102
const edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackAsso_token_
Definition: HitEff.h:77
bool isValid() const
Definition: HandleBase.h:74
bool withinAcceptance
Definition: HitEff.h:107
unsigned int run
Definition: HitEff.h:111
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
size_type size() const
Return the number of contained DetSets.
Definition: DetSetVector.h:283
float ClusterLocX
Definition: HitEff.h:105
ii
Definition: cuy.py:589
int events
Definition: HitEff.h:85
float TrajLocAngleX
Definition: HitEff.h:103
Definition: DetId.h:18
double checkConsistency(const StripClusterParameterEstimator::LocalValues &parameters, double xx, double xerr)
Definition: HitEff.cc:814
T const * product() const
Definition: Handle.h:81
bool useLastMeas_
Definition: HitEff.h:72
unsigned int whatlayer
Definition: HitEff.h:90
virtual int nstrips() const =0
unsigned int trackMultiplicityCut_
Definition: HitEff.h:70
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:21
unsigned int checkLayer(unsigned int iidd, const TrackerTopology *tTopo)
Definition: HitEff.cc:869
unsigned int layers
Definition: HitEff.h:88
float pT
Definition: HitEff.h:110
float ResX
Definition: HitEff.h:106
const edm::EDGetTokenT< reco::TrackCollection > combinatorialTracks_token_
Definition: HitEff.h:75
bool addLumi_
Definition: HitEff.h:67
unsigned int SiStripQualBad
Definition: HitEff.h:107
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
bool addCommonMode_
Definition: HitEff.h:68
size_type size() const
edm::EventID id() const
Definition: EventBase.h:60
fixed size matrix
HLT enums.
float TrajLocErrY
Definition: HitEff.h:104
float TrajLocAngleY
Definition: HitEff.h:103
bool cutOnTracks_
Definition: HitEff.h:69
const_iterator begin() const
first iterator over the map (read only)
bool useFirstMeas_
Definition: HitEff.h:71
float ClusterLocY
Definition: HitEff.h:105
const edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clusters_token_
Definition: HitEff.h:78
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
const edm::EDGetTokenT< LumiScalersCollection > scalerToken_
Definition: HitEff.h:64
TTree * traj
Definition: HitEff.h:84
const edm::EDGetTokenT< std::vector< Trajectory > > trajectories_token_
Definition: HitEff.h:76
T const * product() const
Definition: ESHandle.h:86
const edm::EDGetTokenT< edm::DetSetVector< SiStripRawDigi > > commonModeToken_
Definition: HitEff.h:65
float timeAtIpInOut
Definition: MuonTime.h:14
int EventTrackCKF
Definition: HitEff.h:85
const_iterator begin(bool update=false) const
float TrajLocY
Definition: HitEff.h:103
unsigned int tecSide(const DetId &id) const
void HitEff::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 100 of file HitEff.cc.

References addCommonMode_, addLumi_, bunchx, vertices_cff::chi2, ClusterLocErrX, ClusterLocErrY, ClusterLocX, ClusterLocY, ClusterStoN, commonMode, compSettings, events, EventTrackCKF, TFileService::file(), highPurity, Id, instLumi, TFileService::make(), ModIsBad, nHits, AlCaHLTBitMon_ParallelJobs::p, pT, PU, ResX, ResXSig, run, SiStripQualBad, tquality, traj, TrajGlbX, TrajGlbY, TrajGlbZ, trajHitValid, TrajLocAngleX, TrajLocAngleY, TrajLocErrX, TrajLocErrY, TrajLocX, TrajLocY, whatlayer, and withinAcceptance.

100  {
101 
103  if(compSettings>0){
104  edm::LogInfo("SiStripHitEfficiency:HitEff")<<"the compressions settings are:"<< compSettings << std::endl;
105  fs->file().SetCompressionSettings(compSettings);
106  }
107 
108  traj = fs->make<TTree>("traj","tree of trajectory positions");
109  #ifdef ExtendedCALIBTree
110  traj->Branch("timeDT",&timeDT,"timeDT/F");
111  traj->Branch("timeDTErr",&timeDTErr,"timeDTErr/F");
112  traj->Branch("timeDTDOF",&timeDTDOF,"timeDTDOF/I");
113  traj->Branch("timeECAL",&timeECAL,"timeECAL/F");
114  traj->Branch("dedx",&dedx,"dedx/F");
115  traj->Branch("dedxNOM",&dedxNOM,"dedxNOM/I");
116  traj->Branch("nLostHits",&nLostHits,"nLostHits/I");
117  traj->Branch("chi2",&chi2,"chi2/F");
118  traj->Branch("p",&p,"p/F");
119  #endif
120  traj->Branch("TrajGlbX",&TrajGlbX,"TrajGlbX/F");
121  traj->Branch("TrajGlbY",&TrajGlbY,"TrajGlbY/F");
122  traj->Branch("TrajGlbZ",&TrajGlbZ,"TrajGlbZ/F");
123  traj->Branch("TrajLocX",&TrajLocX,"TrajLocX/F");
124  traj->Branch("TrajLocY",&TrajLocY,"TrajLocY/F");
125  traj->Branch("TrajLocAngleX",&TrajLocAngleX,"TrajLocAngleX/F");
126  traj->Branch("TrajLocAngleY",&TrajLocAngleY,"TrajLocAngleY/F");
127  traj->Branch("TrajLocErrX",&TrajLocErrX,"TrajLocErrX/F");
128  traj->Branch("TrajLocErrY",&TrajLocErrY,"TrajLocErrY/F");
129  traj->Branch("ClusterLocX",&ClusterLocX,"ClusterLocX/F");
130  traj->Branch("ClusterLocY",&ClusterLocY,"ClusterLocY/F");
131  traj->Branch("ClusterLocErrX",&ClusterLocErrX,"ClusterLocErrX/F");
132  traj->Branch("ClusterLocErrY",&ClusterLocErrY,"ClusterLocErrY/F");
133  traj->Branch("ClusterStoN",&ClusterStoN,"ClusterStoN/F");
134  traj->Branch("ResX",&ResX,"ResX/F");
135  traj->Branch("ResXSig",&ResXSig,"ResXSig/F");
136  traj->Branch("ModIsBad",&ModIsBad,"ModIsBad/i");
137  traj->Branch("SiStripQualBad",&SiStripQualBad,"SiStripQualBad/i");
138  traj->Branch("withinAcceptance",&withinAcceptance,"withinAcceptance/O");
139  traj->Branch("nHits",&nHits,"nHits/I");
140  traj->Branch("pT",&pT,"pT/F");
141  traj->Branch("highPurity",&highPurity,"highPurity/O");
142  traj->Branch("trajHitValid", &trajHitValid, "trajHitValid/i");
143  traj->Branch("Id",&Id,"Id/i");
144  traj->Branch("run",&run,"run/i");
145  traj->Branch("event",&event,"event/i");
146  traj->Branch("layer",&whatlayer,"layer/i");
147  traj->Branch("tquality",&tquality,"tquality/I");
148  traj->Branch("bunchx",&bunchx,"bunchx/I");
149  if(addLumi_) {
150  traj->Branch("instLumi",&instLumi,"instLumi/F");
151  traj->Branch("PU",&PU,"PU/F");
152  }
153  if(addCommonMode_) traj->Branch("commonMode",&commonMode,"commonMode/F");
154 
155  events = 0;
156  EventTrackCKF = 0;
157 
158 }
float TrajGlbZ
Definition: HitEff.h:102
float ClusterLocErrX
Definition: HitEff.h:105
float ClusterLocErrY
Definition: HitEff.h:105
float TrajGlbY
Definition: HitEff.h:102
float TrajLocX
Definition: HitEff.h:103
float ResXSig
Definition: HitEff.h:106
int compSettings
Definition: HitEff.h:87
float instLumi
Definition: HitEff.h:113
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
float ClusterStoN
Definition: HitEff.h:105
float TrajLocErrX
Definition: HitEff.h:104
float commonMode
Definition: HitEff.h:114
unsigned int bunchx
Definition: HitEff.h:111
bool highPurity
Definition: HitEff.h:108
float PU
Definition: HitEff.h:113
int nHits
Definition: HitEff.h:109
unsigned int Id
Definition: HitEff.h:107
unsigned int trajHitValid
Definition: HitEff.h:111
unsigned int ModIsBad
Definition: HitEff.h:107
float TrajGlbX
Definition: HitEff.h:102
bool withinAcceptance
Definition: HitEff.h:107
unsigned int run
Definition: HitEff.h:111
float ClusterLocX
Definition: HitEff.h:105
int events
Definition: HitEff.h:85
float TrajLocAngleX
Definition: HitEff.h:103
TFile & file() const
return opened TFile
Definition: TFileService.h:37
unsigned int whatlayer
Definition: HitEff.h:90
float pT
Definition: HitEff.h:110
float ResX
Definition: HitEff.h:106
bool addLumi_
Definition: HitEff.h:67
unsigned int SiStripQualBad
Definition: HitEff.h:107
bool addCommonMode_
Definition: HitEff.h:68
float TrajLocErrY
Definition: HitEff.h:104
float TrajLocAngleY
Definition: HitEff.h:103
float ClusterLocY
Definition: HitEff.h:105
TTree * traj
Definition: HitEff.h:84
int EventTrackCKF
Definition: HitEff.h:85
float TrajLocY
Definition: HitEff.h:103
int tquality
Definition: HitEff.h:112
Definition: event.py:1
bool HitEff::check2DPartner ( unsigned int  iidd,
const std::vector< TrajectoryMeasurement > &  traj 
)

Definition at line 853 of file HitEff.cc.

Referenced by analyze().

853  {
854  unsigned int partner_iidd = 0;
855  bool found2DPartner = false;
856  // first get the id of the other detector
857  if ((iidd & 0x3)==1) partner_iidd = iidd+1;
858  if ((iidd & 0x3)==2) partner_iidd = iidd-1;
859  // next look in the trajectory measurements for a measurement from that detector
860  // loop through trajectory measurements to find the partner_iidd
861  for (std::vector<TrajectoryMeasurement>::const_iterator iTM=traj.begin(); iTM!=traj.end(); ++iTM) {
862  if (iTM->recHit()->geographicalId().rawId()==partner_iidd) {
863  found2DPartner = true;
864  }
865  }
866  return found2DPartner;
867 }
TTree * traj
Definition: HitEff.h:84
double HitEff::checkConsistency ( const StripClusterParameterEstimator::LocalValues parameters,
double  xx,
double  xerr 
)

Definition at line 814 of file HitEff.cc.

References funct::abs(), relativeConstraints::error, mathSSE::sqrt(), and geometryCSVtoXML::xx.

Referenced by analyze().

814  {
815  double error = sqrt(parameters.second.xx() + xerr*xerr);
816  double separation = abs(parameters.first.x() - xx);
817  double consistency = separation/error;
818  return consistency;
819 }
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int HitEff::checkLayer ( unsigned int  iidd,
const TrackerTopology tTopo 
)

Definition at line 869 of file HitEff.cc.

References DEFINE_FWK_MODULE, digitizers_cfi::strip, DetId::subdetId(), StripSubdetector::TEC, TrackerTopology::tecWheel(), StripSubdetector::TIB, TrackerTopology::tibLayer(), StripSubdetector::TID, TrackerTopology::tidWheel(), StripSubdetector::TOB, and TrackerTopology::tobLayer().

Referenced by analyze().

869  {
871  unsigned int subid=strip.subdetId();
872  if (subid == StripSubdetector::TIB) {
873 
874  return tTopo->tibLayer(iidd);
875  }
876  if (subid == StripSubdetector::TOB) {
877 
878  return tTopo->tobLayer(iidd) + 4 ;
879  }
880  if (subid == StripSubdetector::TID) {
881 
882  return tTopo->tidWheel(iidd) + 10;
883  }
884  if (subid == StripSubdetector::TEC) {
885 
886  return tTopo->tecWheel(iidd) + 13 ;
887  }
888  return 0;
889 }
unsigned int tibLayer(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
unsigned int tecWheel(const DetId &id) const
unsigned int tobLayer(const DetId &id) const
void HitEff::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 806 of file HitEff.cc.

References events, EventTrackCKF, LogDebug, and traj.

Referenced by o2olib.O2ORunMgr::executeJob().

806  {
807  traj->GetDirectory()->cd();
808  traj->Write();
809 
810  LogDebug("SiStripHitEfficiency:HitEff") << " Events Analysed " << events << endl;
811  LogDebug("SiStripHitEfficiency:HitEff") << " Number Of Tracked events " << EventTrackCKF << endl;
812 }
#define LogDebug(id)
int events
Definition: HitEff.h:85
TTree * traj
Definition: HitEff.h:84
int EventTrackCKF
Definition: HitEff.h:85
bool HitEff::isDoubleSided ( unsigned int  iidd,
const TrackerTopology tTopo 
) const

Definition at line 821 of file HitEff.cc.

References digitizers_cfi::strip, DetId::subdetId(), StripSubdetector::TEC, TrackerTopology::tecRing(), StripSubdetector::TIB, TrackerTopology::tibLayer(), StripSubdetector::TID, TrackerTopology::tidRing(), StripSubdetector::TOB, and TrackerTopology::tobLayer().

Referenced by analyze().

821  {
823  unsigned int subid=strip.subdetId();
824  unsigned int layer = 0;
825  if (subid == StripSubdetector::TIB) {
826 
827  layer = tTopo->tibLayer(iidd);
828  if (layer == 1 || layer == 2) return true;
829  else return false;
830  }
831  else if (subid == StripSubdetector::TOB) {
832 
833  layer = tTopo->tobLayer(iidd) + 4 ;
834  if (layer == 5 || layer == 6) return true;
835  else return false;
836  }
837  else if (subid == StripSubdetector::TID) {
838 
839  layer = tTopo->tidRing(iidd) + 10;
840  if (layer == 11 || layer == 12) return true;
841  else return false;
842  }
843  else if (subid == StripSubdetector::TEC) {
844 
845  layer = tTopo->tecRing(iidd) + 13 ;
846  if (layer == 14 || layer == 15 || layer == 18) return true;
847  else return false;
848  }
849  else
850  return false;
851 }
unsigned int tibLayer(const DetId &id) const
unsigned int tidRing(const DetId &id) const
unsigned int tecRing(const DetId &id) const
ring id
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
unsigned int tobLayer(const DetId &id) const

Member Data Documentation

bool HitEff::addCommonMode_
private

Definition at line 68 of file HitEff.h.

Referenced by analyze(), beginJob(), and HitEff().

bool HitEff::addLumi_
private

Definition at line 67 of file HitEff.h.

Referenced by analyze(), beginJob(), and HitEff().

unsigned int HitEff::bunchx
private

Definition at line 111 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocErrX
private

Definition at line 105 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocErrY
private

Definition at line 105 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocX
private

Definition at line 105 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocY
private

Definition at line 105 of file HitEff.h.

Referenced by analyze(), and beginJob().

const edm::EDGetTokenT< edmNew::DetSetVector<SiStripCluster> > HitEff::clusters_token_
private

Definition at line 78 of file HitEff.h.

Referenced by analyze().

float HitEff::ClusterStoN
private

Definition at line 105 of file HitEff.h.

Referenced by analyze(), and beginJob().

const edm::EDGetTokenT< reco::TrackCollection > HitEff::combinatorialTracks_token_
private

Definition at line 75 of file HitEff.h.

Referenced by analyze().

float HitEff::commonMode
private

Definition at line 114 of file HitEff.h.

Referenced by analyze(), and beginJob().

const edm::EDGetTokenT<edm::DetSetVector<SiStripRawDigi> > HitEff::commonModeToken_
private

Definition at line 65 of file HitEff.h.

Referenced by analyze().

int HitEff::compSettings
private

Definition at line 87 of file HitEff.h.

Referenced by beginJob(), and HitEff().

edm::ParameterSet HitEff::conf_
private

Definition at line 82 of file HitEff.h.

Referenced by HitEff().

bool HitEff::cutOnTracks_
private

Definition at line 69 of file HitEff.h.

Referenced by analyze(), and HitEff().

bool HitEff::DEBUG
private

Definition at line 89 of file HitEff.h.

Referenced by HitEff().

const edm::EDGetTokenT<DetIdCollection> HitEff::digis_token_
private

Definition at line 79 of file HitEff.h.

Referenced by analyze().

unsigned int HitEff::event
private
int HitEff::events
private
int HitEff::EventTrackCKF
private

Definition at line 85 of file HitEff.h.

Referenced by analyze(), beginJob(), and endJob().

bool HitEff::highPurity
private

Definition at line 108 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::Id
private

Definition at line 107 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::instLumi
private

Definition at line 113 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::layers
private

Definition at line 88 of file HitEff.h.

Referenced by analyze(), and HitEff().

unsigned int HitEff::ModIsBad
private

Definition at line 107 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::nHits
private

Definition at line 109 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::pT
private

Definition at line 110 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::PU
private

Definition at line 113 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ResX
private

Definition at line 106 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ResXSig
private

Definition at line 106 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::run
private
const edm::EDGetTokenT<LumiScalersCollection> HitEff::scalerToken_
private

Definition at line 64 of file HitEff.h.

Referenced by analyze().

unsigned int HitEff::SiStripQualBad
private

Definition at line 107 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::tquality
private

Definition at line 112 of file HitEff.h.

Referenced by beginJob().

const edm::EDGetTokenT<MeasurementTrackerEvent> HitEff::trackerEvent_token_
private

Definition at line 80 of file HitEff.h.

Referenced by analyze().

unsigned int HitEff::trackMultiplicityCut_
private

Definition at line 70 of file HitEff.h.

Referenced by analyze(), and HitEff().

TTree* HitEff::traj
private

Definition at line 84 of file HitEff.h.

Referenced by analyze(), beginJob(), and endJob().

const edm::EDGetTokenT< std::vector<Trajectory> > HitEff::trajectories_token_
private

Definition at line 76 of file HitEff.h.

Referenced by analyze().

float HitEff::TrajGlbX
private

Definition at line 102 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajGlbY
private

Definition at line 102 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajGlbZ
private

Definition at line 102 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::trajHitValid
private

Definition at line 111 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocAngleX
private

Definition at line 103 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocAngleY
private

Definition at line 103 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocErrX
private

Definition at line 104 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocErrY
private

Definition at line 104 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocX
private

Definition at line 103 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocY
private

Definition at line 103 of file HitEff.h.

Referenced by analyze(), and beginJob().

const edm::EDGetTokenT< TrajTrackAssociationCollection > HitEff::trajTrackAsso_token_
private

Definition at line 77 of file HitEff.h.

Referenced by analyze().

bool HitEff::useAllHitsFromTracksWithMissingHits_
private

Definition at line 73 of file HitEff.h.

Referenced by analyze(), and HitEff().

bool HitEff::useFirstMeas_
private

Definition at line 71 of file HitEff.h.

Referenced by analyze(), and HitEff().

bool HitEff::useLastMeas_
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and HitEff().

unsigned int HitEff::whatlayer
private

Definition at line 90 of file HitEff.h.

Referenced by analyze(), and beginJob().

bool HitEff::withinAcceptance
private

Definition at line 107 of file HitEff.h.

Referenced by analyze(), and beginJob().