CMS 3D CMS Logo

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

#include <DQM/TrackingMonitor/src/TrackSplittingMonitor.cc>

Inheritance diagram for TrackSplittingMonitor:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 TrackSplittingMonitor (const edm::ParameterSet &)
 
 ~TrackSplittingMonitor () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

void doProfileX (TH2 *th2, MonitorElement *me)
 
void doProfileX (MonitorElement *th2m, MonitorElement *me)
 

Private Attributes

edm::ParameterSet conf_
 
const CSCGeometrycscGeometry
 
const edm::ESGetToken< CSCGeometry, MuonGeometryRecordcscGeomToken_
 
double d0Cut_
 
MonitorElementdcurvAbsoluteResiduals_global_
 
MonitorElementdcurvAbsoluteResiduals_tracker_
 
MonitorElementdcurvNormalizedResiduals_global_
 
MonitorElementdcurvNormalizedResiduals_tracker_
 
MonitorElementddxyAbsoluteResiduals_global_
 
MonitorElementddxyAbsoluteResiduals_tracker_
 
MonitorElementddxyNormalizedResiduals_global_
 
MonitorElementddxyNormalizedResiduals_tracker_
 
MonitorElementddzAbsoluteResiduals_global_
 
MonitorElementddzAbsoluteResiduals_tracker_
 
MonitorElementddzNormalizedResiduals_global_
 
MonitorElementddzNormalizedResiduals_tracker_
 
MonitorElementdphiAbsoluteResiduals_global_
 
MonitorElementdphiAbsoluteResiduals_tracker_
 
MonitorElementdphiNormalizedResiduals_global_
 
MonitorElementdphiNormalizedResiduals_tracker_
 
MonitorElementdptAbsoluteResiduals_global_
 
MonitorElementdptAbsoluteResiduals_tracker_
 
MonitorElementdptNormalizedResiduals_global_
 
MonitorElementdptNormalizedResiduals_tracker_
 
DQMStoredqmStore_
 
const DTGeometrydtGeometry
 
const edm::ESGetToken< DTGeometry, MuonGeometryRecorddtGeomToken_
 
MonitorElementdthetaAbsoluteResiduals_global_
 
MonitorElementdthetaAbsoluteResiduals_tracker_
 
MonitorElementdthetaNormalizedResiduals_global_
 
MonitorElementdthetaNormalizedResiduals_tracker_
 
double dzCut_
 
std::string histname
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmfToken_
 
double norchiCut_
 
int pixelHitsPerLeg_
 
bool plotMuons_
 
double ptCut_
 
const RPCGeometryrpcGeometry
 
const edm::ESGetToken< RPCGeometry, MuonGeometryRecordrpcGeomToken_
 
edm::InputTag splitMuons_
 
edm::EDGetTokenT< std::vector< reco::Muon > > splitMuonsToken_
 
edm::InputTag splitTracks_
 
edm::EDGetTokenT< std::vector< reco::Track > > splitTracksToken_
 
const TrackerGeometrytheGeometry
 
const MagneticFieldtheMagField
 
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtkGeomToken_
 
int totalHitsPerLeg_
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

Monitoring source for general quantities related to tracks.

Definition at line 40 of file TrackSplittingMonitor.h.

Constructor & Destructor Documentation

◆ TrackSplittingMonitor()

TrackSplittingMonitor::TrackSplittingMonitor ( const edm::ParameterSet iConfig)
explicit

Definition at line 26 of file TrackSplittingMonitor.cc.

References conf_, d0Cut_, dzCut_, edm::ParameterSet::getParameter(), norchiCut_, pixelHitsPerLeg_, plotMuons_, ptCut_, splitMuons_, splitMuonsToken_, splitTracks_, splitTracksToken_, and totalHitsPerLeg_.

28  conf_(iConfig),
34  splitTracks_ = conf_.getParameter<edm::InputTag>("splitTrackCollection");
35  splitMuons_ = conf_.getParameter<edm::InputTag>("splitMuonCollection");
36  splitTracksToken_ = consumes<std::vector<reco::Track> >(splitTracks_);
37  splitMuonsToken_ = mayConsume<std::vector<reco::Muon> >(splitMuons_);
38 
39  plotMuons_ = conf_.getParameter<bool>("ifPlotMuons");
40 
41  // cuts
42  pixelHitsPerLeg_ = conf_.getParameter<int>("pixelHitsPerLeg");
43  totalHitsPerLeg_ = conf_.getParameter<int>("totalHitsPerLeg");
44  d0Cut_ = conf_.getParameter<double>("d0Cut");
45  dzCut_ = conf_.getParameter<double>("dzCut");
46  ptCut_ = conf_.getParameter<double>("ptCut");
47  norchiCut_ = conf_.getParameter<double>("norchiCut");
48 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > cscGeomToken_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
edm::EDGetTokenT< std::vector< reco::Muon > > splitMuonsToken_
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > mfToken_
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
edm::EDGetTokenT< std::vector< reco::Track > > splitTracksToken_
const edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeomToken_

◆ ~TrackSplittingMonitor()

TrackSplittingMonitor::~TrackSplittingMonitor ( )
overridedefault

Member Function Documentation

◆ analyze()

void TrackSplittingMonitor::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 178 of file TrackSplittingMonitor.cc.

References cscGeometry, cscGeomToken_, reco::TrackBase::d0(), d0Cut_, reco::TrackBase::d0Error(), ddxyAbsoluteResiduals_global_, ddxyAbsoluteResiduals_tracker_, ddxyNormalizedResiduals_global_, ddxyNormalizedResiduals_tracker_, dtGeometry, dtGeomToken_, reco::TrackBase::dz(), dzCut_, reco::TrackBase::dzError(), dqm::impl::MonitorElement::Fill(), edm::EventSetup::getData(), reco::Muon::globalTrack(), iEvent, edm::HandleBase::isValid(), mfToken_, norchiCut_, reco::TrackBase::normalizedChi2(), reco::TrackBase::phi(), reco::TrackBase::phiError(), PixelSubdetector::PixelBarrel, pixelHitsPerLeg_, plotMuons_, funct::pow(), reco::TrackBase::pt(), HLT_2022v12_cff::pt1, HLT_2022v12_cff::pt2, ptCut_, reco::TrackBase::ptError(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), rpcGeometry, rpcGeomToken_, RecoMuonCosmics_cff::splitMuons, splitMuonsToken_, splitTracksToken_, mathSSE::sqrt(), theGeometry, theMagField, reco::TrackBase::theta(), reco::TrackBase::thetaError(), tkGeomToken_, and totalHitsPerLeg_.

178  {
179  theMagField = &iSetup.getData(mfToken_);
180  theGeometry = &iSetup.getData(tkGeomToken_);
181  dtGeometry = &iSetup.getData(dtGeomToken_);
182  cscGeometry = &iSetup.getData(cscGeomToken_);
183  rpcGeometry = &iSetup.getData(rpcGeomToken_);
184 
186  if (!splitTracks.isValid())
187  return;
188 
190  if (plotMuons_) {
191  splitMuons = iEvent.getHandle(splitMuonsToken_);
192  }
193 
194  if (splitTracks->size() == 2) {
195  // check that there are 2 tracks in split track collection
196  edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size();
197 
198  // split tracks calculations
199  reco::Track track1 = splitTracks->at(0);
200  reco::Track track2 = splitTracks->at(1);
201 
202  // -------------------------- basic selection ---------------------------
203 
204  // hit counting
205  // looping through the hits for track 1
206  double nRechits1 = 0;
207  double nRechitinBPIX1 = 0;
208  for (trackingRecHit_iterator iHit = track1.recHitsBegin(); iHit != track1.recHitsEnd(); ++iHit) {
209  if ((*iHit)->isValid()) {
210  nRechits1++;
211  int type = (*iHit)->geographicalId().subdetId();
212  if (type == int(PixelSubdetector::PixelBarrel)) {
213  ++nRechitinBPIX1;
214  }
215  }
216  }
217  // looping through the hits for track 2
218  double nRechits2 = 0;
219  double nRechitinBPIX2 = 0;
220  for (trackingRecHit_iterator iHit = track2.recHitsBegin(); iHit != track2.recHitsEnd(); ++iHit) {
221  if ((*iHit)->isValid()) {
222  nRechits2++;
223  int type = (*iHit)->geographicalId().subdetId();
224  if (type == int(PixelSubdetector::PixelBarrel)) {
225  ++nRechitinBPIX2;
226  }
227  }
228  }
229 
230  // DCA of each track
231  double d01 = track1.d0();
232  double dz1 = track1.dz();
233  double d02 = track2.d0();
234  double dz2 = track2.dz();
235 
236  // pT of each track
237  double pt1 = track1.pt();
238  double pt2 = track2.pt();
239 
240  // chi2 of each track
241  double norchi1 = track1.normalizedChi2();
242  double norchi2 = track2.normalizedChi2();
243 
244  // basic selection
245  // pixel hits and total hits
246  if ((nRechitinBPIX1 >= pixelHitsPerLeg_) && (nRechitinBPIX1 >= pixelHitsPerLeg_) &&
247  (nRechits1 >= totalHitsPerLeg_) && (nRechits2 >= totalHitsPerLeg_)) {
248  // dca cut
249  if (((fabs(d01) < d0Cut_)) && (fabs(d02) < d0Cut_) && (fabs(dz2) < dzCut_) && (fabs(dz2) < dzCut_)) {
250  // pt cut
251  if ((pt1 + pt2) / 2 < ptCut_) {
252  // chi2 cut
253  if ((norchi1 < norchiCut_) && (norchi2 < norchiCut_)) {
254  // passed all cuts...
255  edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?";
256 
257  double ddxyVal = d01 - d02;
258  double ddzVal = dz1 - dz2;
259  double dphiVal = track1.phi() - track2.phi();
260  double dthetaVal = track1.theta() - track2.theta();
261  double dptVal = pt1 - pt2;
262  double dcurvVal = (1 / pt1) - (1 / pt2);
263 
264  double d01ErrVal = track1.d0Error();
265  double d02ErrVal = track2.d0Error();
266  double dz1ErrVal = track1.dzError();
267  double dz2ErrVal = track2.dzError();
268  double phi1ErrVal = track1.phiError();
269  double phi2ErrVal = track2.phiError();
270  double theta1ErrVal = track1.thetaError();
271  double theta2ErrVal = track2.thetaError();
272  double pt1ErrVal = track1.ptError();
273  double pt2ErrVal = track2.ptError();
274 
275  ddxyAbsoluteResiduals_tracker_->Fill(10000.0 * ddxyVal / sqrt(2.0));
276  ddxyAbsoluteResiduals_tracker_->Fill(10000.0 * ddzVal / sqrt(2.0));
277  ddxyAbsoluteResiduals_tracker_->Fill(1000.0 * dphiVal / sqrt(2.0));
278  ddxyAbsoluteResiduals_tracker_->Fill(1000.0 * dthetaVal / sqrt(2.0));
279  ddxyAbsoluteResiduals_tracker_->Fill(dptVal / sqrt(2.0));
280  ddxyAbsoluteResiduals_tracker_->Fill(dcurvVal / sqrt(2.0));
281 
282  ddxyNormalizedResiduals_tracker_->Fill(ddxyVal / sqrt(d01ErrVal * d01ErrVal + d02ErrVal * d02ErrVal));
283  ddxyNormalizedResiduals_tracker_->Fill(ddzVal / sqrt(dz1ErrVal * dz1ErrVal + dz2ErrVal * dz2ErrVal));
284  ddxyNormalizedResiduals_tracker_->Fill(dphiVal / sqrt(phi1ErrVal * phi1ErrVal + phi2ErrVal * phi2ErrVal));
286  sqrt(theta1ErrVal * theta1ErrVal + theta2ErrVal * theta2ErrVal));
287  ddxyNormalizedResiduals_tracker_->Fill(dptVal / sqrt(pt1ErrVal * pt1ErrVal + pt2ErrVal * pt2ErrVal));
289  dcurvVal / sqrt(pow(pt1ErrVal, 2) / pow(pt1, 4) + pow(pt2ErrVal, 2) / pow(pt2, 4)));
290 
291  // if do the same for split muons
292  if (plotMuons_ && splitMuons.isValid()) {
293  int gmCtr = 0;
294  bool topGlobalMuonFlag = false;
295  bool bottomGlobalMuonFlag = false;
296  int topGlobalMuon = -1;
297  int bottomGlobalMuon = -1;
298  double topGlobalMuonNorchi2 = 1e10;
299  double bottomGlobalMuonNorchi2 = 1e10;
300 
301  // check if usable split global muons
302  for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++) {
303  if (gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon()) {
304  reco::TrackRef trackerTrackRef1(splitTracks, 0);
305  reco::TrackRef trackerTrackRef2(splitTracks, 1);
306 
307  if (gmI->innerTrack() == trackerTrackRef1) {
308  if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2) {
309  topGlobalMuonFlag = true;
310  topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
311  topGlobalMuon = gmCtr;
312  }
313  }
314  if (gmI->innerTrack() == trackerTrackRef2) {
315  if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2) {
316  bottomGlobalMuonFlag = true;
317  bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
318  bottomGlobalMuon = gmCtr;
319  }
320  }
321  }
322  gmCtr++;
323  }
324 
325  if (bottomGlobalMuonFlag && topGlobalMuonFlag) {
326  reco::Muon muonTop = splitMuons->at(topGlobalMuon);
327  reco::Muon muonBottom = splitMuons->at(bottomGlobalMuon);
328 
329  reco::TrackRef glb1 = muonTop.globalTrack();
330  reco::TrackRef glb2 = muonBottom.globalTrack();
331 
332  double ddxyValGlb = glb1->d0() - glb2->d0();
333  double ddzValGlb = glb1->dz() - glb2->dz();
334  double dphiValGlb = glb1->phi() - glb2->phi();
335  double dthetaValGlb = glb1->theta() - glb2->theta();
336  double dptValGlb = glb1->pt() - glb2->pt();
337  double dcurvValGlb = (1 / glb1->pt()) - (1 / glb2->pt());
338 
339  double d01ErrValGlb = glb1->d0Error();
340  double d02ErrValGlb = glb2->d0Error();
341  double dz1ErrValGlb = glb1->dzError();
342  double dz2ErrValGlb = glb2->dzError();
343  double phi1ErrValGlb = glb1->phiError();
344  double phi2ErrValGlb = glb2->phiError();
345  double theta1ErrValGlb = glb1->thetaError();
346  double theta2ErrValGlb = glb2->thetaError();
347  double pt1ErrValGlb = glb1->ptError();
348  double pt2ErrValGlb = glb2->ptError();
349 
350  ddxyAbsoluteResiduals_global_->Fill(10000.0 * ddxyValGlb / sqrt(2.0));
351  ddxyAbsoluteResiduals_global_->Fill(10000.0 * ddzValGlb / sqrt(2.0));
352  ddxyAbsoluteResiduals_global_->Fill(1000.0 * dphiValGlb / sqrt(2.0));
353  ddxyAbsoluteResiduals_global_->Fill(1000.0 * dthetaValGlb / sqrt(2.0));
354  ddxyAbsoluteResiduals_global_->Fill(dptValGlb / sqrt(2.0));
355  ddxyAbsoluteResiduals_global_->Fill(dcurvValGlb / sqrt(2.0));
356 
358  sqrt(d01ErrValGlb * d01ErrValGlb + d02ErrValGlb * d02ErrValGlb));
360  sqrt(dz1ErrValGlb * dz1ErrValGlb + dz2ErrValGlb * dz2ErrValGlb));
362  dphiValGlb / sqrt(phi1ErrValGlb * phi1ErrValGlb + phi2ErrValGlb * phi2ErrValGlb));
364  dthetaValGlb / sqrt(theta1ErrValGlb * theta1ErrValGlb + theta2ErrValGlb * theta2ErrValGlb));
366  sqrt(pt1ErrValGlb * pt1ErrValGlb + pt2ErrValGlb * pt2ErrValGlb));
368  dcurvValGlb / sqrt(pow(pt1ErrValGlb, 2) / pow(pt1, 4) + pow(pt2ErrValGlb, 2) / pow(pt2, 4)));
369  }
370 
371  } // end of split muons loop
372  }
373  }
374  }
375  }
376  }
377 }
const CSCGeometry * cscGeometry
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:51
double thetaError() const
error on theta
Definition: TrackBase.h:757
MonitorElement * ddxyNormalizedResiduals_tracker_
const DTGeometry * dtGeometry
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > cscGeomToken_
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
edm::EDGetTokenT< std::vector< reco::Muon > > splitMuonsToken_
void Fill(long long x)
const RPCGeometry * rpcGeometry
double pt() const
track transverse momentum
Definition: TrackBase.h:637
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > mfToken_
MonitorElement * ddxyNormalizedResiduals_global_
double dzError() const
error on dz
Definition: TrackBase.h:778
const MagneticField * theMagField
const TrackerGeometry * theGeometry
T sqrt(T t)
Definition: SSEVec.h:19
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:91
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:88
Log< level::Info, false > LogInfo
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:611
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
bool isValid() const
Definition: HandleBase.h:70
double theta() const
polar angle
Definition: TrackBase.h:602
MonitorElement * ddxyAbsoluteResiduals_global_
edm::EDGetTokenT< std::vector< reco::Track > > splitTracksToken_
double d0Error() const
error on d0
Definition: TrackBase.h:772
double phiError() const
error on phi
Definition: TrackBase.h:766
MonitorElement * ddxyAbsoluteResiduals_tracker_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
const edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeomToken_

◆ bookHistograms()

void TrackSplittingMonitor::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Implements DQMEDAnalyzer.

Definition at line 52 of file TrackSplittingMonitor.cc.

References dqm::implementation::IBooker::book1D(), conf_, dcurvAbsoluteResiduals_global_, dcurvAbsoluteResiduals_tracker_, TrackSplittingMonitor_cfi::dcurvBin, TrackSplittingMonitor_cfi::dcurvMax, TrackSplittingMonitor_cfi::dcurvMin, dcurvNormalizedResiduals_global_, dcurvNormalizedResiduals_tracker_, ddxyAbsoluteResiduals_global_, ddxyAbsoluteResiduals_tracker_, TrackSplittingMonitor_cfi::ddxyBin, TrackSplittingMonitor_cfi::ddxyMax, TrackSplittingMonitor_cfi::ddxyMin, ddxyNormalizedResiduals_global_, ddxyNormalizedResiduals_tracker_, ddzAbsoluteResiduals_global_, ddzAbsoluteResiduals_tracker_, TrackSplittingMonitor_cfi::ddzBin, TrackSplittingMonitor_cfi::ddzMax, TrackSplittingMonitor_cfi::ddzMin, ddzNormalizedResiduals_global_, ddzNormalizedResiduals_tracker_, dphiAbsoluteResiduals_global_, dphiAbsoluteResiduals_tracker_, TrackSplittingMonitor_cfi::dphiBin, TrackSplittingMonitor_cfi::dphiMax, TrackSplittingMonitor_cfi::dphiMin, dphiNormalizedResiduals_global_, dphiNormalizedResiduals_tracker_, dptAbsoluteResiduals_global_, dptAbsoluteResiduals_tracker_, TrackSplittingMonitor_cfi::dptBin, TrackSplittingMonitor_cfi::dptMax, TrackSplittingMonitor_cfi::dptMin, dptNormalizedResiduals_global_, dptNormalizedResiduals_tracker_, dthetaAbsoluteResiduals_global_, dthetaAbsoluteResiduals_tracker_, TrackSplittingMonitor_cfi::dthetaBin, TrackSplittingMonitor_cfi::dthetaMax, TrackSplittingMonitor_cfi::dthetaMin, dthetaNormalizedResiduals_global_, dthetaNormalizedResiduals_tracker_, edm::ParameterSet::getParameter(), TrackSplittingMonitor_cfi::normBin, TrackSplittingMonitor_cfi::normMax, TrackSplittingMonitor_cfi::normMin, plotMuons_, dqm::impl::MonitorElement::setAxisTitle(), dqm::implementation::NavigatorBase::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

56 {
57  std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
58  ibooker.setCurrentFolder(MEFolderName);
59 
60  // bin declarations
61  int ddxyBin = conf_.getParameter<int>("ddxyBin");
62  double ddxyMin = conf_.getParameter<double>("ddxyMin");
63  double ddxyMax = conf_.getParameter<double>("ddxyMax");
64 
65  int ddzBin = conf_.getParameter<int>("ddzBin");
66  double ddzMin = conf_.getParameter<double>("ddzMin");
67  double ddzMax = conf_.getParameter<double>("ddzMax");
68 
69  int dphiBin = conf_.getParameter<int>("dphiBin");
70  double dphiMin = conf_.getParameter<double>("dphiMin");
71  double dphiMax = conf_.getParameter<double>("dphiMax");
72 
73  int dthetaBin = conf_.getParameter<int>("dthetaBin");
74  double dthetaMin = conf_.getParameter<double>("dthetaMin");
75  double dthetaMax = conf_.getParameter<double>("dthetaMax");
76 
77  int dptBin = conf_.getParameter<int>("dptBin");
78  double dptMin = conf_.getParameter<double>("dptMin");
79  double dptMax = conf_.getParameter<double>("dptMax");
80 
81  int dcurvBin = conf_.getParameter<int>("dcurvBin");
82  double dcurvMin = conf_.getParameter<double>("dcurvMin");
83  double dcurvMax = conf_.getParameter<double>("dcurvMax");
84 
85  int normBin = conf_.getParameter<int>("normBin");
86  double normMin = conf_.getParameter<double>("normMin");
87  double normMax = conf_.getParameter<double>("normMax");
88 
89  // declare histogram
91  ibooker.book1D("ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax);
93  ibooker.book1D("ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax);
95  ibooker.book1D("dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax);
97  "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax);
99  ibooker.book1D("dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax);
101  ibooker.book1D("dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax);
102 
104  ibooker.book1D("ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax);
106  ibooker.book1D("ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax);
108  ibooker.book1D("dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax);
110  "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax);
112  ibooker.book1D("dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax);
114  ibooker.book1D("dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax);
115 
116  if (plotMuons_) {
118  ibooker.book1D("ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax);
120  ibooker.book1D("ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax);
122  ibooker.book1D("dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax);
124  "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax);
126  ibooker.book1D("dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax);
128  ibooker.book1D("dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax);
129 
131  ibooker.book1D("ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax);
133  ibooker.book1D("ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax);
135  ibooker.book1D("dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax);
137  "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax);
139  ibooker.book1D("dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax);
141  ibooker.book1D("dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax);
142  }
143 
144  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta d_{xy})/#sqrt{2} [#mum]");
145  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta d_{z})/#sqrt{2} [#mum]");
146  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta #phi)/#sqrt{2} [mrad]");
147  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta #theta)/#sqrt{2} [mrad]");
148  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta pT)/#sqrt{2} [GeV]");
149  ddxyAbsoluteResiduals_tracker_->setAxisTitle("(#delta (1/pT))/#sqrt{2} [GeV^{-1}]");
150 
151  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta d_{xy}/#sigma(d_{xy}");
152  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta d_{z}/#sigma(d_{z})");
153  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta #phi/#sigma(d_{#phi})");
154  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta #theta/#sigma(d_{#theta})");
155  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta p_{T}/#sigma(p_{T})");
156  ddxyNormalizedResiduals_tracker_->setAxisTitle("#delta 1/p_{T}/#sigma(1/p_{T})");
157 
158  if (plotMuons_) {
159  ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta d_{xy})/#sqrt{2} [#mum]");
160  ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta d_{z})/#sqrt{2} [#mum]");
161  ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta #phi)/#sqrt{2} [mrad]");
162  ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta #theta)/#sqrt{2} [mrad]");
163  ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta pT)/#sqrt{2} [GeV]");
164  ddxyAbsoluteResiduals_global_->setAxisTitle("(#delta (1/pT))/#sqrt{2} [GeV^{-1}]");
165 
166  ddxyNormalizedResiduals_global_->setAxisTitle("#delta d_{xy}/#sigma(d_{xy}");
167  ddxyNormalizedResiduals_global_->setAxisTitle("#delta d_{z}/#sigma(d_{z})");
168  ddxyNormalizedResiduals_global_->setAxisTitle("#delta #phi/#sigma(d_{#phi})");
169  ddxyNormalizedResiduals_global_->setAxisTitle("#delta #theta/#sigma(d_{#theta})");
170  ddxyNormalizedResiduals_global_->setAxisTitle("#delta p_{T}/#sigma(p_{T})");
171  ddxyNormalizedResiduals_global_->setAxisTitle("#delta 1/p_{T}/#sigma(1/p_{T})");
172  }
173 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
MonitorElement * dthetaAbsoluteResiduals_global_
MonitorElement * dcurvNormalizedResiduals_global_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
MonitorElement * dthetaAbsoluteResiduals_tracker_
MonitorElement * ddzNormalizedResiduals_global_
MonitorElement * ddxyNormalizedResiduals_tracker_
MonitorElement * dcurvAbsoluteResiduals_tracker_
MonitorElement * ddzAbsoluteResiduals_global_
MonitorElement * dphiNormalizedResiduals_global_
MonitorElement * dthetaNormalizedResiduals_tracker_
MonitorElement * ddxyNormalizedResiduals_global_
MonitorElement * dthetaNormalizedResiduals_global_
MonitorElement * dphiNormalizedResiduals_tracker_
MonitorElement * dcurvAbsoluteResiduals_global_
MonitorElement * dphiAbsoluteResiduals_tracker_
MonitorElement * ddzNormalizedResiduals_tracker_
MonitorElement * dptAbsoluteResiduals_global_
MonitorElement * dptAbsoluteResiduals_tracker_
MonitorElement * ddzAbsoluteResiduals_tracker_
MonitorElement * ddxyAbsoluteResiduals_global_
MonitorElement * dcurvNormalizedResiduals_tracker_
MonitorElement * dphiAbsoluteResiduals_global_
MonitorElement * dptNormalizedResiduals_tracker_
MonitorElement * dptNormalizedResiduals_global_
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
MonitorElement * ddxyAbsoluteResiduals_tracker_
virtual void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)

◆ doProfileX() [1/2]

void TrackSplittingMonitor::doProfileX ( TH2 *  th2,
MonitorElement me 
)
private

◆ doProfileX() [2/2]

void TrackSplittingMonitor::doProfileX ( MonitorElement th2m,
MonitorElement me 
)
private

Member Data Documentation

◆ conf_

edm::ParameterSet TrackSplittingMonitor::conf_
private

Definition at line 59 of file TrackSplittingMonitor.h.

Referenced by bookHistograms(), and TrackSplittingMonitor().

◆ cscGeometry

const CSCGeometry* TrackSplittingMonitor::cscGeometry
private

Definition at line 70 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ cscGeomToken_

const edm::ESGetToken<CSCGeometry, MuonGeometryRecord> TrackSplittingMonitor::cscGeomToken_
private

Definition at line 64 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ d0Cut_

double TrackSplittingMonitor::d0Cut_
private

Definition at line 81 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().

◆ dcurvAbsoluteResiduals_global_

MonitorElement* TrackSplittingMonitor::dcurvAbsoluteResiduals_global_
private

Definition at line 106 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dcurvAbsoluteResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dcurvAbsoluteResiduals_tracker_
private

Definition at line 92 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dcurvNormalizedResiduals_global_

MonitorElement* TrackSplittingMonitor::dcurvNormalizedResiduals_global_
private

Definition at line 113 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dcurvNormalizedResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dcurvNormalizedResiduals_tracker_
private

Definition at line 99 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ ddxyAbsoluteResiduals_global_

MonitorElement* TrackSplittingMonitor::ddxyAbsoluteResiduals_global_
private

Definition at line 101 of file TrackSplittingMonitor.h.

Referenced by analyze(), and bookHistograms().

◆ ddxyAbsoluteResiduals_tracker_

MonitorElement* TrackSplittingMonitor::ddxyAbsoluteResiduals_tracker_
private

Definition at line 87 of file TrackSplittingMonitor.h.

Referenced by analyze(), and bookHistograms().

◆ ddxyNormalizedResiduals_global_

MonitorElement* TrackSplittingMonitor::ddxyNormalizedResiduals_global_
private

Definition at line 108 of file TrackSplittingMonitor.h.

Referenced by analyze(), and bookHistograms().

◆ ddxyNormalizedResiduals_tracker_

MonitorElement* TrackSplittingMonitor::ddxyNormalizedResiduals_tracker_
private

Definition at line 94 of file TrackSplittingMonitor.h.

Referenced by analyze(), and bookHistograms().

◆ ddzAbsoluteResiduals_global_

MonitorElement* TrackSplittingMonitor::ddzAbsoluteResiduals_global_
private

Definition at line 102 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ ddzAbsoluteResiduals_tracker_

MonitorElement* TrackSplittingMonitor::ddzAbsoluteResiduals_tracker_
private

Definition at line 88 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ ddzNormalizedResiduals_global_

MonitorElement* TrackSplittingMonitor::ddzNormalizedResiduals_global_
private

Definition at line 109 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ ddzNormalizedResiduals_tracker_

MonitorElement* TrackSplittingMonitor::ddzNormalizedResiduals_tracker_
private

Definition at line 95 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dphiAbsoluteResiduals_global_

MonitorElement* TrackSplittingMonitor::dphiAbsoluteResiduals_global_
private

Definition at line 103 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dphiAbsoluteResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dphiAbsoluteResiduals_tracker_
private

Definition at line 89 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dphiNormalizedResiduals_global_

MonitorElement* TrackSplittingMonitor::dphiNormalizedResiduals_global_
private

Definition at line 110 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dphiNormalizedResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dphiNormalizedResiduals_tracker_
private

Definition at line 96 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dptAbsoluteResiduals_global_

MonitorElement* TrackSplittingMonitor::dptAbsoluteResiduals_global_
private

Definition at line 105 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dptAbsoluteResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dptAbsoluteResiduals_tracker_
private

Definition at line 91 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dptNormalizedResiduals_global_

MonitorElement* TrackSplittingMonitor::dptNormalizedResiduals_global_
private

Definition at line 112 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dptNormalizedResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dptNormalizedResiduals_tracker_
private

Definition at line 98 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dqmStore_

DQMStore* TrackSplittingMonitor::dqmStore_
private

Definition at line 58 of file TrackSplittingMonitor.h.

◆ dtGeometry

const DTGeometry* TrackSplittingMonitor::dtGeometry
private

Definition at line 69 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ dtGeomToken_

const edm::ESGetToken<DTGeometry, MuonGeometryRecord> TrackSplittingMonitor::dtGeomToken_
private

Definition at line 63 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ dthetaAbsoluteResiduals_global_

MonitorElement* TrackSplittingMonitor::dthetaAbsoluteResiduals_global_
private

Definition at line 104 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dthetaAbsoluteResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dthetaAbsoluteResiduals_tracker_
private

Definition at line 90 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dthetaNormalizedResiduals_global_

MonitorElement* TrackSplittingMonitor::dthetaNormalizedResiduals_global_
private

Definition at line 111 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dthetaNormalizedResiduals_tracker_

MonitorElement* TrackSplittingMonitor::dthetaNormalizedResiduals_tracker_
private

Definition at line 97 of file TrackSplittingMonitor.h.

Referenced by bookHistograms().

◆ dzCut_

double TrackSplittingMonitor::dzCut_
private

Definition at line 82 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().

◆ histname

std::string TrackSplittingMonitor::histname
private

Definition at line 56 of file TrackSplittingMonitor.h.

◆ mfToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> TrackSplittingMonitor::mfToken_
private

Definition at line 61 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ norchiCut_

double TrackSplittingMonitor::norchiCut_
private

Definition at line 84 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().

◆ pixelHitsPerLeg_

int TrackSplittingMonitor::pixelHitsPerLeg_
private

Definition at line 79 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().

◆ plotMuons_

bool TrackSplittingMonitor::plotMuons_
private

Definition at line 78 of file TrackSplittingMonitor.h.

Referenced by analyze(), bookHistograms(), and TrackSplittingMonitor().

◆ ptCut_

double TrackSplittingMonitor::ptCut_
private

◆ rpcGeometry

const RPCGeometry* TrackSplittingMonitor::rpcGeometry
private

Definition at line 71 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ rpcGeomToken_

const edm::ESGetToken<RPCGeometry, MuonGeometryRecord> TrackSplittingMonitor::rpcGeomToken_
private

Definition at line 65 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ splitMuons_

edm::InputTag TrackSplittingMonitor::splitMuons_
private

Definition at line 74 of file TrackSplittingMonitor.h.

Referenced by TrackSplittingMonitor().

◆ splitMuonsToken_

edm::EDGetTokenT<std::vector<reco::Muon> > TrackSplittingMonitor::splitMuonsToken_
private

Definition at line 76 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().

◆ splitTracks_

edm::InputTag TrackSplittingMonitor::splitTracks_
private

Definition at line 73 of file TrackSplittingMonitor.h.

Referenced by TrackSplittingMonitor().

◆ splitTracksToken_

edm::EDGetTokenT<std::vector<reco::Track> > TrackSplittingMonitor::splitTracksToken_
private

Definition at line 75 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().

◆ theGeometry

const TrackerGeometry* TrackSplittingMonitor::theGeometry
private

Definition at line 67 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ theMagField

const MagneticField* TrackSplittingMonitor::theMagField
private

Definition at line 68 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ tkGeomToken_

const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> TrackSplittingMonitor::tkGeomToken_
private

Definition at line 62 of file TrackSplittingMonitor.h.

Referenced by analyze().

◆ totalHitsPerLeg_

int TrackSplittingMonitor::totalHitsPerLeg_
private

Definition at line 80 of file TrackSplittingMonitor.h.

Referenced by analyze(), and TrackSplittingMonitor().