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
DTChamberEfficiency Class Reference

#include <DTChamberEfficiency.h>

Inheritance diagram for DTChamberEfficiency:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 
void beginJob ()
 
void beginRun (const edm::Run &, const edm::EventSetup &)
 
 DTChamberEfficiency (const edm::ParameterSet &pset)
 
void endJob ()
 
 ~DTChamberEfficiency ()
 
- 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
 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
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void bookHistos ()
 
bool chamberSelection (const DetId &idDetLay, reco::TransientTrack &trans_track) const
 
std::vector< const DetLayer * > compatibleLayers (const NavigationSchool &navigationSchool, const DetLayer *initialLayer, const FreeTrajectoryState &fts, PropagationDirection propDir)
 
edm::ESHandle< Propagatorpropagator () const
 
MeasurementContainer segQualityCut (const MeasurementContainer &seg_list) const
 

Private Attributes

bool debug
 
edm::ESHandle< DTGeometrydtGeom
 
std::vector< std::vector
< MonitorElement * > > 
histosPerW
 
edm::InputTag labelRPCRecHits
 
edm::ESHandle< MagneticFieldmagfield
 
edm::InputTag thecscSegments
 
DQMStoretheDbe
 
edm::InputTag thedt4DSegments
 
Chi2MeasurementEstimatortheEstimator
 
double theMaxChi2
 
MuonDetLayerMeasurementstheMeasurementExtractor
 
int theMinNrec
 
std::string theNavigationType
 
double theNSigma
 
MuonServiceProxytheService
 
edm::ESHandle
< GlobalTrackingGeometry
theTrackingGeometry
 
edm::InputTag theTracksLabel_
 
edm::EDGetTokenT
< reco::TrackCollection
theTracksToken_
 

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

Description:

This class provides the histograms for the calculation of the efficiency of muons reconstruction in the DTs. It is applicable both in presence or absence of a magnetic field. Histos are 2D Sector vs Chamber plots for each wheel

Author
: Mario Pelliccioni - INFN Torino pelli.nosp@m.cci@.nosp@m.cern..nosp@m.ch $date : 05/12/2008 16:51:04 CET $

Modification:

Definition at line 55 of file DTChamberEfficiency.h.

Constructor & Destructor Documentation

DTChamberEfficiency::DTChamberEfficiency ( const edm::ParameterSet pset)

Definition at line 55 of file DTChamberEfficiency.cc.

References Chi2MeasurementEstimatorESProducer_cfi::Chi2MeasurementEstimator, debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), MuonServiceProxy_cff::MuonServiceProxy, cppFunctionSkipper::operator, and createJobs::theNSigma.

56 {
57  // Get the debug parameter for verbose output
58  debug = pSet.getUntrackedParameter<bool>("debug",false);
59 
60  LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiency")
61  << "DTChamberEfficiency: constructor called";
62 
63  // Get the DQM needed services
65  theDbe->setCurrentFolder("DT/05-ChamberEff/Task");
66 
67  // service parameters
68  ParameterSet serviceParameters = pSet.getParameter<ParameterSet>("ServiceParameters");
69  theService = new MuonServiceProxy(serviceParameters);
70 
71  theTracksLabel_ = pSet.getParameter<InputTag>("TrackCollection");
72  theTracksToken_ = consumes<reco::TrackCollection>(theTracksLabel_);
73 
74  theMaxChi2 = static_cast<unsigned int>(pSet.getParameter<double>("theMaxChi2"));
75  theNSigma = pSet.getParameter<double>("theNSigma");
76  theMinNrec = static_cast<int>(pSet.getParameter<double>("theMinNrec"));
77 
78  labelRPCRecHits = pSet.getParameter<InputTag>("theRPCRecHits");
79 
80  thedt4DSegments = pSet.getParameter<InputTag>("dt4DSegments");
81  thecscSegments = pSet.getParameter<InputTag>("cscSegments");
82 
84 
85  theMeasurementExtractor = new MuonDetLayerMeasurements(thedt4DSegments,thecscSegments,
86  labelRPCRecHits,iC,true,false,false);
87 
88  theNavigationType = pSet.getParameter<string>("NavigationType");
89 
91 }
MuonDetLayerMeasurements * theMeasurementExtractor
T getParameter(std::string const &) const
MuonServiceProxy * theService
edm::InputTag thecscSegments
edm::InputTag labelRPCRecHits
edm::EDGetTokenT< reco::TrackCollection > theTracksToken_
edm::InputTag thedt4DSegments
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
edm::InputTag theTracksLabel_
Chi2MeasurementEstimator * theEstimator
DTChamberEfficiency::~DTChamberEfficiency ( )

Definition at line 94 of file DTChamberEfficiency.cc.

References LogTrace.

95 {
96  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
97  << "DTChamberEfficiency: destructor called";
98 
99  // free memory
100  delete theService;
102  delete theEstimator;
103 }
MuonDetLayerMeasurements * theMeasurementExtractor
MuonServiceProxy * theService
#define LogTrace(id)
Chi2MeasurementEstimator * theEstimator

Member Function Documentation

void DTChamberEfficiency::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 165 of file DTChamberEfficiency.cc.

References alongMomentum, debug, DetId::det(), MuonSubdetId::DT, DTChamberId, event(), HcalObjRepresent::Fill(), TrajectoryStateOnSurface::freeState(), mergeVDriftHistosByStation::histos, i, TrajectoryStateOnSurface::isValid(), edm::HandleBase::isValid(), LogTrace, HLT_ES_cff::magfield, DetId::Muon, oppositeToMomentum, DetId::rawId(), reco::TransientTrack::recHitsSize(), DTChamberId::sector(), DTChamberId::station(), DetId::subdetId(), TrajectoryStateOnSurface::surface(), testEve_cfg::tracks, and DTChamberId::wheel().

167 {
168 
169  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") <<
170  "--- [DTChamberEfficiency] Event analysed #Run: " <<
171  event.id().run() << " #Event: " << event.id().event() << endl;
172 
173  theService->update(eventSetup);
175 
176  //Read tracks from event
178  event.getByToken(theTracksToken_, tracks);
179 
180  if(tracks.isValid()) { // check the validity of the collection
181 
182  //loop over the muons
183  for(reco::TrackCollection::const_iterator track = tracks->begin(); track!=tracks->end(); ++track) {
184 
185 
187  const int recHitsize = (int)trans_track.recHitsSize();
188  if(recHitsize < theMinNrec) continue;
189 
190  // printout the DT rechits used by the track
191  if (debug) {
192  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "--- New track" << endl;
193  set<DTChamberId> chAlrUsed;
194  for(trackingRecHit_iterator rHit = trans_track.recHitsBegin();
195  rHit != trans_track.recHitsEnd(); ++rHit) {
196  DetId rHitid = (*rHit)->geographicalId();
197  if(!(rHitid.det() == DetId::Muon && rHitid.subdetId() == MuonSubdetId::DT)) continue;
198  DTChamberId wId(rHitid.rawId());
199  if(chAlrUsed.find(wId) != chAlrUsed.end()) continue;
200  chAlrUsed.insert(wId);
201  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " " << wId << endl;
202  }
203  }
204 
205  // Get the layer on which the seed relies
206  DetId id = trans_track.track().innerDetId();
207  const DetLayer *initialLayer = theService->detLayerGeometry()->idToLayer(id);
208 
209  TrajectoryStateOnSurface init_fs = trans_track.innermostMeasurementState();
210  const FreeTrajectoryState *init_fs_free = init_fs.freeState();
211 
212  //get the list of compatible layers
213  vector<const DetLayer*> layer_list = compatibleLayers(*theService->muonNavigationSchool(), initialLayer,*init_fs_free,alongMomentum);
214  vector<const DetLayer*> layer_list_2 = compatibleLayers(*theService->muonNavigationSchool(), initialLayer,*init_fs_free,oppositeToMomentum);
215 
216  layer_list.insert(layer_list.end(),layer_list_2.begin(),layer_list_2.end());
217 
218  set<DTChamberId> alreadyCheckedCh;
219 
220  //loop over the list of compatible layers
221  for(int i=0;i< (int)layer_list.size();i++) {
222 
223  //propagate the track to the i-th layer
224  TrajectoryStateOnSurface tsos = propagator()->propagate(init_fs,layer_list.at(i)->surface());
225  if(!tsos.isValid()) continue;
226 
227  //determine the chambers kinematically compatible with the track on the i-th layer
228  vector<DetWithState> dss = layer_list.at(i)->compatibleDets(tsos, *propagator(), *theEstimator);
229 
230  if(dss.size() == 0) continue;
231 
232  // get the first det (it's the most compatible)
233  const DetWithState detWithState = dss.front();
234  const DetId idDetLay = detWithState.first->geographicalId();
235 
236  // check if this is a DT and the track has the needed quality
237  if(!chamberSelection(idDetLay,trans_track)) continue;
238 
239  DTChamberId DTid = (DTChamberId) idDetLay;
240 
241  // check if the chamber has already been counted
242  if(alreadyCheckedCh.find(DTid) != alreadyCheckedCh.end()) continue;
243  alreadyCheckedCh.insert(DTid);
244 
245  // get the compatible measurements
246  MeasurementContainer detMeasurements_initial = theMeasurementExtractor->measurements(layer_list.at(i),
247  detWithState.first,
248  detWithState.second,
249  *theEstimator, event);
250  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " chamber: " << DTid
251  << " has: " << detMeasurements_initial.size() << " comp. meas." << endl;
252 
253 
254  //we want to be more picky about the quality of the segments:
255  //exclude the segments with less than 12 hits
256  MeasurementContainer detMeasurements = segQualityCut(detMeasurements_initial);
257 
258  // get the histos for this chamber
259  vector<MonitorElement *> histos = histosPerW[DTid.wheel()+2];
260  // fill them
261  if (detMeasurements_initial.size() != 0) histos[0]->Fill(DTid.sector(),DTid.station(),1.);
262  if (detMeasurements.size() != 0) histos[1]->Fill(DTid.sector(),DTid.station(),1.);
263  histos[2]->Fill(DTid.sector(),DTid.station(),1.);
264 
265 
266  }
267  }
268  } else {
269  LogInfo("DTDQM|DTMonitorModule|DTChamberEfficiency") << "[DTChamberEfficiency] Collection: " << theTracksLabel_
270  << " is not valid!" << endl;
271  }
272  return;
273 }
MuonDetLayerMeasurements * theMeasurementExtractor
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
int i
Definition: DBlmapReader.cc:9
MuonServiceProxy * theService
MeasurementContainer segQualityCut(const MeasurementContainer &seg_list) const
std::vector< std::vector< MonitorElement * > > histosPerW
size_t recHitsSize() const
number of RecHits
edm::EDGetTokenT< reco::TrackCollection > theTracksToken_
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const SurfaceType & surface() const
bool chamberSelection(const DetId &idDetLay, reco::TransientTrack &trans_track) const
edm::ESHandle< Propagator > propagator() const
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
FreeTrajectoryState const * freeState(bool withErrors=true) const
std::vector< const DetLayer * > compatibleLayers(const NavigationSchool &navigationSchool, const DetLayer *initialLayer, const FreeTrajectoryState &fts, PropagationDirection propDir)
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:76
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
MeasurementContainer measurements(const DetLayer *layer, const GeomDet *det, const TrajectoryStateOnSurface &stateOnDet, const MeasurementEstimator &est, const edm::Event &iEvent)
edm::ESHandle< MagneticField > magfield
#define LogTrace(id)
Definition: DetId.h:18
tuple tracks
Definition: testEve_cfg.py:39
edm::InputTag theTracksLabel_
std::vector< TrajectoryMeasurement > MeasurementContainer
T const * product() const
Definition: ESHandle.h:86
void setEvent(const edm::Event &)
set event
int sector() const
Definition: DTChamberId.h:61
Chi2MeasurementEstimator * theEstimator
static const int DT
Definition: MuonSubdetId.h:12
int station() const
Return the station number.
Definition: DTChamberId.h:51
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTChamberEfficiency::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 105 of file DTChamberEfficiency.cc.

References bookHistos(), and LogTrace.

105  {
106 
107  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
108  << "DTChamberEfficiency: beginOfJob";
109 
110  bookHistos();
111 
112  return;
113 }
#define LogTrace(id)
void DTChamberEfficiency::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 115 of file DTChamberEfficiency.cc.

References edm::EventSetup::get(), and HLT_ES_cff::magfield.

116 {
117  // Get the DT Geometry
118  setup.get<MuonGeometryRecord>().get(dtGeom);
119 
120  setup.get<IdealMagneticFieldRecord>().get(magfield);
122 
123  return;
124 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
edm::ESHandle< DTGeometry > dtGeom
edm::ESHandle< MagneticField > magfield
const T & get() const
Definition: EventSetup.h:55
void DTChamberEfficiency::bookHistos ( )
private

Definition at line 135 of file DTChamberEfficiency.cc.

References mergeVDriftHistosByStation::histos, and LogTrace.

136 {
137  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
138  << "DTChamberEfficiency: booking histos";
139 
140  // Create the monitor elements
141  theDbe->setCurrentFolder("DT/05-ChamberEff/Task");
142 
143  for(int wheel=-2;wheel<=2;wheel++){
144 
145  vector<MonitorElement *> histos;
146 
147  stringstream wheel_str; wheel_str << wheel;
148 
149  histos.push_back(theDbe->book2D("hCountSectVsChamb_All_W"+ wheel_str.str(),
150  "Countings for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.));
151 
152  histos.push_back(theDbe->book2D("hCountSectVsChamb_Qual_W"+ wheel_str.str(),
153  "Countings for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.));
154 
155 
156  histos.push_back(theDbe->book2D("hExtrapSectVsChamb_W"+ wheel_str.str(),
157  "Extrapolations for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.));
158 
159  histosPerW.push_back(histos);
160  }
161 
162  return;
163 }
std::vector< std::vector< MonitorElement * > > histosPerW
#define LogTrace(id)
bool DTChamberEfficiency::chamberSelection ( const DetId idDetLay,
reco::TransientTrack trans_track 
) const
private

Definition at line 275 of file DTChamberEfficiency.cc.

References DetId::det(), MuonSubdetId::DT, DetId::Muon, reco::TransientTrack::recHit(), reco::TransientTrack::recHitsSize(), and DetId::subdetId().

276 {
277  //check that we have a muon and that is a DT detector
278  if(!(idDetLay.det() == DetId::Muon && idDetLay.subdetId() == MuonSubdetId::DT)) return false;
279 
280  if(trans_track.recHitsSize() == 2)
281  if(trans_track.recHit(0)->geographicalId() == idDetLay ||
282  trans_track.recHit(1)->geographicalId() == idDetLay) return false;
283 
284  return true;
285 }
size_t recHitsSize() const
number of RecHits
TrackingRecHitRef recHit(size_t i) const
get n-th recHit
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
static const int DT
Definition: MuonSubdetId.h:12
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
vector< const DetLayer * > DTChamberEfficiency::compatibleLayers ( const NavigationSchool navigationSchool,
const DetLayer initialLayer,
const FreeTrajectoryState fts,
PropagationDirection  propDir 
)
private

Definition at line 317 of file DTChamberEfficiency.cc.

References DirectMuonNavigation::compatibleLayers(), and NavigationSchool::compatibleLayers().

319 {
320 
321 vector<const DetLayer*> detLayers;
322 
323 if(theNavigationType == "Standard"){
324  // ask for compatible layers
325  detLayers = navigationSchool.compatibleLayers(*initialLayer,fts,propDir);
326  // I have to fit by hand the first layer until the seedTSOS is defined on the first rechit layer
327  // In fact the first layer is not returned by initialLayer->compatibleLayers.
328 
329  detLayers.insert(detLayers.begin(),initialLayer);
330 
331  }
332  else if (theNavigationType == "Direct"){
333 
334  DirectMuonNavigation navigation(ESHandle<MuonDetLayerGeometry>(&*theService->detLayerGeometry()));
335  detLayers = navigation.compatibleLayers(fts,propDir);
336  }
337  else
338  LogError("DTDQM|DTMonitorModule|DTChamberEfficiency") << "No Properly Navigation Selected!!"<<endl;
339 
340  return detLayers;
341 }
MuonServiceProxy * theService
std::vector< const DetLayer * > compatibleLayers(const DetLayer &detLayer, Args &&...args) const
Returns all layers compatible.
void DTChamberEfficiency::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 126 of file DTChamberEfficiency.cc.

References LogTrace.

127 {
128  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
129  << "DTChamberEfficiency: endOfJob";
130 
131  return;
132 }
#define LogTrace(id)
ESHandle< Propagator > DTChamberEfficiency::propagator ( void  ) const
inlineprivate

Definition at line 343 of file DTChamberEfficiency.cc.

343  {
344  return theService->propagator("SteppingHelixPropagatorAny");
345 }
MuonServiceProxy * theService
MeasurementContainer DTChamberEfficiency::segQualityCut ( const MeasurementContainer seg_list) const
private

Definition at line 289 of file DTChamberEfficiency.cc.

References DTChamberId, query::result, and DTChamberId::station().

290 {
291 
293 
294  for(MeasurementContainer::const_iterator mescont_Itr = seg_list.begin();
295  mescont_Itr != seg_list.end(); ++mescont_Itr){
296 
297  //get the rechits of the segment
298  TransientTrackingRecHit::ConstRecHitContainer recHit_list = mescont_Itr->recHit()->transientHits();
299 
300  //loop over the rechits and get the number of hits
301  int nhit_seg(0);
302  for(TransientTrackingRecHit::ConstRecHitContainer::const_iterator recList_Itr = recHit_list.begin();
303  recList_Itr != recHit_list.end(); ++recList_Itr){
304 
305  nhit_seg += (int)(*recList_Itr)->transientHits().size();
306  }
307 
308  DTChamberId tmpId = (DTChamberId)mescont_Itr->recHit()->hit()->geographicalId();
309 
310  if(tmpId.station() < 4 && nhit_seg >= 12) result.push_back(*mescont_Itr);
311  if(tmpId.station() == 4 && nhit_seg >= 8) result.push_back(*mescont_Itr);
312  }
313 
314  return result;
315 }
tuple result
Definition: query.py:137
std::vector< ConstRecHitPointer > ConstRecHitContainer
std::vector< TrajectoryMeasurement > MeasurementContainer
int station() const
Return the station number.
Definition: DTChamberId.h:51

Member Data Documentation

bool DTChamberEfficiency::debug
private
edm::ESHandle<DTGeometry> DTChamberEfficiency::dtGeom
private

Definition at line 99 of file DTChamberEfficiency.h.

std::vector<std::vector<MonitorElement*> > DTChamberEfficiency::histosPerW
private

Definition at line 112 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::labelRPCRecHits
private

Definition at line 89 of file DTChamberEfficiency.h.

edm::ESHandle<MagneticField> DTChamberEfficiency::magfield
private

Definition at line 107 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::thecscSegments
private

Definition at line 91 of file DTChamberEfficiency.h.

DQMStore* DTChamberEfficiency::theDbe
private

Definition at line 101 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::thedt4DSegments
private

Definition at line 90 of file DTChamberEfficiency.h.

Chi2MeasurementEstimator* DTChamberEfficiency::theEstimator
private

Definition at line 105 of file DTChamberEfficiency.h.

double DTChamberEfficiency::theMaxChi2
private

Definition at line 93 of file DTChamberEfficiency.h.

MuonDetLayerMeasurements* DTChamberEfficiency::theMeasurementExtractor
private

Definition at line 104 of file DTChamberEfficiency.h.

int DTChamberEfficiency::theMinNrec
private

Definition at line 95 of file DTChamberEfficiency.h.

std::string DTChamberEfficiency::theNavigationType
private

Definition at line 97 of file DTChamberEfficiency.h.

double DTChamberEfficiency::theNSigma
private

Definition at line 94 of file DTChamberEfficiency.h.

MuonServiceProxy* DTChamberEfficiency::theService
private

Definition at line 103 of file DTChamberEfficiency.h.

edm::ESHandle<GlobalTrackingGeometry> DTChamberEfficiency::theTrackingGeometry
private

Definition at line 108 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::theTracksLabel_
private

Definition at line 86 of file DTChamberEfficiency.h.

edm::EDGetTokenT<reco::TrackCollection> DTChamberEfficiency::theTracksToken_
private

Definition at line 87 of file DTChamberEfficiency.h.