CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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
 
virtual ~HitEff ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

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

Private Attributes

unsigned int bunchx
 
float chi2
 
float ClusterLocErrX
 
float ClusterLocErrY
 
float ClusterLocX
 
float ClusterLocY
 
float ClusterStoN
 
edm::ParameterSet conf_
 
bool DEBUG
 
float dedx
 
int dedxNOM
 
unsigned int event
 
int events
 
int EventTrackCKF
 
unsigned int Id
 
int istep
 
unsigned int layers
 
unsigned int ModIsBad
 
int nHits
 
int nLostHits
 
float p
 
float pT
 
float ResX
 
float ResXSig
 
unsigned int run
 
unsigned int SiStripQualBad
 
float timeDT
 
int timeDTDOF
 
float timeDTErr
 
float timeECAL
 
int tquality
 
TTree * traj
 
float TrajGlbX
 
float TrajGlbY
 
float TrajGlbZ
 
unsigned int trajHitValid
 
float TrajLocAngleX
 
float TrajLocAngleY
 
float TrajLocErrX
 
float TrajLocErrY
 
float TrajLocX
 
float TrajLocY
 
unsigned int whatlayer
 
bool withinAcceptance
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- 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 44 of file HitEff.h.

Constructor & Destructor Documentation

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

Definition at line 73 of file HitEff.cc.

References conf_, DEBUG, edm::ParameterSet::getParameter(), and layers.

73  :
74  conf_(conf)
75 {
76  layers =conf_.getParameter<int>("Layer");
77  DEBUG = conf_.getParameter<bool>("Debug");
78 }
T getParameter(std::string const &) const
bool DEBUG
Definition: HitEff.h:66
unsigned int layers
Definition: HitEff.h:65
edm::ParameterSet conf_
Definition: HitEff.h:60
HitEff::~HitEff ( )
virtual

Definition at line 81 of file HitEff.cc.

81 { }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 133 of file HitEff.cc.

References funct::abs(), AnalyticalPropagator_cfi::AnalyticalPropagator, anyDirection, edm::EventBase::bunchCrossing(), bunchx, check2DPartner(), checkConsistency(), checkLayer(), chi2, ClusterLocErrX, ClusterLocErrY, ClusterLocX, ClusterLocY, ClusterStoN, conf_, gather_cfg::cout, DEBUG, HLT_25ns14e33_v1_cff::estimator, edm::EventID::event(), events, EventTrackCKF, edm::EventSetup::get(), edm::Event::getByLabel(), edm::ParameterSet::getParameter(), i, edm::EventBase::id(), Id, cuy::ii, input, ires, isDoubleSided(), getDQMSummary::iter, layers, StripClusterParameterEstimator::localParameters(), LayerMeasurements::measurements(), ModIsBad, reco::MuonTime::nDof, nHits, nLostHits, p, Parameters::parameters, edm::Handle< T >::product(), edm::ESHandle< class >::product(), HLT_25ns14e33_v1_cff::propagator, pT, TrajectoryMeasurement::recHit(), dt_dqm_sourceclient_common_cff::reco, ResX, ResXSig, edm::EventID::run(), run, SiStripClusterInfo::signalOverNoise(), SiStripQualBad, mathSSE::sqrt(), StripSubdetector::TEC, TrackerTopology::tecSide(), reco::MuonTime::timeAtIpInOut, reco::MuonTime::timeAtIpInOutErr, timeDT, timeDTDOF, timeDTErr, timeECAL, tmp, patCandidatesForDimuonsSequences_cff::tracker, traj, TrajGlbX, TrajGlbY, TrajGlbZ, trajHitValid, TrajLocAngleX, TrajLocAngleY, TrajLocErrX, TrajLocErrY, TrajLocX, TrajLocY, whatlayer, and withinAcceptance.

133  {
134  //Retrieve tracker topology from geometry
135  edm::ESHandle<TrackerTopology> tTopoHandle;
136  es.get<IdealGeometryRecord>().get(tTopoHandle);
137  const TrackerTopology* const tTopo = tTopoHandle.product();
138 
139  // bool DEBUG = false;
140 
141  if (DEBUG) cout << "beginning analyze from HitEff" << endl;
142 
143  using namespace edm;
144  using namespace reco;
145  // Step A: Get Inputs
146 
147  int run_nr = e.id().run();
148  int ev_nr = e.id().event();
149  int bunch_nr = e.bunchCrossing();
150 
151  //CombinatoriaTrack
152  edm::Handle<reco::TrackCollection> trackCollectionCKF;
153  edm::InputTag TkTagCKF = conf_.getParameter<edm::InputTag>("combinatorialTracks");
154  e.getByLabel(TkTagCKF,trackCollectionCKF);
155 
156  edm::Handle<std::vector<Trajectory> > TrajectoryCollectionCKF;
157  edm::InputTag TkTrajCKF = conf_.getParameter<edm::InputTag>("trajectories");
158  e.getByLabel(TkTrajCKF,TrajectoryCollectionCKF);
159 
160  // Clusters
161  // get the SiStripClusters from the event
163  e.getByLabel("siStripClusters", theClusters);
164 
165  //get tracker geometry
167  es.get<TrackerDigiGeometryRecord>().get(tracker);
168  const TrackerGeometry * tkgeom=&(* tracker);
169 
170  //get Cluster Parameter Estimator
171  //std::string cpe = conf_.getParameter<std::string>("StripCPE");
173  es.get<TkStripCPERecord>().get("StripCPEfromTrackAngle", parameterestimator);
174  const StripClusterParameterEstimator &stripcpe(*parameterestimator);
175 
176  // get the SiStripQuality records
177  edm::ESHandle<SiStripQuality> SiStripQuality_;
178 //LQ commenting the try/catch that causes problem in 74X calibTree production
179 // try {
180 // es.get<SiStripQualityRcd>().get("forCluster",SiStripQuality_);
181 // }
182 // catch (...) {
183  es.get<SiStripQualityRcd>().get(SiStripQuality_);
184 // }
185 
186  edm::ESHandle<MagneticField> magFieldHandle;
187  es.get<IdealMagneticFieldRecord>().get(magFieldHandle);
188  const MagneticField* magField_ = magFieldHandle.product();
189 
190  // get the list of module IDs with FED-detected errors
191  edm::Handle< DetIdCollection > fedErrorIds;
192  e.getByLabel("siStripDigis", fedErrorIds );
193 
194  ESHandle<MeasurementTracker> measurementTrackerHandle;
195  es.get<CkfComponentsRecord>().get(measurementTrackerHandle);
196 
197  edm::Handle<MeasurementTrackerEvent> measurementTrackerEvent;
198  e.getByLabel("MeasurementTrackerEvent", measurementTrackerEvent);
199 
201  es.get<TrackingComponentsRecord>().get("Chi2",est);
202 
204  es.get<TrackingComponentsRecord>().get("PropagatorWithMaterial",prop);
205  const Propagator* thePropagator = prop.product();
206 
207  events++;
208 
209  // *************** SiStripCluster Collection
210  const edmNew::DetSetVector<SiStripCluster>& input = *theClusters;
211 
212  //go through clusters to write out global position of good clusters for the layer understudy for comparison
213  // Loop through clusters just to print out locations
214  // Commented out to avoid discussion, should really be deleted.
215  /*
216  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
217  // DSViter is a vector of SiStripClusters located on a single module
218  unsigned int ClusterId = DSViter->id();
219  DetId ClusterDetId(ClusterId);
220  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
221 
222  edmNew::DetSet<SiStripCluster>::const_iterator begin=DSViter->begin();
223  edmNew::DetSet<SiStripCluster>::const_iterator end =DSViter->end();
224  for(edmNew::DetSet<SiStripCluster>::const_iterator iter=begin;iter!=end;++iter) {
225  //iter is a single SiStripCluster
226  StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
227 
228  const Surface* surface;
229  surface = &(tracker->idToDet(ClusterDetId)->surface());
230  LocalPoint lp = parameters.first;
231  GlobalPoint gp = surface->toGlobal(lp);
232  unsigned int layer = checkLayer(ClusterId, tTopo);
233  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;
234  }
235  }
236  */
237 
238  // Tracking
239  const reco::TrackCollection *tracksCKF=trackCollectionCKF.product();
240  if (DEBUG) cout << "number ckf tracks found = " << tracksCKF->size() << endl;
241  //if (tracksCKF->size() == 1 ){
242  if (tracksCKF->size() > 0 && tracksCKF->size()<100) {
243  if (DEBUG) cout << "starting checking good event with < 100 tracks" << endl;
244 
245  EventTrackCKF++;
246 
247  /*
248 
249  //get dEdx info if available
250  Handle<ValueMap<DeDxData> > dEdxUncalibHandle;
251  if (e.getByLabel("dedxMedianCTF", dEdxUncalibHandle)) {
252  const ValueMap<DeDxData> dEdxTrackUncalib = *dEdxUncalibHandle.product();
253 
254  reco::TrackRef itTrack = reco::TrackRef( trackCollectionCKF, 0 );
255  dedx = dEdxTrackUncalib[itTrack].dEdx();
256  dedxNOM = dEdxTrackUncalib[itTrack].numberOfMeasurements();
257  } else {
258  dedx = -999.0; dedxNOM = -999;
259  }
260 
261  */
262 
263  //get muon and ecal timing info if available
265  if(e.getByLabel("muonsWitht0Correction",muH)){
266  const MuonCollection & muonsT0 = *muH.product();
267  if(muonsT0.size()!=0) {
268  MuonTime mt0 = muonsT0[0].time();
269  timeDT = mt0.timeAtIpInOut;
271  timeDTDOF = mt0.nDof;
272 
273  bool hasCaloEnergyInfo = muonsT0[0].isEnergyValid();
274  if (hasCaloEnergyInfo) timeECAL = muonsT0[0].calEnergy().ecal_time;
275  }
276  } else {
277  timeDT = -999.0; timeDTErr = -999.0; timeDTDOF = -999; timeECAL = -999.0;
278  }
279 
280  // actually should do a loop over all the tracks in the event here
281 
282  for (vector<Trajectory>::const_iterator itraj = TrajectoryCollectionCKF.product()->begin();
283  itraj != TrajectoryCollectionCKF.product()->end();
284  itraj++) {
285 
286  // for each track, fill some variables such as number of hits and momentum
287  nHits = itraj->foundHits();
288  nLostHits = itraj->lostHits();
289  chi2 = (itraj->chiSquared()/itraj->ndof());
290  pT = sqrt( ( itraj->lastMeasurement().updatedState().globalMomentum().x() *
291  itraj->lastMeasurement().updatedState().globalMomentum().x()) +
292  ( itraj->lastMeasurement().updatedState().globalMomentum().y() *
293  itraj->lastMeasurement().updatedState().globalMomentum().y()) );
294  p = itraj->lastMeasurement().updatedState().globalMomentum().mag();
295 
296  //Put in code to check track quality
297 
298 
299  std::vector<TrajectoryMeasurement> TMeas=itraj->measurements();
300  vector<TrajectoryMeasurement>::iterator itm;
301  double xloc = 0.;
302  double yloc = 0.;
303  double xErr = 0.;
304  double yErr = 0.;
305  double angleX = -999.;
306  double angleY = -999.;
307  double xglob,yglob,zglob;
308 
309  for (itm=TMeas.begin();itm!=TMeas.end();itm++){
310  auto theInHit = (*itm).recHit();
311 
312  if(DEBUG) cout << "theInHit is valid = " << theInHit->isValid() << endl;
313 
314  unsigned int iidd = theInHit->geographicalId().rawId();
315 
316  unsigned int TKlayers = checkLayer(iidd, tTopo);
317  if (DEBUG) cout << "TKlayer from trajectory: " << TKlayers << " from module = " << iidd << " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
318 
319  // If Trajectory measurement from TOB 6 or TEC 9, skip it because it's always valid they are filled later
320  if ( TKlayers == 10 || TKlayers == 22 ) {
321  if (DEBUG) cout << "skipping original TM for TOB 6 or TEC 9" << endl;
322  continue;
323  }
324 
325  // Make vector of TrajectoryAtInvalidHits to hold the trajectories
326  std::vector<TrajectoryAtInvalidHit> TMs;
327 
328  // Make AnalyticalPropagator to use in TAVH constructor
330 
331  // for double sided layers check both sensors--if no hit was found on either sensor surface,
332  // the trajectory measurements only have one invalid hit entry on the matched surface
333  // so get the TrajectoryAtInvalidHit for both surfaces and include them in the study
334  if (isDoubleSided(iidd, tTopo) && ((iidd & 0x3)==0) ) {
335  // do hit eff check twice--once for each sensor
336  //add a TM for each surface
337  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 1));
338  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 2));
339  } else if ( isDoubleSided(iidd, tTopo) && (!check2DPartner(iidd, TMeas)) ) {
340  // if only one hit was found the trajectory measurement is on that sensor surface, and the other surface from
341  // the matched layer should be added to the study as well
342  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 1));
343  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator, 2));
344  if (DEBUG) cout << " found a hit with a missing partner" << endl;
345  } else {
346  //only add one TM for the single surface and the other will be added in the next iteration
347  TMs.push_back(TrajectoryAtInvalidHit(*itm, tTopo, tkgeom, propagator));
348  }
349 
351  //Now check for tracks at TOB6 and TEC9
352 
353  // to make sure we only propagate on the last TOB5 hit check the next entry isn't also in TOB5
354  // to avoid bias, make sure the TOB5 hit is valid (an invalid hit on TOB5 could only exist with a valid hit on TOB6)
355 
356  bool isValid = theInHit->isValid();
357  bool isLast = (itm==(TMeas.end()-1));
358  bool isLastTOB5 = true;
359  if (!isLast) {
360  if ( checkLayer((++itm)->recHit()->geographicalId().rawId(), tTopo) == 9 ) isLastTOB5 = false;
361  else isLastTOB5 = true;
362  --itm;
363  }
364 
365  if ( TKlayers==9 && isValid && isLastTOB5 ) {
366  // if ( TKlayers==9 && itm==TMeas.rbegin()) {
367  // if ( TKlayers==9 && (itm==TMeas.back()) ) { // to check for only the last entry in the trajectory for propagation
368  std::vector< BarrelDetLayer const*> barrelTOBLayers = measurementTrackerHandle->geometricSearchTracker()->tobLayers() ;
369  const DetLayer* tob6 = barrelTOBLayers[barrelTOBLayers.size()-1];
370  const MeasurementEstimator* estimator = est.product();
371  const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(*measurementTrackerHandle, *measurementTrackerEvent);
372  const TrajectoryStateOnSurface tsosTOB5 = itm->updatedState();
373  vector<TrajectoryMeasurement> tmp = theLayerMeasurements->measurements(*tob6, tsosTOB5, *thePropagator, *estimator);
374 
375  if ( !tmp.empty()) {
376  if (DEBUG) cout << "size of TM from propagation = " << tmp.size() << endl;
377 
378  // take the last of the TMs, which is always an invalid hit
379  // if no detId is available, ie detId==0, then no compatible layer was crossed
380  // otherwise, use that TM for the efficiency measurement
381  TrajectoryMeasurement tob6TM(tmp.back());
382  auto tob6Hit = tob6TM.recHit();
383 
384  if (tob6Hit->geographicalId().rawId()!=0) {
385  if (DEBUG) cout << "tob6 hit actually being added to TM vector" << endl;
386  TMs.push_back(TrajectoryAtInvalidHit(tob6TM, tTopo, tkgeom, propagator));
387  }
388  }
389  }
390 
391  bool isLastTEC8 = true;
392  if (!isLast) {
393  if ( checkLayer((++itm)->recHit()->geographicalId().rawId(), tTopo) == 21 ) isLastTEC8 = false;
394  else isLastTEC8 = true;
395  --itm;
396  }
397 
398  if ( TKlayers==21 && isValid && isLastTEC8 ) {
399 
400  std::vector< const ForwardDetLayer*> posTecLayers = measurementTrackerHandle->geometricSearchTracker()->posTecLayers() ;
401  const DetLayer* tec9pos = posTecLayers[posTecLayers.size()-1];
402  std::vector< const ForwardDetLayer*> negTecLayers = measurementTrackerHandle->geometricSearchTracker()->negTecLayers() ;
403  const DetLayer* tec9neg = negTecLayers[negTecLayers.size()-1];
404 
405  const MeasurementEstimator* estimator = est.product();
406  const LayerMeasurements* theLayerMeasurements = new LayerMeasurements(*measurementTrackerHandle, *measurementTrackerEvent);
407  const TrajectoryStateOnSurface tsosTEC9 = itm->updatedState();
408 
409  // check if track on positive or negative z
410  if (!iidd == StripSubdetector::TEC) cout << "there is a problem with TEC 9 extrapolation" << endl;
411 
412  //cout << " tec9 id = " << iidd << " and side = " << tTopo->tecSide(iidd) << endl;
413  vector<TrajectoryMeasurement> tmp;
414  if ( tTopo->tecSide(iidd) == 1 ) {
415  tmp = theLayerMeasurements->measurements(*tec9neg, tsosTEC9, *thePropagator, *estimator);
416  //cout << "on negative side" << endl;
417  }
418  if ( tTopo->tecSide(iidd) == 2 ) {
419  tmp = theLayerMeasurements->measurements(*tec9pos, tsosTEC9, *thePropagator, *estimator);
420  //cout << "on positive side" << endl;
421  }
422 
423  if ( !tmp.empty()) {
424  // take the last of the TMs, which is always an invalid hit
425  // if no detId is available, ie detId==0, then no compatible layer was crossed
426  // otherwise, use that TM for the efficiency measurement
427  TrajectoryMeasurement tec9TM(tmp.back());
428  auto tec9Hit = tec9TM.recHit();
429 
430  unsigned int tec9id = tec9Hit->geographicalId().rawId();
431  if (DEBUG) cout << "tec9id = " << tec9id << " is Double sided = " << isDoubleSided(tec9id, tTopo) << " and 0x3 = " << (tec9id & 0x3) << endl;
432 
433  if (tec9Hit->geographicalId().rawId()!=0) {
434  if (DEBUG) cout << "tec9 hit actually being added to TM vector" << endl;
435  // in tec the hit can be single or doubled sided. whenever the invalid hit at the end of vector of TMs is
436  // double sided it is always on the matched surface, so we need to split it into the true sensor surfaces
437  if (isDoubleSided(tec9id, tTopo)) {
438  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator, 1));
439  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator, 2));
440  } else
441  TMs.push_back(TrajectoryAtInvalidHit(tec9TM, tTopo, tkgeom, propagator));
442  }
443  } //else cout << "tec9 tmp empty" << endl;
444  }
445 
447 
448  // Modules Constraints
449 
450  for(std::vector<TrajectoryAtInvalidHit>::const_iterator TM=TMs.begin();TM!=TMs.end();++TM) {
451 
452  // --> Get trajectory from combinatedState
453  iidd = TM->monodet_id();
454  if (DEBUG) cout << "setting iidd = " << iidd << " before checking efficiency and ";
455 
456  xloc = TM->localX();
457  yloc = TM->localY();
458 
459  xErr = TM->localErrorX();
460  yErr = TM->localErrorY();
461 
462  angleX = atan( TM->localDxDz() );
463  angleY = atan( TM->localDyDz() );
464 
465  xglob = TM->globalX();
466  yglob = TM->globalY();
467  zglob = TM->globalZ();
468  withinAcceptance = TM->withinAcceptance();
469 
470  trajHitValid = TM->validHit();
471 
472  // reget layer from iidd here, to account for TOB 6 and TEC 9 TKlayers being off
473  TKlayers = checkLayer(iidd, tTopo);
474 
475  if ( (layers == TKlayers) || (layers == 0) ) { // Look at the layer not used to reconstruct the track
476  whatlayer = TKlayers;
477  if (DEBUG) cout << "Looking at layer under study" << endl;
478  TrajGlbX = 0.0; TrajGlbY = 0.0; TrajGlbZ = 0.0; ModIsBad = 2; Id = 0; SiStripQualBad = 0;
479  run = 0; event = 0; TrajLocX = 0.0; TrajLocY = 0.0; TrajLocErrX = 0.0; TrajLocErrY = 0.0;
480  TrajLocAngleX = -999.0; TrajLocAngleY = -999.0; ResX = 0.0; ResXSig = 0.0;
481  ClusterLocX = 0.0; ClusterLocY = 0.0; ClusterLocErrX = 0.0; ClusterLocErrY = 0.0; ClusterStoN = 0.0;
482  bunchx = 0;
483 
484  // RPhi RecHit Efficiency
485 
486  if (input.size() > 0 ) {
487  if (DEBUG) cout << "Checking clusters with size = " << input.size() << endl;
488  int nClusters = 0;
489  std::vector< std::vector<float> > VCluster_info; //fill with X residual, X residual pull, local X, sig(X), local Y, sig(Y), StoN
490  for (edmNew::DetSetVector<SiStripCluster>::const_iterator DSViter = input.begin(); DSViter != input.end(); DSViter++) {
491  // DSViter is a vector of SiStripClusters located on a single module
492  //if (DEBUG) cout << "the ID from the DSViter = " << DSViter->id() << endl;
493  unsigned int ClusterId = DSViter->id();
494  if (ClusterId == iidd) {
495  if (DEBUG) cout << "found (ClusterId == iidd) with ClusterId = " << ClusterId << " and iidd = " << iidd << endl;
496  DetId ClusterDetId(ClusterId);
497  const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tkgeom->idToDetUnit(ClusterDetId);
498 
499  for(edmNew::DetSet<SiStripCluster>::const_iterator iter=DSViter->begin();iter!=DSViter->end();++iter) {
500  //iter is a single SiStripCluster
501  StripClusterParameterEstimator::LocalValues parameters=stripcpe.localParameters(*iter,*stripdet);
502  float res = (parameters.first.x() - xloc);
503  float sigma = checkConsistency(parameters , xloc, xErr);
504  // The consistency is probably more accurately measured with the Chi2MeasurementEstimator. To use it
505  // you need a TransientTrackingRecHit instead of the cluster
506  //theEstimator= new Chi2MeasurementEstimator(30);
507  //const Chi2MeasurementEstimator *theEstimator(100);
508  //theEstimator->estimate(TM->tsos(), TransientTrackingRecHit);
509 
510  SiStripClusterInfo clusterInfo = SiStripClusterInfo(*iter, es, ClusterId);
511  // signal to noise from SiStripClusterInfo not working in 225. I'll fix this after the interface
512  // redesign in 300 -ku
513  //float cluster_info[7] = {res, sigma, parameters.first.x(), sqrt(parameters.second.xx()), parameters.first.y(), sqrt(parameters.second.yy()), signal_to_noise};
514  std::vector< float > cluster_info;
515  cluster_info.push_back(res);
516  cluster_info.push_back(sigma);
517  cluster_info.push_back(parameters.first.x());
518  cluster_info.push_back(sqrt(parameters.second.xx()));
519  cluster_info.push_back(parameters.first.y());
520  cluster_info.push_back(sqrt(parameters.second.yy()));
521  cluster_info.push_back( clusterInfo.signalOverNoise() );
522  //cluster_info.push_back( clusterInfo.getSignalOverNoise() );
523  VCluster_info.push_back(cluster_info);
524  nClusters++;
525  if (DEBUG) cout << "Have ID match. residual = " << VCluster_info.back()[0] << " res sigma = " << VCluster_info.back()[1] << endl;
526  if (DEBUG) cout << "trajectory measurement compatability estimate = " << (*itm).estimate() << endl;
527  if (DEBUG) cout << "hit position = " << parameters.first.x() << " hit error = " << sqrt(parameters.second.xx()) << " trajectory position = " << xloc << " traj error = " << xErr << endl;
528  }
529  }
530  }
531  float FinalResSig = 1000.0;
532  float FinalCluster[7]= {1000.0, 1000.0, 0.0, 0.0, 0.0, 0.0, 0.0};
533  if (nClusters > 0) {
534  if (DEBUG) cout << "found clusters > 0" << endl;
535  if (nClusters > 1) {
536  //get the smallest one
537  vector< vector<float> >::iterator ires;
538  for (ires=VCluster_info.begin(); ires!=VCluster_info.end(); ires++){
539  if ( abs((*ires)[1]) < abs(FinalResSig)) {
540  FinalResSig = (*ires)[1];
541  for (unsigned int i = 0; i<ires->size(); i++) {
542  if (DEBUG) cout << "filling final cluster. i = " << i << " before fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
543  FinalCluster[i] = (*ires)[i];
544  if (DEBUG) cout << "filling final cluster. i = " << i << " after fill FinalCluster[i]=" << FinalCluster[i] << " and (*ires)[i] =" << (*ires)[i] << endl;
545  }
546  }
547  if (DEBUG) cout << "iresidual = " << (*ires)[0] << " isigma = " << (*ires)[1] << " and FinalRes = " << FinalCluster[0] << endl;
548  }
549  }
550  else {
551  FinalResSig = VCluster_info.at(0)[1];
552  for (unsigned int i = 0; i<VCluster_info.at(0).size(); i++) {
553  FinalCluster[i] = VCluster_info.at(0)[i];
554  }
555  }
556  nClusters=0;
557  VCluster_info.clear();
558  }
559 
560  if (DEBUG) cout << "Final residual in X = " << FinalCluster[0] << "+-" << (FinalCluster[0]/FinalResSig) << endl;
561  if (DEBUG) cout << "Checking location of trajectory: abs(yloc) = " << abs(yloc) << " abs(xloc) = " << abs(xloc) << endl;
562  if (DEBUG) cout << "Checking location of cluster hit: yloc = " << FinalCluster[4] << "+-" << FinalCluster[5] << " xloc = " << FinalCluster[2] << "+-" << FinalCluster[3] << endl;
563  if (DEBUG) cout << "Final cluster signal to noise = " << FinalCluster[6] << endl;
564 
565  float exclusionWidth = 0.4;
566  float TOBexclusion = 0.0;
567  float TECexRing5 = -0.89;
568  float TECexRing6 = -0.56;
569  float TECexRing7 = 0.60;
570  //Added by Chris Edelmaier to do TEC bonding exclusion
571  int subdetector = ((iidd>>25) & 0x7);
572  int ringnumber = ((iidd>>5) & 0x7);
573 
574  //New TOB and TEC bonding region exclusion zone
575  if((TKlayers >= 5 && TKlayers < 11)||((subdetector == 6)&&( (ringnumber >= 5)&&(ringnumber <=7) ))) {
576  //There are only 2 cases that we need to exclude for
577  float highzone = 0.0;
578  float lowzone = 0.0;
579  float higherr = yloc + 5.0*yErr;
580  float lowerr = yloc - 5.0*yErr;
581  if(TKlayers >= 5 && TKlayers < 11) {
582  //TOB zone
583  highzone = TOBexclusion + exclusionWidth;
584  lowzone = TOBexclusion - exclusionWidth;
585  } else if (ringnumber == 5) {
586  //TEC ring 5
587  highzone = TECexRing5 + exclusionWidth;
588  lowzone = TECexRing5 - exclusionWidth;
589  } else if (ringnumber == 6) {
590  //TEC ring 6
591  highzone = TECexRing6 + exclusionWidth;
592  lowzone = TECexRing6 - exclusionWidth;
593  } else if (ringnumber == 7) {
594  //TEC ring 7
595  highzone = TECexRing7 + exclusionWidth;
596  lowzone = TECexRing7 - exclusionWidth;
597  }
598  //Now that we have our exclusion region, we just have to properly identify it
599  if((highzone <= higherr)&&(highzone >= lowerr)) withinAcceptance = false;
600  if((lowzone >= lowerr)&&(lowzone <= higherr)) withinAcceptance = false;
601  if((higherr <= highzone)&&(higherr >= lowzone)) withinAcceptance = false;
602  if((lowerr >= lowzone)&&(lowerr <= highzone)) withinAcceptance = false;
603  }
604 
605  // fill ntuple varibles
606  //get global position from module id number iidd
607  TrajGlbX = xglob;
608  TrajGlbY = yglob;
609  TrajGlbZ = zglob;
610  Id = iidd;
611  run = run_nr;
612  event = ev_nr;
613  bunchx = bunch_nr;
614  //if ( SiStripQuality_->IsModuleBad(iidd) ) {
615  if ( SiStripQuality_->getBadApvs(iidd)!=0 ) {
616  SiStripQualBad = 1;
617  if(DEBUG) cout << "strip is bad from SiStripQuality" << endl;
618  } else {
619  SiStripQualBad = 0;
620  if(DEBUG) cout << "strip is good from SiStripQuality" << endl;
621  }
622 
623  //check for FED-detected errors and include those in SiStripQualBad
624  for (unsigned int ii=0;ii< fedErrorIds->size();ii++) {
625  if (iidd == (*fedErrorIds)[ii].rawId() )
626  SiStripQualBad = 1;
627  }
628 
629  TrajLocX = xloc;
630  TrajLocY = yloc;
631  TrajLocErrX = xErr;
632  TrajLocErrY = yErr;
633  TrajLocAngleX = angleX;
634  TrajLocAngleY = angleY;
635  ResX = FinalCluster[0];
636  ResXSig = FinalResSig;
637  if (FinalResSig != FinalCluster[1]) if (DEBUG) cout << "Problem with best cluster selection because FinalResSig = " << FinalResSig << " and FinalCluster[1] = " << FinalCluster[1] << endl;
638  ClusterLocX = FinalCluster[2];
639  ClusterLocY = FinalCluster[4];
640  ClusterLocErrX = FinalCluster[3];
641  ClusterLocErrY = FinalCluster[5];
642  ClusterStoN = FinalCluster[6];
643 
644  if (DEBUG) cout << "before check good" << endl;
645 
646  if ( FinalResSig < 999.0) { //could make requirement on track/hit consistency, but for
647  //now take anything with a hit on the module
648  if (DEBUG) cout << "hit being counted as good " << FinalCluster[0] << " FinalRecHit " <<
649  iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd <<
650  " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
651  ModIsBad = 0;
652  traj->Fill();
653  }
654  else {
655  if (DEBUG) cout << "hit being counted as bad ######### Invalid RPhi FinalResX " << FinalCluster[0] << " FinalRecHit " <<
656  iidd << " TKlayers " << TKlayers << " xloc " << xloc << " yloc " << yloc << " module " << iidd <<
657  " matched/stereo/rphi = " << ((iidd & 0x3)==0) << "/" << ((iidd & 0x3)==1) << "/" << ((iidd & 0x3)==2) << endl;
658  ModIsBad = 1;
659  traj->Fill();
660 
661  if (DEBUG) cout << " RPhi Error " << sqrt(xErr*xErr + yErr*yErr) << " ErrorX " << xErr << " yErr " << yErr << endl;
662  } if (DEBUG) cout << "after good location check" << endl;
663  } if (DEBUG) cout << "after list of clusters" << endl;
664  } if (DEBUG) cout << "After layers=TKLayers if" << endl;
665  } if (DEBUG) cout << "After looping over TrajAtValidHit list" << endl;
666  } if (DEBUG) cout << "end TMeasurement loop" << endl;
667  } if (DEBUG) cout << "end of trajectories loop" << endl;
668  }
669 }
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
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
int i
Definition: DBlmapReader.cc:9
float TrajGlbZ
Definition: HitEff.h:71
dictionary parameters
Definition: Parameters.py:2
bool check2DPartner(unsigned int iidd, const std::vector< TrajectoryMeasurement > &traj)
Definition: HitEff.cc:718
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
float ClusterLocErrX
Definition: HitEff.h:73
float ClusterLocErrY
Definition: HitEff.h:73
ConstRecHitPointer const & recHit() const
float TrajGlbY
Definition: HitEff.h:71
float TrajLocX
Definition: HitEff.h:72
float ResXSig
Definition: HitEff.h:74
float chi2
Definition: HitEff.h:77
int ires[2]
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
int bunchCrossing() const
Definition: EventBase.h:66
float ClusterStoN
Definition: HitEff.h:73
data_type const * const_iterator
Definition: DetSetNew.h:30
float TrajLocErrX
Definition: HitEff.h:72
int ii
Definition: cuy.py:588
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
bool DEBUG
Definition: HitEff.h:66
unsigned int bunchx
Definition: HitEff.h:78
static std::string const input
Definition: EdmProvDump.cc:43
float p
Definition: HitEff.h:77
float timeECAL
Definition: HitEff.h:81
float signalOverNoise() const
bool isDoubleSided(unsigned int iidd, const TrackerTopology *tTopo) const
Definition: HitEff.cc:686
float timeDTErr
Definition: HitEff.h:79
int nHits
Definition: HitEff.h:76
unsigned int Id
Definition: HitEff.h:75
T sqrt(T t)
Definition: SSEVec.h:48
int timeDTDOF
Definition: HitEff.h:80
float timeAtIpInOutErr
Definition: MuonTime.h:15
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:78
unsigned int ModIsBad
Definition: HitEff.h:75
float TrajGlbX
Definition: HitEff.h:71
bool withinAcceptance
Definition: HitEff.h:75
unsigned int run
Definition: HitEff.h:78
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
float ClusterLocX
Definition: HitEff.h:73
int events
Definition: HitEff.h:63
float TrajLocAngleX
Definition: HitEff.h:72
Definition: DetId.h:18
double checkConsistency(const StripClusterParameterEstimator::LocalValues &parameters, double xx, double xerr)
Definition: HitEff.cc:679
T const * product() const
Definition: Handle.h:81
unsigned int whatlayer
Definition: HitEff.h:67
unsigned int checkLayer(unsigned int iidd, const TrackerTopology *tTopo)
Definition: HitEff.cc:734
T const * product() const
Definition: ESHandle.h:86
float timeDT
Definition: HitEff.h:79
unsigned int layers
Definition: HitEff.h:65
float pT
Definition: HitEff.h:77
float ResX
Definition: HitEff.h:74
std::pair< LocalPoint, LocalError > LocalValues
unsigned int SiStripQualBad
Definition: HitEff.h:75
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
edm::EventID id() const
Definition: EventBase.h:60
float TrajLocErrY
Definition: HitEff.h:72
float TrajLocAngleY
Definition: HitEff.h:72
tuple cout
Definition: gather_cfg.py:121
float ClusterLocY
Definition: HitEff.h:73
TTree * traj
Definition: HitEff.h:62
edm::ParameterSet conf_
Definition: HitEff.h:60
float timeAtIpInOut
Definition: MuonTime.h:14
int EventTrackCKF
Definition: HitEff.h:63
int nLostHits
Definition: HitEff.h:76
float TrajLocY
Definition: HitEff.h:72
unsigned int tecSide(const DetId &id) const
void HitEff::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 83 of file HitEff.cc.

References bunchx, chi2, ClusterLocErrX, ClusterLocErrY, ClusterLocX, ClusterLocY, ClusterStoN, dedx, dedxNOM, event, events, EventTrackCKF, Id, istep, TFileService::make(), ModIsBad, nHits, nLostHits, p, pT, ResX, ResXSig, run, SiStripQualBad, timeDT, timeDTDOF, timeDTErr, timeECAL, tquality, traj, TrajGlbX, TrajGlbY, TrajGlbZ, trajHitValid, TrajLocAngleX, TrajLocAngleY, TrajLocErrX, TrajLocErrY, TrajLocX, TrajLocY, whatlayer, and withinAcceptance.

83  {
84 
86 
87  traj = fs->make<TTree>("traj","tree of trajectory positions");
88  traj->Branch("TrajGlbX",&TrajGlbX,"TrajGlbX/F");
89  traj->Branch("TrajGlbY",&TrajGlbY,"TrajGlbY/F");
90  traj->Branch("TrajGlbZ",&TrajGlbZ,"TrajGlbZ/F");
91  traj->Branch("TrajLocX",&TrajLocX,"TrajLocX/F");
92  traj->Branch("TrajLocY",&TrajLocY,"TrajLocY/F");
93  traj->Branch("TrajLocErrX",&TrajLocErrX,"TrajLocErrX/F");
94  traj->Branch("TrajLocErrY",&TrajLocErrY,"TrajLocErrY/F");
95  traj->Branch("TrajLocAngleX",&TrajLocAngleX,"TrajLocAngleX/F");
96  traj->Branch("TrajLocAngleY",&TrajLocAngleY,"TrajLocAngleY/F");
97  traj->Branch("ClusterLocX",&ClusterLocX,"ClusterLocX/F");
98  traj->Branch("ClusterLocY",&ClusterLocY,"ClusterLocY/F");
99  traj->Branch("ClusterLocErrX",&ClusterLocErrX,"ClusterLocErrX/F");
100  traj->Branch("ClusterLocErrY",&ClusterLocErrY,"ClusterLocErrY/F");
101  traj->Branch("ClusterStoN",&ClusterStoN,"ClusterStoN/F");
102  traj->Branch("ResX",&ResX,"ResX/F");
103  traj->Branch("ResXSig",&ResXSig,"ResXSig/F");
104  traj->Branch("ModIsBad",&ModIsBad,"ModIsBad/i");
105  traj->Branch("SiStripQualBad",&SiStripQualBad,"SiStripQualBad/i");
106  traj->Branch("withinAcceptance",&withinAcceptance,"withinAcceptance/O");
107  traj->Branch("nHits",&nHits,"nHits/I");
108  traj->Branch("nLostHits",&nLostHits,"nLostHits/I");
109  traj->Branch("chi2",&chi2,"chi2/F");
110  traj->Branch("p",&p,"p/F");
111  traj->Branch("pT",&pT,"pT/F");
112  traj->Branch("trajHitValid", &trajHitValid, "trajHitValid/i");
113  traj->Branch("Id",&Id,"Id/i");
114  traj->Branch("run",&run,"run/i");
115  traj->Branch("event",&event,"event/i");
116  traj->Branch("layer",&whatlayer,"layer/i");
117  traj->Branch("timeDT",&timeDT,"timeDT/F");
118  traj->Branch("timeDTErr",&timeDTErr,"timeDTErr/F");
119  traj->Branch("timeDTDOF",&timeDTDOF,"timeDTDOF/I");
120  traj->Branch("timeECAL",&timeECAL,"timeECAL/F");
121  traj->Branch("dedx",&dedx,"dedx/F");
122  traj->Branch("dedxNOM",&dedxNOM,"dedxNOM/I");
123  traj->Branch("tquality",&tquality,"tquality/I");
124  traj->Branch("istep",&istep,"istep/I");
125  traj->Branch("bunchx",&bunchx,"bunchx/I");
126 
127  events = 0;
128  EventTrackCKF = 0;
129 
130 }
float TrajGlbZ
Definition: HitEff.h:71
float ClusterLocErrX
Definition: HitEff.h:73
float ClusterLocErrY
Definition: HitEff.h:73
float TrajGlbY
Definition: HitEff.h:71
float TrajLocX
Definition: HitEff.h:72
float ResXSig
Definition: HitEff.h:74
float chi2
Definition: HitEff.h:77
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
float ClusterStoN
Definition: HitEff.h:73
float TrajLocErrX
Definition: HitEff.h:72
unsigned int bunchx
Definition: HitEff.h:78
float p
Definition: HitEff.h:77
float timeECAL
Definition: HitEff.h:81
float timeDTErr
Definition: HitEff.h:79
int nHits
Definition: HitEff.h:76
unsigned int Id
Definition: HitEff.h:75
int timeDTDOF
Definition: HitEff.h:80
unsigned int trajHitValid
Definition: HitEff.h:78
unsigned int ModIsBad
Definition: HitEff.h:75
int istep
Definition: HitEff.h:84
float TrajGlbX
Definition: HitEff.h:71
bool withinAcceptance
Definition: HitEff.h:75
unsigned int run
Definition: HitEff.h:78
float ClusterLocX
Definition: HitEff.h:73
int events
Definition: HitEff.h:63
float TrajLocAngleX
Definition: HitEff.h:72
int dedxNOM
Definition: HitEff.h:82
unsigned int event
Definition: HitEff.h:78
unsigned int whatlayer
Definition: HitEff.h:67
float timeDT
Definition: HitEff.h:79
float pT
Definition: HitEff.h:77
float ResX
Definition: HitEff.h:74
unsigned int SiStripQualBad
Definition: HitEff.h:75
float TrajLocErrY
Definition: HitEff.h:72
float TrajLocAngleY
Definition: HitEff.h:72
float dedx
Definition: HitEff.h:81
float ClusterLocY
Definition: HitEff.h:73
TTree * traj
Definition: HitEff.h:62
int EventTrackCKF
Definition: HitEff.h:63
int nLostHits
Definition: HitEff.h:76
float TrajLocY
Definition: HitEff.h:72
int tquality
Definition: HitEff.h:83
bool HitEff::check2DPartner ( unsigned int  iidd,
const std::vector< TrajectoryMeasurement > &  traj 
)

Definition at line 718 of file HitEff.cc.

Referenced by analyze().

718  {
719  unsigned int partner_iidd = 0;
720  bool found2DPartner = false;
721  // first get the id of the other detector
722  if ((iidd & 0x3)==1) partner_iidd = iidd+1;
723  if ((iidd & 0x3)==2) partner_iidd = iidd-1;
724  // next look in the trajectory measurements for a measurement from that detector
725  // loop through trajectory measurements to find the partner_iidd
726  for (std::vector<TrajectoryMeasurement>::const_iterator iTM=traj.begin(); iTM!=traj.end(); ++iTM) {
727  if (iTM->recHit()->geographicalId().rawId()==partner_iidd) {
728  found2DPartner = true;
729  }
730  }
731  return found2DPartner;
732 }
TTree * traj
Definition: HitEff.h:62
double HitEff::checkConsistency ( const StripClusterParameterEstimator::LocalValues parameters,
double  xx,
double  xerr 
)

Definition at line 679 of file HitEff.cc.

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

Referenced by analyze().

679  {
680  double error = sqrt(parameters.second.xx() + xerr*xerr);
681  double separation = abs(parameters.first.x() - xx);
682  double consistency = separation/error;
683  return consistency;
684 }
dictionary parameters
Definition: Parameters.py:2
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int HitEff::checkLayer ( unsigned int  iidd,
const TrackerTopology tTopo 
)

Definition at line 734 of file HitEff.cc.

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

Referenced by analyze().

734  {
736  unsigned int subid=strip.subdetId();
737  if (subid == StripSubdetector::TIB) {
738 
739  return tTopo->tibLayer(iidd);
740  }
741  if (subid == StripSubdetector::TOB) {
742 
743  return tTopo->tobLayer(iidd) + 4 ;
744  }
745  if (subid == StripSubdetector::TID) {
746 
747  return tTopo->tidWheel(iidd) + 10;
748  }
749  if (subid == StripSubdetector::TEC) {
750 
751  return tTopo->tecWheel(iidd) + 13 ;
752  }
753  return 0;
754 }
unsigned int tibLayer(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
unsigned int tecWheel(const DetId &id) const
unsigned int tobLayer(const DetId &id) const
void HitEff::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 671 of file HitEff.cc.

References gather_cfg::cout, DEBUG, events, EventTrackCKF, and traj.

671  {
672  traj->GetDirectory()->cd();
673  traj->Write();
674 
675  if (DEBUG) cout << " Events Analysed " << events << endl;
676  if (DEBUG) cout << " Number Of Tracked events " << EventTrackCKF << endl;
677 }
bool DEBUG
Definition: HitEff.h:66
int events
Definition: HitEff.h:63
tuple cout
Definition: gather_cfg.py:121
TTree * traj
Definition: HitEff.h:62
int EventTrackCKF
Definition: HitEff.h:63
bool HitEff::isDoubleSided ( unsigned int  iidd,
const TrackerTopology tTopo 
) const

Definition at line 686 of file HitEff.cc.

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

Referenced by analyze().

686  {
688  unsigned int subid=strip.subdetId();
689  unsigned int layer = 0;
690  if (subid == StripSubdetector::TIB) {
691 
692  layer = tTopo->tibLayer(iidd);
693  if (layer == 1 || layer == 2) return true;
694  else return false;
695  }
696  else if (subid == StripSubdetector::TOB) {
697 
698  layer = tTopo->tobLayer(iidd) + 4 ;
699  if (layer == 5 || layer == 6) return true;
700  else return false;
701  }
702  else if (subid == StripSubdetector::TID) {
703 
704  layer = tTopo->tidRing(iidd) + 10;
705  if (layer == 11 || layer == 12) return true;
706  else return false;
707  }
708  else if (subid == StripSubdetector::TEC) {
709 
710  layer = tTopo->tecRing(iidd) + 13 ;
711  if (layer == 14 || layer == 15 || layer == 18) return true;
712  else return false;
713  }
714  else
715  return false;
716 }
unsigned int tibLayer(const DetId &id) const
unsigned int tidRing(const DetId &id) const
unsigned int tecRing(const DetId &id) const
ring id
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
unsigned int tobLayer(const DetId &id) const

Member Data Documentation

unsigned int HitEff::bunchx
private

Definition at line 78 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::chi2
private

Definition at line 77 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocErrX
private

Definition at line 73 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocErrY
private

Definition at line 73 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocX
private

Definition at line 73 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterLocY
private

Definition at line 73 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ClusterStoN
private

Definition at line 73 of file HitEff.h.

Referenced by analyze(), and beginJob().

edm::ParameterSet HitEff::conf_
private

Definition at line 60 of file HitEff.h.

Referenced by analyze(), and HitEff().

bool HitEff::DEBUG
private

Definition at line 66 of file HitEff.h.

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

float HitEff::dedx
private

Definition at line 81 of file HitEff.h.

Referenced by beginJob().

int HitEff::dedxNOM
private

Definition at line 82 of file HitEff.h.

Referenced by beginJob().

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

Definition at line 63 of file HitEff.h.

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

unsigned int HitEff::Id
private

Definition at line 75 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::istep
private

Definition at line 84 of file HitEff.h.

Referenced by beginJob().

unsigned int HitEff::layers
private

Definition at line 65 of file HitEff.h.

Referenced by analyze(), and HitEff().

unsigned int HitEff::ModIsBad
private

Definition at line 75 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::nHits
private

Definition at line 76 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::nLostHits
private

Definition at line 76 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::p
private
float HitEff::pT
private

Definition at line 77 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ResX
private

Definition at line 74 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::ResXSig
private

Definition at line 74 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::run
private
unsigned int HitEff::SiStripQualBad
private

Definition at line 75 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::timeDT
private

Definition at line 79 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::timeDTDOF
private

Definition at line 80 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::timeDTErr
private

Definition at line 79 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::timeECAL
private

Definition at line 81 of file HitEff.h.

Referenced by analyze(), and beginJob().

int HitEff::tquality
private

Definition at line 83 of file HitEff.h.

Referenced by beginJob().

TTree* HitEff::traj
private

Definition at line 62 of file HitEff.h.

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

float HitEff::TrajGlbX
private

Definition at line 71 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajGlbY
private

Definition at line 71 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajGlbZ
private

Definition at line 71 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::trajHitValid
private

Definition at line 78 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocAngleX
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocAngleY
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocErrX
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocErrY
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocX
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and beginJob().

float HitEff::TrajLocY
private

Definition at line 72 of file HitEff.h.

Referenced by analyze(), and beginJob().

unsigned int HitEff::whatlayer
private

Definition at line 67 of file HitEff.h.

Referenced by analyze(), and beginJob().

bool HitEff::withinAcceptance
private

Definition at line 75 of file HitEff.h.

Referenced by analyze(), and beginJob().