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

Referenced by analyze().

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

Definition at line 673 of file HitEff.cc.

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

Referenced by analyze().

673  {
674  double error = sqrt(parameters.second.xx() + xerr*xerr);
675  double separation = abs(parameters.first.x() - xx);
676  double consistency = separation/error;
677  return consistency;
678 }
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 728 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().

728  {
730  unsigned int subid=strip.subdetId();
731  if (subid == StripSubdetector::TIB) {
732 
733  return tTopo->tibLayer(iidd);
734  }
735  if (subid == StripSubdetector::TOB) {
736 
737  return tTopo->tobLayer(iidd) + 4 ;
738  }
739  if (subid == StripSubdetector::TID) {
740 
741  return tTopo->tidWheel(iidd) + 10;
742  }
743  if (subid == StripSubdetector::TEC) {
744 
745  return tTopo->tecWheel(iidd) + 13 ;
746  }
747  return 0;
748 }
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 665 of file HitEff.cc.

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

665  {
666  traj->GetDirectory()->cd();
667  traj->Write();
668 
669  if (DEBUG) cout << " Events Analysed " << events << endl;
670  if (DEBUG) cout << " Number Of Tracked events " << EventTrackCKF << endl;
671 }
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 680 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().

680  {
682  unsigned int subid=strip.subdetId();
683  unsigned int layer = 0;
684  if (subid == StripSubdetector::TIB) {
685 
686  layer = tTopo->tibLayer(iidd);
687  if (layer == 1 || layer == 2) return true;
688  else return false;
689  }
690  else if (subid == StripSubdetector::TOB) {
691 
692  layer = tTopo->tobLayer(iidd) + 4 ;
693  if (layer == 5 || layer == 6) return true;
694  else return false;
695  }
696  else if (subid == StripSubdetector::TID) {
697 
698  layer = tTopo->tidRing(iidd) + 10;
699  if (layer == 11 || layer == 12) return true;
700  else return false;
701  }
702  else if (subid == StripSubdetector::TEC) {
703 
704  layer = tTopo->tecRing(iidd) + 13 ;
705  if (layer == 14 || layer == 15 || layer == 18) return true;
706  else return false;
707  }
708  else
709  return false;
710 }
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().