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

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

Inheritance diagram for TrackSplittingMonitor:
edm::EDAnalyzer

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
 
virtual void beginJob (void)
 
virtual void endJob (void)
 
 TrackSplittingMonitor (const edm::ParameterSet &)
 
 ~TrackSplittingMonitor ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

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

Private Attributes

edm::ParameterSet conf_
 
edm::ESHandle< CSCGeometrycscGeometry
 
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_
 
edm::ESHandle< DTGeometrydtGeometry
 
MonitorElementdthetaAbsoluteResiduals_global_
 
MonitorElementdthetaAbsoluteResiduals_tracker_
 
MonitorElementdthetaNormalizedResiduals_global_
 
MonitorElementdthetaNormalizedResiduals_tracker_
 
double dzCut_
 
std::string histname
 
double norchiCut_
 
int pixelHitsPerLeg_
 
bool plotMuons_
 
double ptCut_
 
edm::ESHandle< RPCGeometryrpcGeometry
 
edm::InputTag splitMuons_
 
edm::InputTag splitTracks_
 
edm::ESHandle< TrackerGeometrytheGeometry
 
edm::ESHandle< MagneticFieldtheMagField
 
int totalHitsPerLeg_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- 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::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

Monitoring source for general quantities related to tracks.

Definition at line 36 of file TrackSplittingMonitor.h.

Constructor & Destructor Documentation

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

Definition at line 42 of file TrackSplittingMonitor.cc.

42  {
43 }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 159 of file TrackSplittingMonitor.cc.

References cscGeometry, reco::TrackBase::d0(), d0Cut_, reco::TrackBase::d0Error(), ddxyAbsoluteResiduals_global_, ddxyAbsoluteResiduals_tracker_, ddxyNormalizedResiduals_global_, ddxyNormalizedResiduals_tracker_, dtGeometry, reco::TrackBase::dz(), dzCut_, reco::TrackBase::dzError(), MonitorElement::Fill(), edm::EventSetup::get(), edm::Event::getByLabel(), reco::Muon::globalTrack(), edm::HandleBase::isValid(), norchiCut_, reco::TrackBase::normalizedChi2(), reco::TrackBase::phi(), reco::TrackBase::phiError(), PixelSubdetector::PixelBarrel, pixelHitsPerLeg_, plotMuons_, funct::pow(), reco::TrackBase::pt(), ptCut_, reco::TrackBase::ptError(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), rpcGeometry, RecoMuonCosmics_cff::splitMuons, splitMuons_, splitTracks_, mathSSE::sqrt(), theGeometry, theMagField, reco::TrackBase::theta(), reco::TrackBase::thetaError(), and totalHitsPerLeg_.

159  {
160 
161 
162  iSetup.get<IdealMagneticFieldRecord>().get(theMagField);
164  iSetup.get<MuonGeometryRecord>().get(dtGeometry);
165  iSetup.get<MuonGeometryRecord>().get(cscGeometry);
166  iSetup.get<MuonGeometryRecord>().get(rpcGeometry);
167 
169  iEvent.getByLabel(splitTracks_, splitTracks);
170  if (!splitTracks.isValid()) return;
171 
173  if (plotMuons_){
174  iEvent.getByLabel(splitMuons_, splitMuons);
175  }
176 
177  if (splitTracks->size() == 2){
178  // check that there are 2 tracks in split track collection
179  edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size();
180 
181  // split tracks calculations
182  reco::Track track1 = splitTracks->at(0);
183  reco::Track track2 = splitTracks->at(1);
184 
185 
186  // -------------------------- basic selection ---------------------------
187 
188  // hit counting
189  // looping through the hits for track 1
190  double nRechits1 =0;
191  double nRechitinBPIX1 =0;
192  for (trackingRecHit_iterator iHit = track1.recHitsBegin(); iHit != track1.recHitsEnd(); ++iHit) {
193  if((*iHit)->isValid()) {
194  nRechits1++;
195  int type =(*iHit)->geographicalId().subdetId();
196  if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX1;}
197  }
198  }
199  // looping through the hits for track 2
200  double nRechits2 =0;
201  double nRechitinBPIX2 =0;
202  for (trackingRecHit_iterator iHit = track2.recHitsBegin(); iHit != track2.recHitsEnd(); ++iHit) {
203  if((*iHit)->isValid()) {
204  nRechits2++;
205  int type =(*iHit)->geographicalId().subdetId();
206  if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX2;}
207  }
208  }
209 
210  // DCA of each track
211  double d01 = track1.d0();
212  double dz1 = track1.dz();
213  double d02 = track2.d0();
214  double dz2 = track2.dz();
215 
216  // pT of each track
217  double pt1 = track1.pt();
218  double pt2 = track2.pt();
219 
220  // chi2 of each track
221  double norchi1 = track1.normalizedChi2();
222  double norchi2 = track2.normalizedChi2();
223 
224  // basic selection
225  // pixel hits and total hits
226  if ((nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechits1 >= totalHitsPerLeg_)&&(nRechits2 >= totalHitsPerLeg_)){
227  // dca cut
228  if ( ((fabs(d01) < d0Cut_))&&(fabs(d02) < d0Cut_)&&(fabs(dz2) < dzCut_)&&(fabs(dz2) < dzCut_) ){
229  // pt cut
230  if ( (pt1+pt2)/2 < ptCut_){
231  // chi2 cut
232  if ( (norchi1 < norchiCut_)&&(norchi2 < norchiCut_) ){
233 
234  // passed all cuts...
235  edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?";
236 
237  double ddxyVal = d01 - d02;
238  double ddzVal = dz1 - dz2;
239  double dphiVal = track1.phi() - track2.phi();
240  double dthetaVal = track1.theta() - track2.theta();
241  double dptVal = pt1 - pt2;
242  double dcurvVal = (1/pt1) - (1/pt2);
243 
244  double d01ErrVal = track1.d0Error();
245  double d02ErrVal = track2.d0Error();
246  double dz1ErrVal = track1.dzError();
247  double dz2ErrVal = track2.dzError();
248  double phi1ErrVal = track1.phiError();
249  double phi2ErrVal = track2.phiError();
250  double theta1ErrVal = track1.thetaError();
251  double theta2ErrVal = track2.thetaError();
252  double pt1ErrVal = track1.ptError();
253  double pt2ErrVal = track2.ptError();
254 
255  ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddxyVal/sqrt(2.0) );
256  ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddzVal/sqrt(2.0) );
257  ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dphiVal/sqrt(2.0) );
258  ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dthetaVal/sqrt(2.0) );
259  ddxyAbsoluteResiduals_tracker_->Fill( dptVal/sqrt(2.0) );
260  ddxyAbsoluteResiduals_tracker_->Fill( dcurvVal/sqrt(2.0) );
261 
262  ddxyNormalizedResiduals_tracker_->Fill( ddxyVal/sqrt( d01ErrVal*d01ErrVal + d02ErrVal*d02ErrVal ) );
263  ddxyNormalizedResiduals_tracker_->Fill( ddzVal/sqrt( dz1ErrVal*dz1ErrVal + dz2ErrVal*dz2ErrVal ) );
264  ddxyNormalizedResiduals_tracker_->Fill( dphiVal/sqrt( phi1ErrVal*phi1ErrVal + phi2ErrVal*phi2ErrVal ) );
265  ddxyNormalizedResiduals_tracker_->Fill( dthetaVal/sqrt( theta1ErrVal*theta1ErrVal + theta2ErrVal*theta2ErrVal ) );
266  ddxyNormalizedResiduals_tracker_->Fill( dptVal/sqrt( pt1ErrVal*pt1ErrVal + pt2ErrVal*pt2ErrVal ) );
267  ddxyNormalizedResiduals_tracker_->Fill( dcurvVal/sqrt( pow(pt1ErrVal,2)/pow(pt1,4) + pow(pt2ErrVal,2)/pow(pt2,4) ) );
268 
269  // if do the same for split muons
270  if (plotMuons_ && splitMuons.isValid()){
271 
272  int gmCtr = 0;
273  bool topGlobalMuonFlag = false;
274  bool bottomGlobalMuonFlag = false;
275  int topGlobalMuon = -1;
276  int bottomGlobalMuon = -1;
277  double topGlobalMuonNorchi2 = 1e10;
278  double bottomGlobalMuonNorchi2 = 1e10;
279 
280  // check if usable split global muons
281  for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++){
282  if ( gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon() ){
283 
284  reco::TrackRef trackerTrackRef1( splitTracks, 0 );
285  reco::TrackRef trackerTrackRef2( splitTracks, 1 );
286 
287  if (gmI->innerTrack() == trackerTrackRef1){
288  if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2){
289  topGlobalMuonFlag = true;
290  topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
291  topGlobalMuon = gmCtr;
292  }
293  }
294  if (gmI->innerTrack() == trackerTrackRef2){
295  if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2){
296  bottomGlobalMuonFlag = true;
297  bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
298  bottomGlobalMuon = gmCtr;
299  }
300  }
301  }
302  gmCtr++;
303  }
304 
305  if (bottomGlobalMuonFlag && topGlobalMuonFlag) {
306 
307  reco::Muon muonTop = splitMuons->at( topGlobalMuon );
308  reco::Muon muonBottom = splitMuons->at( bottomGlobalMuon );
309 
310  reco::TrackRef glb1 = muonTop.globalTrack();
311  reco::TrackRef glb2 = muonBottom.globalTrack();
312 
313  double ddxyValGlb = glb1->d0() - glb2->d0();
314  double ddzValGlb = glb1->dz() - glb2->dz();
315  double dphiValGlb = glb1->phi() - glb2->phi();
316  double dthetaValGlb = glb1->theta() - glb2->theta();
317  double dptValGlb = glb1->pt() - glb2->pt();
318  double dcurvValGlb = (1/glb1->pt()) - (1/glb2->pt());
319 
320  double d01ErrValGlb = glb1->d0Error();
321  double d02ErrValGlb = glb2->d0Error();
322  double dz1ErrValGlb = glb1->dzError();
323  double dz2ErrValGlb = glb2->dzError();
324  double phi1ErrValGlb = glb1->phiError();
325  double phi2ErrValGlb = glb2->phiError();
326  double theta1ErrValGlb = glb1->thetaError();
327  double theta2ErrValGlb = glb2->thetaError();
328  double pt1ErrValGlb = glb1->ptError();
329  double pt2ErrValGlb = glb2->ptError();
330 
331  ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddxyValGlb/sqrt(2.0) );
332  ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddzValGlb/sqrt(2.0) );
333  ddxyAbsoluteResiduals_global_->Fill( 1000.0*dphiValGlb/sqrt(2.0) );
334  ddxyAbsoluteResiduals_global_->Fill( 1000.0*dthetaValGlb/sqrt(2.0) );
335  ddxyAbsoluteResiduals_global_->Fill( dptValGlb/sqrt(2.0) );
336  ddxyAbsoluteResiduals_global_->Fill( dcurvValGlb/sqrt(2.0) );
337 
338  ddxyNormalizedResiduals_global_->Fill( ddxyValGlb/sqrt( d01ErrValGlb*d01ErrValGlb + d02ErrValGlb*d02ErrValGlb ) );
339  ddxyNormalizedResiduals_global_->Fill( ddzValGlb/sqrt( dz1ErrValGlb*dz1ErrValGlb + dz2ErrValGlb*dz2ErrValGlb ) );
340  ddxyNormalizedResiduals_global_->Fill( dphiValGlb/sqrt( phi1ErrValGlb*phi1ErrValGlb + phi2ErrValGlb*phi2ErrValGlb ) );
341  ddxyNormalizedResiduals_global_->Fill( dthetaValGlb/sqrt( theta1ErrValGlb*theta1ErrValGlb + theta2ErrValGlb*theta2ErrValGlb ) );
342  ddxyNormalizedResiduals_global_->Fill( dptValGlb/sqrt( pt1ErrValGlb*pt1ErrValGlb + pt2ErrValGlb*pt2ErrValGlb ) );
343  ddxyNormalizedResiduals_global_->Fill( dcurvValGlb/sqrt( pow(pt1ErrValGlb,2)/pow(pt1,4) + pow(pt2ErrValGlb,2)/pow(pt2,4) ) );
344 
345  }
346 
347 
348  } // end of split muons loop
349  }
350  }
351  }
352  }
353  }
354 }
type
Definition: HCALResponse.h:22
edm::ESHandle< TrackerGeometry > theGeometry
double d0Error() const
error on d0
Definition: TrackBase.h:211
double d0() const
dxy parameter in perigee convention (d0 = - dxy)
Definition: TrackBase.h:123
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:111
double theta() const
polar angle
Definition: TrackBase.h:117
MonitorElement * ddxyNormalizedResiduals_tracker_
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:139
edm::ESHandle< RPCGeometry > rpcGeometry
void Fill(long long x)
MonitorElement * ddxyNormalizedResiduals_global_
T sqrt(T t)
Definition: SSEVec.h:46
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:194
double phiError() const
error on phi
Definition: TrackBase.h:207
edm::ESHandle< MagneticField > theMagField
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
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:127
double dzError() const
error on dz
Definition: TrackBase.h:215
const T & get() const
Definition: EventSetup.h:55
MonitorElement * ddxyAbsoluteResiduals_global_
edm::ESHandle< DTGeometry > dtGeometry
edm::ESHandle< CSCGeometry > cscGeometry
MonitorElement * ddxyAbsoluteResiduals_tracker_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double thetaError() const
error on theta
Definition: TrackBase.h:201
virtual TrackRef globalTrack() const
reference to Track reconstructed in both tracked and muon detector
Definition: Muon.h:55
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
void TrackSplittingMonitor::beginJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 45 of file TrackSplittingMonitor.cc.

References DQMStore::book1D(), conf_, d0Cut_, dcurvAbsoluteResiduals_global_, dcurvAbsoluteResiduals_tracker_, dcurvNormalizedResiduals_global_, dcurvNormalizedResiduals_tracker_, ddxyAbsoluteResiduals_global_, ddxyAbsoluteResiduals_tracker_, ddxyNormalizedResiduals_global_, ddxyNormalizedResiduals_tracker_, ddzAbsoluteResiduals_global_, ddzAbsoluteResiduals_tracker_, ddzNormalizedResiduals_global_, ddzNormalizedResiduals_tracker_, dphiAbsoluteResiduals_global_, dphiAbsoluteResiduals_tracker_, dphiNormalizedResiduals_global_, dphiNormalizedResiduals_tracker_, dptAbsoluteResiduals_global_, dptAbsoluteResiduals_tracker_, dptNormalizedResiduals_global_, dptNormalizedResiduals_tracker_, dqmStore_, dthetaAbsoluteResiduals_global_, dthetaAbsoluteResiduals_tracker_, dthetaNormalizedResiduals_global_, dthetaNormalizedResiduals_tracker_, dzCut_, edm::ParameterSet::getParameter(), norchiCut_, pixelHitsPerLeg_, plotMuons_, ptCut_, MonitorElement::setAxisTitle(), DQMStore::setCurrentFolder(), splitMuons_, splitTracks_, and totalHitsPerLeg_.

45  {
46 
47  std::string MEFolderName = conf_.getParameter<std::string>("FolderName");
48  dqmStore_->setCurrentFolder(MEFolderName);
49 
50  //get input tags
51  splitTracks_ = conf_.getParameter< edm::InputTag >("splitTrackCollection");
52  splitMuons_ = conf_.getParameter< edm::InputTag >("splitMuonCollection");
53  plotMuons_ = conf_.getParameter<bool>("ifPlotMuons");
54 
55  // cuts
56  pixelHitsPerLeg_ = conf_.getParameter<int>("pixelHitsPerLeg");
57  totalHitsPerLeg_ = conf_.getParameter<int>("totalHitsPerLeg");
58  d0Cut_ = conf_.getParameter<double>("d0Cut");
59  dzCut_ = conf_.getParameter<double>("dzCut");
60  ptCut_ = conf_.getParameter<double>("ptCut");
61  norchiCut_ = conf_.getParameter<double>("norchiCut");
62 
63 
64  // bin declarations
65  int ddxyBin = conf_.getParameter<int>("ddxyBin");
66  double ddxyMin = conf_.getParameter<double>("ddxyMin");
67  double ddxyMax = conf_.getParameter<double>("ddxyMax");
68 
69  int ddzBin = conf_.getParameter<int>("ddzBin");
70  double ddzMin = conf_.getParameter<double>("ddzMin");
71  double ddzMax = conf_.getParameter<double>("ddzMax");
72 
73  int dphiBin = conf_.getParameter<int>("dphiBin");
74  double dphiMin = conf_.getParameter<double>("dphiMin");
75  double dphiMax = conf_.getParameter<double>("dphiMax");
76 
77  int dthetaBin = conf_.getParameter<int>("dthetaBin");
78  double dthetaMin = conf_.getParameter<double>("dthetaMin");
79  double dthetaMax = conf_.getParameter<double>("dthetaMax");
80 
81  int dptBin = conf_.getParameter<int>("dptBin");
82  double dptMin = conf_.getParameter<double>("dptMin");
83  double dptMax = conf_.getParameter<double>("dptMax");
84 
85  int dcurvBin = conf_.getParameter<int>("dcurvBin");
86  double dcurvMin = conf_.getParameter<double>("dcurvMin");
87  double dcurvMax = conf_.getParameter<double>("dcurvMax");
88 
89  int normBin = conf_.getParameter<int>("normBin");
90  double normMin = conf_.getParameter<double>("normMin");
91  double normMax = conf_.getParameter<double>("normMax");
92 
93  // declare histogram
94  ddxyAbsoluteResiduals_tracker_ = dqmStore_->book1D( "ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax );
95  ddzAbsoluteResiduals_tracker_ = dqmStore_->book1D( "ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax );
96  dphiAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax );
97  dthetaAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax );
98  dptAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax );
99  dcurvAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax );
100 
101  ddxyNormalizedResiduals_tracker_ = dqmStore_->book1D( "ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax );
102  ddzNormalizedResiduals_tracker_ = dqmStore_->book1D( "ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax );
103  dphiNormalizedResiduals_tracker_ = dqmStore_->book1D( "dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax );
104  dthetaNormalizedResiduals_tracker_ = dqmStore_->book1D( "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax );
105  dptNormalizedResiduals_tracker_ = dqmStore_->book1D( "dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax );
106  dcurvNormalizedResiduals_tracker_ = dqmStore_->book1D( "dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax );
107 
108  if (plotMuons_){
109  ddxyAbsoluteResiduals_global_ = dqmStore_->book1D( "ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax );
110  ddzAbsoluteResiduals_global_ = dqmStore_->book1D( "ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax );
111  dphiAbsoluteResiduals_global_ = dqmStore_->book1D( "dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax );
112  dthetaAbsoluteResiduals_global_ = dqmStore_->book1D( "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax );
113  dptAbsoluteResiduals_global_ = dqmStore_->book1D( "dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax );
114  dcurvAbsoluteResiduals_global_ = dqmStore_->book1D( "dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax );
115 
116  ddxyNormalizedResiduals_global_ = dqmStore_->book1D( "ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax );
117  ddzNormalizedResiduals_global_ = dqmStore_->book1D( "ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax );
118  dphiNormalizedResiduals_global_ = dqmStore_->book1D( "dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax );
119  dthetaNormalizedResiduals_global_ = dqmStore_->book1D( "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax );
120  dptNormalizedResiduals_global_ = dqmStore_->book1D( "dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax );
121  dcurvNormalizedResiduals_global_ = dqmStore_->book1D( "dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax );
122  }
123 
124  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" );
125  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" );
126  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" );
127  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" );
128  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" );
129  ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" );
130 
131  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" );
132  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" );
133  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" );
134  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" );
135  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" );
136  ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" );
137 
138  if (plotMuons_){
139  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" );
140  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" );
141  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" );
142  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" );
143  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" );
144  ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" );
145 
146  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" );
147  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" );
148  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" );
149  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" );
150  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" );
151  ddxyNormalizedResiduals_global_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" );
152  }
153 
154 }
T getParameter(std::string const &) const
MonitorElement * dthetaAbsoluteResiduals_global_
MonitorElement * dcurvNormalizedResiduals_global_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
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 * ddxyAbsoluteResiduals_tracker_
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
void TrackSplittingMonitor::doProfileX ( TH2 *  th2,
MonitorElement me 
)
private
void TrackSplittingMonitor::doProfileX ( MonitorElement th2m,
MonitorElement me 
)
private
void TrackSplittingMonitor::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 357 of file TrackSplittingMonitor.cc.

References conf_, dqmStore_, edm::ParameterSet::getParameter(), dumpDBToFile_GT_ttrig_cfg::outputFileName, DQMStore::save(), and DQMStore::showDirStructure().

357  {
358  bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile");
359  std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
360  if(outputMEsInRootFile){
362  dqmStore_->save(outputFileName);
363  }
364 }
T getParameter(std::string const &) const
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2113
void showDirStructure(void) const
Definition: DQMStore.cc:2761

Member Data Documentation

edm::ParameterSet TrackSplittingMonitor::conf_
private

Definition at line 57 of file TrackSplittingMonitor.h.

Referenced by beginJob(), endJob(), and TrackSplittingMonitor().

edm::ESHandle<CSCGeometry> TrackSplittingMonitor::cscGeometry
private

Definition at line 62 of file TrackSplittingMonitor.h.

Referenced by analyze().

double TrackSplittingMonitor::d0Cut_
private

Definition at line 72 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

MonitorElement* TrackSplittingMonitor::dcurvAbsoluteResiduals_global_
private

Definition at line 98 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dcurvAbsoluteResiduals_tracker_
private

Definition at line 84 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dcurvNormalizedResiduals_global_
private

Definition at line 105 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dcurvNormalizedResiduals_tracker_
private

Definition at line 91 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::ddxyAbsoluteResiduals_global_
private

Definition at line 93 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

MonitorElement* TrackSplittingMonitor::ddxyAbsoluteResiduals_tracker_
private

Definition at line 79 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

MonitorElement* TrackSplittingMonitor::ddxyNormalizedResiduals_global_
private

Definition at line 100 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

MonitorElement* TrackSplittingMonitor::ddxyNormalizedResiduals_tracker_
private

Definition at line 86 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

MonitorElement* TrackSplittingMonitor::ddzAbsoluteResiduals_global_
private

Definition at line 94 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::ddzAbsoluteResiduals_tracker_
private

Definition at line 80 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::ddzNormalizedResiduals_global_
private

Definition at line 101 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::ddzNormalizedResiduals_tracker_
private

Definition at line 87 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dphiAbsoluteResiduals_global_
private

Definition at line 95 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dphiAbsoluteResiduals_tracker_
private

Definition at line 81 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dphiNormalizedResiduals_global_
private

Definition at line 102 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dphiNormalizedResiduals_tracker_
private

Definition at line 88 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dptAbsoluteResiduals_global_
private

Definition at line 97 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dptAbsoluteResiduals_tracker_
private

Definition at line 83 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dptNormalizedResiduals_global_
private

Definition at line 104 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dptNormalizedResiduals_tracker_
private

Definition at line 90 of file TrackSplittingMonitor.h.

Referenced by beginJob().

DQMStore* TrackSplittingMonitor::dqmStore_
private

Definition at line 56 of file TrackSplittingMonitor.h.

Referenced by beginJob(), endJob(), and TrackSplittingMonitor().

edm::ESHandle<DTGeometry> TrackSplittingMonitor::dtGeometry
private

Definition at line 61 of file TrackSplittingMonitor.h.

Referenced by analyze().

MonitorElement* TrackSplittingMonitor::dthetaAbsoluteResiduals_global_
private

Definition at line 96 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dthetaAbsoluteResiduals_tracker_
private

Definition at line 82 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dthetaNormalizedResiduals_global_
private

Definition at line 103 of file TrackSplittingMonitor.h.

Referenced by beginJob().

MonitorElement* TrackSplittingMonitor::dthetaNormalizedResiduals_tracker_
private

Definition at line 89 of file TrackSplittingMonitor.h.

Referenced by beginJob().

double TrackSplittingMonitor::dzCut_
private

Definition at line 73 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

std::string TrackSplittingMonitor::histname
private

Definition at line 54 of file TrackSplittingMonitor.h.

double TrackSplittingMonitor::norchiCut_
private

Definition at line 75 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

int TrackSplittingMonitor::pixelHitsPerLeg_
private

Definition at line 70 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

bool TrackSplittingMonitor::plotMuons_
private

Definition at line 69 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

double TrackSplittingMonitor::ptCut_
private

Definition at line 74 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

edm::ESHandle<RPCGeometry> TrackSplittingMonitor::rpcGeometry
private

Definition at line 63 of file TrackSplittingMonitor.h.

Referenced by analyze().

edm::InputTag TrackSplittingMonitor::splitMuons_
private

Definition at line 66 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

edm::InputTag TrackSplittingMonitor::splitTracks_
private

Definition at line 65 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().

edm::ESHandle<TrackerGeometry> TrackSplittingMonitor::theGeometry
private

Definition at line 59 of file TrackSplittingMonitor.h.

Referenced by analyze().

edm::ESHandle<MagneticField> TrackSplittingMonitor::theMagField
private

Definition at line 60 of file TrackSplittingMonitor.h.

Referenced by analyze().

int TrackSplittingMonitor::totalHitsPerLeg_
private

Definition at line 71 of file TrackSplittingMonitor.h.

Referenced by analyze(), and beginJob().