CMS 3D CMS Logo

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

#include <DTChamberEfficiency.h>

Inheritance diagram for DTChamberEfficiency:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup) override
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 
 DTChamberEfficiency (const edm::ParameterSet &pset)
 
 ~DTChamberEfficiency () override
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 

Private Member Functions

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
 
edm::InputTag thedt4DSegments
 
Chi2MeasurementEstimatortheEstimator
 
double theMaxChi2
 
MuonDetLayerMeasurementstheMeasurementExtractor
 
int theMinNrec
 
std::string theNavigationType
 
double theNSigma
 
MuonServiceProxytheService
 
edm::ESHandle< GlobalTrackingGeometrytheTrackingGeometry
 
edm::InputTag theTracksLabel_
 
edm::EDGetTokenT< reco::TrackCollectiontheTracksToken_
 

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 60 of file DTChamberEfficiency.h.

Constructor & Destructor Documentation

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

Definition at line 55 of file DTChamberEfficiency.cc.

References Chi2MeasurementEstimator_cfi::Chi2MeasurementEstimator, debug, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), MuonServiceProxy_cff::MuonServiceProxy, 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  // service parameters
64  ParameterSet serviceParameters = pSet.getParameter<ParameterSet>("ServiceParameters");
65  theService = new MuonServiceProxy(serviceParameters);
66 
67  theTracksLabel_ = pSet.getParameter<InputTag>("TrackCollection");
68  theTracksToken_ = consumes<reco::TrackCollection>(theTracksLabel_);
69 
70  theMaxChi2 = static_cast<unsigned int>(pSet.getParameter<double>("theMaxChi2"));
71  theNSigma = pSet.getParameter<double>("theNSigma");
72  theMinNrec = static_cast<int>(pSet.getParameter<double>("theMinNrec"));
73 
74  labelRPCRecHits = pSet.getParameter<InputTag>("theRPCRecHits");
75 
76  thedt4DSegments = pSet.getParameter<InputTag>("dt4DSegments");
77  thecscSegments = pSet.getParameter<InputTag>("cscSegments");
78 
79  edm::ConsumesCollector iC = consumesCollector();
80 
81  theMeasurementExtractor = new MuonDetLayerMeasurements(thedt4DSegments,thecscSegments,
82  labelRPCRecHits,InputTag(),InputTag(),iC,true,false,false,false);
83 
84  theNavigationType = pSet.getParameter<string>("NavigationType");
85 
87 }
MuonDetLayerMeasurements * theMeasurementExtractor
T getParameter(std::string const &) const
MuonServiceProxy * theService
edm::InputTag thecscSegments
edm::InputTag labelRPCRecHits
edm::EDGetTokenT< reco::TrackCollection > theTracksToken_
edm::InputTag thedt4DSegments
edm::InputTag theTracksLabel_
Chi2MeasurementEstimator * theEstimator
DTChamberEfficiency::~DTChamberEfficiency ( )
override

Definition at line 90 of file DTChamberEfficiency.cc.

References LogTrace.

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

Member Function Documentation

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

Definition at line 143 of file DTChamberEfficiency.cc.

References alongMomentum, debug, DetId::det(), MuonSubdetId::DT, event(), HcalObjRepresent::Fill(), TrajectoryStateOnSurface::freeState(), plotFactory::histos, mps_fire::i, createfilelist::int, TrajectoryStateOnSurface::isValid(), edm::HandleBase::isValid(), LogTrace, DetId::Muon, oppositeToMomentum, PhotonConversionTrajectorySeedProducerFromQuadruplets_cfi::propagator, DetId::rawId(), reco::TransientTrack::recHitsSize(), DTChamberId::sector(), DTChamberId::station(), DetId::subdetId(), TrajectoryStateOnSurface::surface(), HiIsolationCommonParameters_cff::track, l1t::tracks, and DTChamberId::wheel().

145 {
146 
147  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") <<
148  "--- [DTChamberEfficiency] Event analysed #Run: " <<
149  event.id().run() << " #Event: " << event.id().event() << endl;
150 
151  theService->update(eventSetup);
153 
154  //Read tracks from event
156  event.getByToken(theTracksToken_, tracks);
157 
158  if(tracks.isValid()) { // check the validity of the collection
159 
160  //loop over the muons
161  for(reco::TrackCollection::const_iterator track = tracks->begin(); track!=tracks->end(); ++track) {
162 
163 
165  const int recHitsize = (int)trans_track.recHitsSize();
166  if(recHitsize < theMinNrec) continue;
167 
168  // printout the DT rechits used by the track
169  if (debug) {
170  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "--- New track" << endl;
171  set<DTChamberId> chAlrUsed;
172  for(trackingRecHit_iterator rHit = trans_track.recHitsBegin();
173  rHit != trans_track.recHitsEnd(); ++rHit) {
174  DetId rHitid = (*rHit)->geographicalId();
175  if(!(rHitid.det() == DetId::Muon && rHitid.subdetId() == MuonSubdetId::DT)) continue;
176  DTChamberId wId(rHitid.rawId());
177  if(chAlrUsed.find(wId) != chAlrUsed.end()) continue;
178  chAlrUsed.insert(wId);
179  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " " << wId << endl;
180  }
181  }
182 
183  // Get the layer on which the seed relies
184  DetId id = trans_track.track().innerDetId();
185  const DetLayer *initialLayer = theService->detLayerGeometry()->idToLayer(id);
186 
187  TrajectoryStateOnSurface init_fs = trans_track.innermostMeasurementState();
188  const FreeTrajectoryState *init_fs_free = init_fs.freeState();
189 
190  //get the list of compatible layers
191  vector<const DetLayer*> layer_list = compatibleLayers(*theService->muonNavigationSchool(), initialLayer,*init_fs_free,alongMomentum);
192  vector<const DetLayer*> layer_list_2 = compatibleLayers(*theService->muonNavigationSchool(), initialLayer,*init_fs_free,oppositeToMomentum);
193 
194  layer_list.insert(layer_list.end(),layer_list_2.begin(),layer_list_2.end());
195 
196  set<DTChamberId> alreadyCheckedCh;
197 
198  //loop over the list of compatible layers
199  for(int i=0;i< (int)layer_list.size();i++) {
200 
201  //propagate the track to the i-th layer
202  TrajectoryStateOnSurface tsos = propagator()->propagate(init_fs,layer_list.at(i)->surface());
203  if(!tsos.isValid()) continue;
204 
205  //determine the chambers kinematically compatible with the track on the i-th layer
206  vector<DetWithState> dss = layer_list.at(i)->compatibleDets(tsos, *propagator(), *theEstimator);
207 
208  if(dss.empty()) continue;
209 
210  // get the first det (it's the most compatible)
211  const DetWithState detWithState = dss.front();
212  const DetId idDetLay = detWithState.first->geographicalId();
213 
214  // check if this is a DT and the track has the needed quality
215  if(!chamberSelection(idDetLay,trans_track)) continue;
216 
217  DTChamberId DTid = (DTChamberId) idDetLay;
218 
219  // check if the chamber has already been counted
220  if(alreadyCheckedCh.find(DTid) != alreadyCheckedCh.end()) continue;
221  alreadyCheckedCh.insert(DTid);
222 
223  // get the compatible measurements
224  MeasurementContainer detMeasurements_initial = theMeasurementExtractor->measurements(layer_list.at(i),
225  detWithState.first,
226  detWithState.second,
227  *theEstimator, event);
228  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " chamber: " << DTid
229  << " has: " << detMeasurements_initial.size() << " comp. meas." << endl;
230 
231 
232  //we want to be more picky about the quality of the segments:
233  //exclude the segments with less than 12 hits
234  MeasurementContainer detMeasurements = segQualityCut(detMeasurements_initial);
235 
236  // get the histos for this chamber
237  vector<MonitorElement *> histos = histosPerW[DTid.wheel()+2];
238  // fill them
239  if (!detMeasurements_initial.empty()) histos[0]->Fill(DTid.sector(),DTid.station(),1.);
240  if (!detMeasurements.empty()) histos[1]->Fill(DTid.sector(),DTid.station(),1.);
241  histos[2]->Fill(DTid.sector(),DTid.station(),1.);
242 
243 
244  }
245  }
246  } else {
247  LogInfo("DTDQM|DTMonitorModule|DTChamberEfficiency") << "[DTChamberEfficiency] Collection: " << theTracksLabel_
248  << " is not valid!" << endl;
249  }
250  return;
251 }
MuonDetLayerMeasurements * theMeasurementExtractor
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
MuonServiceProxy * theService
MeasurementContainer segQualityCut(const MeasurementContainer &seg_list) const
std::vector< std::vector< MonitorElement * > > histosPerW
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
size_t recHitsSize() const
number of RecHits
edm::EDGetTokenT< reco::TrackCollection > theTracksToken_
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
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
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:74
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
edm::InputTag theTracksLabel_
std::vector< TrajectoryMeasurement > MeasurementContainer
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:53
void setEvent(const edm::Event &)
set event
int sector() const
Definition: DTChamberId.h:61
Chi2MeasurementEstimator * theEstimator
static constexpr int DT
Definition: MuonSubdetId.h:12
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
T const * product() const
Definition: ESHandle.h:86
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39
void DTChamberEfficiency::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  context 
)
overrideprotected

Definition at line 113 of file DTChamberEfficiency.cc.

References DQMStore::IBooker::book2D(), plotFactory::histos, LogTrace, DQMStore::IBooker::setCurrentFolder(), and makeMuonMisalignmentScenario::wheel.

113  {
114 
115  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
116  << "DTChamberEfficiency: booking histos";
117 
118  // Create the monitor elements
119  ibooker.setCurrentFolder("DT/05-ChamberEff/Task");
120 
121  for(int wheel=-2;wheel<=2;wheel++){
122 
123  vector<MonitorElement *> histos;
124 
125  stringstream wheel_str; wheel_str << wheel;
126 
127  histos.push_back(ibooker.book2D("hCountSectVsChamb_All_W"+ wheel_str.str(),
128  "Countings for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.));
129 
130  histos.push_back(ibooker.book2D("hCountSectVsChamb_Qual_W"+ wheel_str.str(),
131  "Countings for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.));
132 
133 
134  histos.push_back(ibooker.book2D("hExtrapSectVsChamb_W"+ wheel_str.str(),
135  "Extrapolations for wheel " + wheel_str.str(),14,1.,15.,4,1.,5.));
136 
137  histosPerW.push_back(histos);
138  }
139 
140  return;
141 }
std::vector< std::vector< MonitorElement * > > histosPerW
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
#define LogTrace(id)
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
bool DTChamberEfficiency::chamberSelection ( const DetId idDetLay,
reco::TransientTrack trans_track 
) const
private

Definition at line 253 of file DTChamberEfficiency.cc.

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

254 {
255  //check that we have a muon and that is a DT detector
256  if(!(idDetLay.det() == DetId::Muon && idDetLay.subdetId() == MuonSubdetId::DT)) return false;
257 
258  if(trans_track.recHitsSize() == 2)
259  if(trans_track.recHit(0)->geographicalId() == idDetLay ||
260  trans_track.recHit(1)->geographicalId() == idDetLay) return false;
261 
262  return true;
263 }
size_t recHitsSize() const
number of RecHits
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
TrackingRecHitRef recHit(size_t i) const
get n-th recHit
static constexpr int DT
Definition: MuonSubdetId.h:12
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39
vector< const DetLayer * > DTChamberEfficiency::compatibleLayers ( const NavigationSchool navigationSchool,
const DetLayer initialLayer,
const FreeTrajectoryState fts,
PropagationDirection  propDir 
)
private

Definition at line 293 of file DTChamberEfficiency.cc.

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

295 {
296 
297 vector<const DetLayer*> detLayers;
298 
299 if(theNavigationType == "Standard"){
300  // ask for compatible layers
301  detLayers = navigationSchool.compatibleLayers(*initialLayer,fts,propDir);
302  // I have to fit by hand the first layer until the seedTSOS is defined on the first rechit layer
303  // In fact the first layer is not returned by initialLayer->compatibleLayers.
304 
305  detLayers.insert(detLayers.begin(),initialLayer);
306 
307  }
308  else if (theNavigationType == "Direct"){
309 
310  DirectMuonNavigation navigation(ESHandle<MuonDetLayerGeometry>(&*theService->detLayerGeometry()));
311  detLayers = navigation.compatibleLayers(fts,propDir);
312  }
313  else
314  LogError("DTDQM|DTMonitorModule|DTChamberEfficiency") << "No Properly Navigation Selected!!"<<endl;
315 
316  return detLayers;
317 }
MuonServiceProxy * theService
std::vector< const DetLayer * > compatibleLayers(const DetLayer &detLayer, Args &&...args) const
Returns all layers compatible.
void DTChamberEfficiency::dqmBeginRun ( const edm::Run run,
const edm::EventSetup setup 
)
override

Definition at line 101 of file DTChamberEfficiency.cc.

References edm::EventSetup::get().

102 {
103  // Get the DT Geometry
104  setup.get<MuonGeometryRecord>().get(dtGeom);
105 
106  setup.get<IdealMagneticFieldRecord>().get(magfield);
108 
109  return;
110 }
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry
edm::ESHandle< DTGeometry > dtGeom
edm::ESHandle< MagneticField > magfield
T get() const
Definition: EventSetup.h:71
ESHandle< Propagator > DTChamberEfficiency::propagator ( ) const
inlineprivate

Definition at line 319 of file DTChamberEfficiency.cc.

319  {
320  return theService->propagator("SteppingHelixPropagatorAny");
321 }
MuonServiceProxy * theService
MeasurementContainer DTChamberEfficiency::segQualityCut ( const MeasurementContainer seg_list) const
private

Definition at line 265 of file DTChamberEfficiency.cc.

References createfilelist::int, mps_fire::result, and DTChamberId::station().

266 {
267 
269 
270  for(MeasurementContainer::const_iterator mescont_Itr = seg_list.begin();
271  mescont_Itr != seg_list.end(); ++mescont_Itr){
272 
273  //get the rechits of the segment
274  TransientTrackingRecHit::ConstRecHitContainer recHit_list = mescont_Itr->recHit()->transientHits();
275 
276  //loop over the rechits and get the number of hits
277  int nhit_seg(0);
278  for(TransientTrackingRecHit::ConstRecHitContainer::const_iterator recList_Itr = recHit_list.begin();
279  recList_Itr != recHit_list.end(); ++recList_Itr){
280 
281  nhit_seg += (int)(*recList_Itr)->transientHits().size();
282  }
283 
284  DTChamberId tmpId = (DTChamberId)mescont_Itr->recHit()->hit()->geographicalId();
285 
286  if(tmpId.station() < 4 && nhit_seg >= 12) result.push_back(*mescont_Itr);
287  if(tmpId.station() == 4 && nhit_seg >= 8) result.push_back(*mescont_Itr);
288  }
289 
290  return result;
291 }
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 104 of file DTChamberEfficiency.h.

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

Definition at line 113 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::labelRPCRecHits
private

Definition at line 94 of file DTChamberEfficiency.h.

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

Definition at line 110 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::thecscSegments
private

Definition at line 96 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::thedt4DSegments
private

Definition at line 95 of file DTChamberEfficiency.h.

Chi2MeasurementEstimator* DTChamberEfficiency::theEstimator
private

Definition at line 108 of file DTChamberEfficiency.h.

double DTChamberEfficiency::theMaxChi2
private

Definition at line 98 of file DTChamberEfficiency.h.

MuonDetLayerMeasurements* DTChamberEfficiency::theMeasurementExtractor
private

Definition at line 107 of file DTChamberEfficiency.h.

int DTChamberEfficiency::theMinNrec
private

Definition at line 100 of file DTChamberEfficiency.h.

std::string DTChamberEfficiency::theNavigationType
private

Definition at line 102 of file DTChamberEfficiency.h.

double DTChamberEfficiency::theNSigma
private

Definition at line 99 of file DTChamberEfficiency.h.

MuonServiceProxy* DTChamberEfficiency::theService
private

Definition at line 106 of file DTChamberEfficiency.h.

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

Definition at line 111 of file DTChamberEfficiency.h.

edm::InputTag DTChamberEfficiency::theTracksLabel_
private

Definition at line 91 of file DTChamberEfficiency.h.

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

Definition at line 92 of file DTChamberEfficiency.h.