CMS 3D CMS Logo

DTChamberEfficiency.cc
Go to the documentation of this file.
1 /******* \class DTEffAnalyzer *******
2  *
3  * Description:
4  *
5  * detailed description
6  *
7  * \author : Mario Pelliccioni, pellicci@cern.ch
8  * $date : 20/11/2008 16:50:57 CET $
9  *
10  * Modification:
11  *
12  *********************************/
13 
14 #include "DTChamberEfficiency.h"
15 
17 
21 
24 
27 
29 
31 
33 
38 
43 
44 #include <cmath>
45 
46 using namespace std;
47 using namespace edm;
48 
50  // Get the debug parameter for verbose output
51  debug = pSet.getUntrackedParameter<bool>("debug", false);
52 
53  LogVerbatim("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: constructor called";
54 
55  // service parameters
56  ParameterSet serviceParameters = pSet.getParameter<ParameterSet>("ServiceParameters");
57  theService = new MuonServiceProxy(serviceParameters, consumesCollector());
58 
59  theTracksLabel_ = pSet.getUntrackedParameter<InputTag>("TrackCollection");
60  theTracksToken_ = consumes<reco::TrackCollection>(theTracksLabel_);
61 
62  theMaxChi2 = static_cast<unsigned int>(pSet.getParameter<double>("theMaxChi2"));
63  theNSigma = pSet.getParameter<double>("theNSigma");
64  theMinNrec = static_cast<int>(pSet.getParameter<double>("theMinNrec"));
65 
66  labelRPCRecHits = pSet.getUntrackedParameter<InputTag>("theRPCRecHits");
67 
68  thedt4DSegments = pSet.getUntrackedParameter<InputTag>("dt4DSegments");
69  thecscSegments = pSet.getUntrackedParameter<InputTag>("cscSegments");
70 
71  edm::ConsumesCollector iC = consumesCollector();
72 
73  theMeasurementExtractor = new MuonDetLayerMeasurements(
74  thedt4DSegments, thecscSegments, labelRPCRecHits, InputTag(), InputTag(), iC, true, false, false, false);
75 
76  theNavigationType = pSet.getParameter<string>("NavigationType");
77 
78  theEstimator = new Chi2MeasurementEstimator(theMaxChi2, theNSigma);
79 }
80 
82  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: destructor called";
83 
84  // free memory
85  delete theService;
86  delete theMeasurementExtractor;
87  delete theEstimator;
88 }
89 
91 
93  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "DTChamberEfficiency: booking histos";
94 
95  // Create the monitor elements
96  ibooker.setCurrentFolder("DT/05-ChamberEff/Task");
97 
98  for (int wheel = -2; wheel <= 2; wheel++) {
99  vector<MonitorElement*> histos;
100 
101  stringstream wheel_str;
102  wheel_str << wheel;
103 
104  histos.push_back(ibooker.book2D(
105  "hCountSectVsChamb_All_W" + wheel_str.str(), "Countings for wheel " + wheel_str.str(), 14, 1., 15., 4, 1., 5.));
106 
107  histos.push_back(ibooker.book2D(
108  "hCountSectVsChamb_Qual_W" + wheel_str.str(), "Countings for wheel " + wheel_str.str(), 14, 1., 15., 4, 1., 5.));
109 
110  histos.push_back(ibooker.book2D(
111  "hExtrapSectVsChamb_W" + wheel_str.str(), "Extrapolations for wheel " + wheel_str.str(), 14, 1., 15., 4, 1., 5.));
112 
113  histosPerW.push_back(histos);
114  }
115 
116  return;
117 }
118 
120  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
121  << "--- [DTChamberEfficiency] Event analysed #Run: " << event.id().run() << " #Event: " << event.id().event()
122  << endl;
123 
124  theService->update(eventSetup);
125  theMeasurementExtractor->setEvent(event);
126 
127  //Read tracks from event
129  event.getByToken(theTracksToken_, tracks);
130 
131  if (tracks.isValid()) { // check the validity of the collection
132 
133  const edm::ESHandle<GlobalTrackingGeometry>& globalTrackingGeometry = theService->trackingGeometry();
134  const MagneticField* magneticField = theService->magneticField().product();
135 
136  //loop over the muons
137  for (reco::TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); ++track) {
138  reco::TransientTrack trans_track(*track, magneticField, globalTrackingGeometry);
139  const int recHitsize = (int)trans_track.recHitsSize();
140  if (recHitsize < theMinNrec)
141  continue;
142 
143  // printout the DT rechits used by the track
144  if (debug) {
145  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << "--- New track" << endl;
146  set<DTChamberId> chAlrUsed;
147  for (trackingRecHit_iterator rHit = trans_track.recHitsBegin(); rHit != trans_track.recHitsEnd(); ++rHit) {
148  DetId rHitid = (*rHit)->geographicalId();
149  if (!(rHitid.det() == DetId::Muon && rHitid.subdetId() == MuonSubdetId::DT))
150  continue;
151  DTChamberId wId(rHitid.rawId());
152  if (chAlrUsed.find(wId) != chAlrUsed.end())
153  continue;
154  chAlrUsed.insert(wId);
155  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency") << " " << wId << endl;
156  }
157  }
158 
159  // Get the layer on which the seed relies
160  DetId id = trans_track.track().innerDetId();
161  const DetLayer* initialLayer = theService->detLayerGeometry()->idToLayer(id);
162 
163  TrajectoryStateOnSurface init_fs = trans_track.innermostMeasurementState();
164  const FreeTrajectoryState* init_fs_free = init_fs.freeState();
165 
166  //get the list of compatible layers
167  vector<const DetLayer*> layer_list =
168  compatibleLayers(*theService->muonNavigationSchool(), initialLayer, *init_fs_free, alongMomentum);
169  vector<const DetLayer*> layer_list_2 =
170  compatibleLayers(*theService->muonNavigationSchool(), initialLayer, *init_fs_free, oppositeToMomentum);
171 
172  layer_list.insert(layer_list.end(), layer_list_2.begin(), layer_list_2.end());
173 
174  set<DTChamberId> alreadyCheckedCh;
175 
176  //loop over the list of compatible layers
177  for (int i = 0; i < (int)layer_list.size(); i++) {
178  //propagate the track to the i-th layer
179  TrajectoryStateOnSurface tsos = propagator()->propagate(init_fs, layer_list.at(i)->surface());
180  if (!tsos.isValid())
181  continue;
182 
183  //determine the chambers kinematically compatible with the track on the i-th layer
184  vector<DetWithState> dss = layer_list.at(i)->compatibleDets(tsos, *propagator(), *theEstimator);
185 
186  if (dss.empty())
187  continue;
188 
189  // get the first det (it's the most compatible)
190  const DetWithState detWithState = dss.front();
191  const DetId idDetLay = detWithState.first->geographicalId();
192 
193  // check if this is a DT and the track has the needed quality
194  if (!chamberSelection(idDetLay, trans_track))
195  continue;
196 
197  DTChamberId DTid = (DTChamberId)idDetLay;
198 
199  // check if the chamber has already been counted
200  if (alreadyCheckedCh.find(DTid) != alreadyCheckedCh.end())
201  continue;
202  alreadyCheckedCh.insert(DTid);
203 
204  // get the compatible measurements
205  MeasurementContainer detMeasurements_initial = theMeasurementExtractor->measurements(
206  layer_list.at(i), detWithState.first, detWithState.second, *theEstimator, event);
207  LogTrace("DTDQM|DTMonitorModule|DTChamberEfficiency")
208  << " chamber: " << DTid << " has: " << detMeasurements_initial.size() << " comp. meas." << endl;
209 
210  //we want to be more picky about the quality of the segments:
211  //exclude the segments with less than 12 hits
212  MeasurementContainer detMeasurements = segQualityCut(detMeasurements_initial);
213 
214  // get the histos for this chamber
215  vector<MonitorElement*> histos = histosPerW[DTid.wheel() + 2];
216  // fill them
217  if (!detMeasurements_initial.empty())
218  histos[0]->Fill(DTid.sector(), DTid.station(), 1.);
219  if (!detMeasurements.empty())
220  histos[1]->Fill(DTid.sector(), DTid.station(), 1.);
221  histos[2]->Fill(DTid.sector(), DTid.station(), 1.);
222  }
223  }
224  } else {
225  LogInfo("DTDQM|DTMonitorModule|DTChamberEfficiency")
226  << "[DTChamberEfficiency] Collection: " << theTracksLabel_ << " is not valid!" << endl;
227  }
228  return;
229 }
230 
231 bool DTChamberEfficiency::chamberSelection(const DetId& idDetLay, reco::TransientTrack& trans_track) const {
232  //check that we have a muon and that is a DT detector
233  if (!(idDetLay.det() == DetId::Muon && idDetLay.subdetId() == MuonSubdetId::DT))
234  return false;
235 
236  if (trans_track.recHitsSize() == 2)
237  if (trans_track.recHit(0)->geographicalId() == idDetLay || trans_track.recHit(1)->geographicalId() == idDetLay)
238  return false;
239 
240  return true;
241 }
242 
245 
246  for (MeasurementContainer::const_iterator mescont_Itr = seg_list.begin(); mescont_Itr != seg_list.end();
247  ++mescont_Itr) {
248  //get the rechits of the segment
249  TransientTrackingRecHit::ConstRecHitContainer recHit_list = mescont_Itr->recHit()->transientHits();
250 
251  //loop over the rechits and get the number of hits
252  int nhit_seg(0);
253  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator recList_Itr = recHit_list.begin();
254  recList_Itr != recHit_list.end();
255  ++recList_Itr) {
256  nhit_seg += (int)(*recList_Itr)->transientHits().size();
257  }
258 
259  DTChamberId tmpId = (DTChamberId)mescont_Itr->recHit()->hit()->geographicalId();
260 
261  if (tmpId.station() < 4 && nhit_seg >= 12)
262  result.push_back(*mescont_Itr);
263  if (tmpId.station() == 4 && nhit_seg >= 8)
264  result.push_back(*mescont_Itr);
265  }
266 
267  return result;
268 }
269 
270 vector<const DetLayer*> DTChamberEfficiency::compatibleLayers(const NavigationSchool& navigationSchool,
271  const DetLayer* initialLayer,
272  const FreeTrajectoryState& fts,
273  PropagationDirection propDir) {
274  vector<const DetLayer*> detLayers;
275 
276  if (theNavigationType == "Standard") {
277  // ask for compatible layers
278  detLayers = navigationSchool.compatibleLayers(*initialLayer, fts, propDir);
279  // I have to fit by hand the first layer until the seedTSOS is defined on the first rechit layer
280  // In fact the first layer is not returned by initialLayer->compatibleLayers.
281 
282  detLayers.insert(detLayers.begin(), initialLayer);
283 
284  } else if (theNavigationType == "Direct") {
285  DirectMuonNavigation navigation(ESHandle<MuonDetLayerGeometry>(&*theService->detLayerGeometry()));
286  detLayers = navigation.compatibleLayers(fts, propDir);
287  } else
288  LogError("DTDQM|DTMonitorModule|DTChamberEfficiency") << "No Properly Navigation Selected!!" << endl;
289 
290  return detLayers;
291 }
292 
294  return theService->propagator("SteppingHelixPropagatorAny");
295 }
296 
297 // Local Variables:
298 // show-trailing-whitespace: t
299 // truncate-lines: t
300 // End:
Log< level::Info, true > LogVerbatim
int station() const
Return the station number.
Definition: DTChamberId.h:42
const Track & track() const
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< const DetLayer * > compatibleLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
trackingRecHit_iterator recHitsEnd() const
last iterator to RecHits
DTChamberEfficiency(const edm::ParameterSet &pset)
bool chamberSelection(const DetId &idDetLay, reco::TransientTrack &trans_track) const
MeasurementContainer segQualityCut(const MeasurementContainer &seg_list) const
PropagationDirection
trackingRecHit_iterator recHitsBegin() const
first iterator to RecHits
Log< level::Error, false > LogError
const SurfaceType & surface() const
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
#define LogTrace(id)
T getUntrackedParameter(std::string const &, T const &) const
edm::ESHandle< Propagator > propagator() const
std::vector< const DetLayer * > compatibleLayers(const DetLayer &detLayer, Args &&... args) const
Returns all layers compatible.
std::pair< const GeomDet *, TrajectoryStateOnSurface > DetWithState
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
TrackingRecHitRef recHit(size_t i) const
get n-th recHit
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:48
TrajectoryStateOnSurface innermostMeasurementState() const
std::vector< const DetLayer * > compatibleLayers(const NavigationSchool &navigationSchool, const DetLayer *initialLayer, const FreeTrajectoryState &fts, PropagationDirection propDir)
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:82
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
std::vector< ConstRecHitPointer > ConstRecHitContainer
Log< level::Info, false > LogInfo
Definition: DetId.h:17
#define debug
Definition: HDRShower.cc:19
auto const & tracks
cannot be loose
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
std::vector< TrajectoryMeasurement > MeasurementContainer
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:212
histos
Definition: combine.py:4
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
HLT enums.
int sector() const
Definition: DTChamberId.h:49
FreeTrajectoryState const * freeState(bool withErrors=true) const
static constexpr int DT
Definition: MuonSubdetId.h:11
Definition: event.py:1
Definition: Run.h:45