CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/DQM/TrackingMonitor/plugins/TrackSplittingMonitor.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2010/03/27 10:38:10 $
00005  *  $Revision: 1.5 $
00006  *  \author Suchandra Dutta , Giorgia Mila
00007  */
00008 
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00011 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00012 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h" 
00013 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00014 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00015 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00016 #include "MagneticField/Engine/interface/MagneticField.h"
00017 
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 #include "FWCore/Utilities/interface/InputTag.h"
00020 #include "DQMServices/Core/interface/DQMStore.h"
00021 #include "DQM/TrackingMonitor/interface/TrackAnalyzer.h"
00022 #include "DQM/TrackingMonitor/plugins/TrackSplittingMonitor.h"
00023 
00024 #include "FWCore/Framework/interface/ESHandle.h"
00025 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h" 
00026 
00027 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00028 #include "TrackingTools/Records/interface/TransientRecHitRecord.h" 
00029 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00030 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00031 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00032 #include <string>
00033 
00034 #include "DataFormats/MuonReco/interface/Muon.h"
00035 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00036 
00037 TrackSplittingMonitor::TrackSplittingMonitor(const edm::ParameterSet& iConfig) {
00038         dqmStore_ = edm::Service<DQMStore>().operator->();
00039         conf_ = iConfig;
00040 }
00041 
00042 TrackSplittingMonitor::~TrackSplittingMonitor() { 
00043 }
00044 
00045 void TrackSplittingMonitor::beginJob(void) {
00046         
00047   std::string MEFolderName = conf_.getParameter<std::string>("FolderName"); 
00048   dqmStore_->setCurrentFolder(MEFolderName);
00049   
00050   //get input tags
00051   splitTracks_ = conf_.getParameter< edm::InputTag >("splitTrackCollection");
00052   splitMuons_ = conf_.getParameter< edm::InputTag >("splitMuonCollection");
00053   plotMuons_ = conf_.getParameter<bool>("ifPlotMuons");
00054   
00055   // cuts
00056   pixelHitsPerLeg_ = conf_.getParameter<int>("pixelHitsPerLeg");
00057   totalHitsPerLeg_ = conf_.getParameter<int>("totalHitsPerLeg");
00058   d0Cut_ = conf_.getParameter<double>("d0Cut");
00059   dzCut_ = conf_.getParameter<double>("dzCut");
00060   ptCut_ = conf_.getParameter<double>("ptCut");
00061   norchiCut_ = conf_.getParameter<double>("norchiCut");
00062   
00063   
00064   // bin declarations
00065   int    ddxyBin = conf_.getParameter<int>("ddxyBin");
00066   double ddxyMin = conf_.getParameter<double>("ddxyMin");
00067   double ddxyMax = conf_.getParameter<double>("ddxyMax");
00068   
00069   int    ddzBin = conf_.getParameter<int>("ddzBin");
00070   double ddzMin = conf_.getParameter<double>("ddzMin");
00071   double ddzMax = conf_.getParameter<double>("ddzMax");
00072   
00073   int    dphiBin = conf_.getParameter<int>("dphiBin");
00074   double dphiMin = conf_.getParameter<double>("dphiMin");
00075   double dphiMax = conf_.getParameter<double>("dphiMax");
00076   
00077   int    dthetaBin = conf_.getParameter<int>("dthetaBin");
00078   double dthetaMin = conf_.getParameter<double>("dthetaMin");
00079   double dthetaMax = conf_.getParameter<double>("dthetaMax");
00080   
00081   int    dptBin = conf_.getParameter<int>("dptBin");
00082   double dptMin = conf_.getParameter<double>("dptMin");
00083   double dptMax = conf_.getParameter<double>("dptMax");
00084   
00085   int    dcurvBin = conf_.getParameter<int>("dcurvBin");
00086   double dcurvMin = conf_.getParameter<double>("dcurvMin");
00087   double dcurvMax = conf_.getParameter<double>("dcurvMax");
00088   
00089   int    normBin = conf_.getParameter<int>("normBin");
00090   double normMin = conf_.getParameter<double>("normMin");
00091   double normMax = conf_.getParameter<double>("normMax");       
00092         
00093   // declare histogram
00094   ddxyAbsoluteResiduals_tracker_ = dqmStore_->book1D( "ddxyAbsoluteResiduals_tracker", "ddxyAbsoluteResiduals_tracker", ddxyBin, ddxyMin, ddxyMax );
00095   ddzAbsoluteResiduals_tracker_ = dqmStore_->book1D( "ddzAbsoluteResiduals_tracker", "ddzAbsoluteResiduals_tracker", ddzBin, ddzMin, ddzMax );
00096   dphiAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dphiAbsoluteResiduals_tracker", "dphiAbsoluteResiduals_tracker", dphiBin, dphiMin, dphiMax );
00097   dthetaAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dthetaAbsoluteResiduals_tracker", "dthetaAbsoluteResiduals_tracker", dthetaBin, dthetaMin, dthetaMax );
00098   dptAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dptAbsoluteResiduals_tracker", "dptAbsoluteResiduals_tracker", dptBin, dptMin, dptMax );
00099   dcurvAbsoluteResiduals_tracker_ = dqmStore_->book1D( "dcurvAbsoluteResiduals_tracker", "dcurvAbsoluteResiduals_tracker", dcurvBin, dcurvMin, dcurvMax );
00100   
00101   ddxyNormalizedResiduals_tracker_ = dqmStore_->book1D( "ddxyNormalizedResiduals_tracker", "ddxyNormalizedResiduals_tracker", normBin, normMin, normMax );
00102   ddzNormalizedResiduals_tracker_ = dqmStore_->book1D( "ddzNormalizedResiduals_tracker", "ddzNormalizedResiduals_tracker", normBin, normMin, normMax );
00103   dphiNormalizedResiduals_tracker_ = dqmStore_->book1D( "dphiNormalizedResiduals_tracker", "dphiNormalizedResiduals_tracker", normBin, normMin, normMax );
00104   dthetaNormalizedResiduals_tracker_ = dqmStore_->book1D( "dthetaNormalizedResiduals_tracker", "dthetaNormalizedResiduals_tracker", normBin, normMin, normMax );
00105   dptNormalizedResiduals_tracker_ = dqmStore_->book1D( "dptNormalizedResiduals_tracker", "dptNormalizedResiduals_tracker", normBin, normMin, normMax );
00106   dcurvNormalizedResiduals_tracker_ = dqmStore_->book1D( "dcurvNormalizedResiduals_tracker", "dcurvNormalizedResiduals_tracker", normBin, normMin, normMax );
00107         
00108   if (plotMuons_){
00109     ddxyAbsoluteResiduals_global_ = dqmStore_->book1D( "ddxyAbsoluteResiduals_global", "ddxyAbsoluteResiduals_global", ddxyBin, ddxyMin, ddxyMax );
00110     ddzAbsoluteResiduals_global_ = dqmStore_->book1D( "ddzAbsoluteResiduals_global", "ddzAbsoluteResiduals_global", ddzBin, ddzMin, ddzMax );
00111     dphiAbsoluteResiduals_global_ = dqmStore_->book1D( "dphiAbsoluteResiduals_global", "dphiAbsoluteResiduals_global", dphiBin, dphiMin, dphiMax );
00112     dthetaAbsoluteResiduals_global_ = dqmStore_->book1D( "dthetaAbsoluteResiduals_global", "dthetaAbsoluteResiduals_global", dthetaBin, dthetaMin, dthetaMax );
00113     dptAbsoluteResiduals_global_ = dqmStore_->book1D( "dptAbsoluteResiduals_global", "dptAbsoluteResiduals_global", dptBin, dptMin, dptMax );
00114     dcurvAbsoluteResiduals_global_ = dqmStore_->book1D( "dcurvAbsoluteResiduals_global", "dcurvAbsoluteResiduals_global", dcurvBin, dcurvMin, dcurvMax );
00115     
00116     ddxyNormalizedResiduals_global_ = dqmStore_->book1D( "ddxyNormalizedResiduals_global", "ddxyNormalizedResiduals_global", normBin, normMin, normMax );
00117     ddzNormalizedResiduals_global_ = dqmStore_->book1D( "ddzNormalizedResiduals_global", "ddzNormalizedResiduals_global", normBin, normMin, normMax );
00118     dphiNormalizedResiduals_global_ = dqmStore_->book1D( "dphiNormalizedResiduals_global", "dphiNormalizedResiduals_global", normBin, normMin, normMax );
00119     dthetaNormalizedResiduals_global_ = dqmStore_->book1D( "dthetaNormalizedResiduals_global", "dthetaNormalizedResiduals_global", normBin, normMin, normMax );
00120     dptNormalizedResiduals_global_ = dqmStore_->book1D( "dptNormalizedResiduals_global", "dptNormalizedResiduals_global", normBin, normMin, normMax );
00121     dcurvNormalizedResiduals_global_ = dqmStore_->book1D( "dcurvNormalizedResiduals_global", "dcurvNormalizedResiduals_global", normBin, normMin, normMax );
00122   }
00123   
00124   ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" );
00125   ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" );
00126   ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" );
00127   ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" );
00128   ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" );
00129   ddxyAbsoluteResiduals_tracker_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" );
00130   
00131   ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" );
00132   ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" );
00133   ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" );
00134   ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" );
00135   ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" );
00136   ddxyNormalizedResiduals_tracker_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" );
00137   
00138   if (plotMuons_){
00139     ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{xy})/#sqrt{2} [#mum]" );
00140     ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta d_{z})/#sqrt{2} [#mum]" );
00141     ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #phi)/#sqrt{2} [mrad]" );
00142     ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta #theta)/#sqrt{2} [mrad]" );
00143     ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta pT)/#sqrt{2} [GeV]" );
00144     ddxyAbsoluteResiduals_global_->setAxisTitle( "(#delta (1/pT))/#sqrt{2} [GeV^{-1}]" );
00145     
00146     ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{xy}/#sigma(d_{xy}" );
00147     ddxyNormalizedResiduals_global_->setAxisTitle( "#delta d_{z}/#sigma(d_{z})" );
00148     ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #phi/#sigma(d_{#phi})" );
00149     ddxyNormalizedResiduals_global_->setAxisTitle( "#delta #theta/#sigma(d_{#theta})" );
00150     ddxyNormalizedResiduals_global_->setAxisTitle( "#delta p_{T}/#sigma(p_{T})" );
00151     ddxyNormalizedResiduals_global_->setAxisTitle( "#delta 1/p_{T}/#sigma(1/p_{T})" );
00152   }
00153           
00154 }
00155 
00156 //
00157 // -- Analyse
00158 //
00159 void TrackSplittingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
00160         
00161   
00162   iSetup.get<IdealMagneticFieldRecord>().get(theMagField);   
00163   iSetup.get<TrackerDigiGeometryRecord>().get(theGeometry);
00164   iSetup.get<MuonGeometryRecord>().get(dtGeometry);
00165   iSetup.get<MuonGeometryRecord>().get(cscGeometry);
00166   iSetup.get<MuonGeometryRecord>().get(rpcGeometry);
00167   
00168   edm::Handle<std::vector<reco::Track> > splitTracks;
00169   iEvent.getByLabel(splitTracks_, splitTracks);
00170   if (!splitTracks.isValid()) return;
00171 
00172   edm::Handle<std::vector<reco::Muon> > splitMuons;
00173   if (plotMuons_){
00174     iEvent.getByLabel(splitMuons_, splitMuons);
00175   }
00176   
00177   if (splitTracks->size() == 2){
00178     // check that there are 2 tracks in split track collection
00179     edm::LogInfo("TrackSplittingMonitor") << "Split Track size: " << splitTracks->size();
00180     
00181     // split tracks calculations
00182     reco::Track track1 = splitTracks->at(0);
00183     reco::Track track2 = splitTracks->at(1);
00184     
00185     
00186     // -------------------------- basic selection ---------------------------
00187     
00188     // hit counting
00189     // looping through the hits for track 1
00190     double nRechits1 =0;
00191     double nRechitinBPIX1 =0;
00192     for (trackingRecHit_iterator iHit = track1.recHitsBegin(); iHit != track1.recHitsEnd(); ++iHit) {
00193       if((*iHit)->isValid()) {       
00194         nRechits1++;
00195         int type =(*iHit)->geographicalId().subdetId();
00196         if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX1;}
00197       }
00198     }    
00199     // looping through the hits for track 2
00200     double nRechits2 =0;
00201     double nRechitinBPIX2 =0;
00202     for (trackingRecHit_iterator iHit = track2.recHitsBegin(); iHit != track2.recHitsEnd(); ++iHit) {
00203       if((*iHit)->isValid()) {       
00204         nRechits2++;
00205         int type =(*iHit)->geographicalId().subdetId();
00206         if(type==int(PixelSubdetector::PixelBarrel)){++nRechitinBPIX2;}
00207       }
00208     }    
00209     
00210     // DCA of each track
00211     double d01 = track1.d0();
00212     double dz1 = track1.dz();
00213     double d02 = track2.d0();
00214     double dz2 = track2.dz();
00215     
00216     // pT of each track
00217     double pt1 = track1.pt();
00218     double pt2 = track2.pt();
00219     
00220     // chi2 of each track
00221     double norchi1 = track1.normalizedChi2();
00222     double norchi2 = track2.normalizedChi2();           
00223     
00224     // basic selection
00225     // pixel hits and total hits
00226     if ((nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechitinBPIX1 >= pixelHitsPerLeg_)&&(nRechits1 >= totalHitsPerLeg_)&&(nRechits2 >= totalHitsPerLeg_)){
00227       // dca cut
00228       if ( ((fabs(d01) < d0Cut_))&&(fabs(d02) < d0Cut_)&&(fabs(dz2) < dzCut_)&&(fabs(dz2) < dzCut_) ){
00229         // pt cut
00230         if ( (pt1+pt2)/2 < ptCut_){
00231           // chi2 cut
00232           if ( (norchi1 < norchiCut_)&&(norchi2 < norchiCut_) ){
00233             
00234             // passed all cuts...
00235             edm::LogInfo("TrackSplittingMonitor") << " Setected after all cuts ?";
00236 
00237             double ddxyVal = d01 - d02;
00238             double ddzVal = dz1 - dz2;                                          
00239             double dphiVal = track1.phi() - track2.phi();
00240             double dthetaVal = track1.theta() - track2.theta();                                         
00241             double dptVal = pt1 - pt2;
00242             double dcurvVal = (1/pt1) - (1/pt2);
00243             
00244             double d01ErrVal = track1.d0Error();
00245             double d02ErrVal = track2.d0Error();
00246             double dz1ErrVal = track1.dzError();
00247             double dz2ErrVal = track2.dzError();
00248             double phi1ErrVal = track1.phiError();
00249             double phi2ErrVal = track2.phiError();
00250             double theta1ErrVal = track1.thetaError();
00251             double theta2ErrVal = track2.thetaError();
00252             double pt1ErrVal = track1.ptError();
00253             double pt2ErrVal = track2.ptError();
00254             
00255             ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddxyVal/sqrt(2.0) );
00256             ddxyAbsoluteResiduals_tracker_->Fill( 10000.0*ddzVal/sqrt(2.0) );
00257             ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dphiVal/sqrt(2.0) );
00258             ddxyAbsoluteResiduals_tracker_->Fill( 1000.0*dthetaVal/sqrt(2.0) );
00259             ddxyAbsoluteResiduals_tracker_->Fill( dptVal/sqrt(2.0) );
00260             ddxyAbsoluteResiduals_tracker_->Fill( dcurvVal/sqrt(2.0) );
00261             
00262             ddxyNormalizedResiduals_tracker_->Fill( ddxyVal/sqrt( d01ErrVal*d01ErrVal + d02ErrVal*d02ErrVal ) );
00263             ddxyNormalizedResiduals_tracker_->Fill( ddzVal/sqrt( dz1ErrVal*dz1ErrVal + dz2ErrVal*dz2ErrVal ) );
00264             ddxyNormalizedResiduals_tracker_->Fill( dphiVal/sqrt( phi1ErrVal*phi1ErrVal + phi2ErrVal*phi2ErrVal ) );
00265             ddxyNormalizedResiduals_tracker_->Fill( dthetaVal/sqrt( theta1ErrVal*theta1ErrVal + theta2ErrVal*theta2ErrVal ) );
00266             ddxyNormalizedResiduals_tracker_->Fill( dptVal/sqrt( pt1ErrVal*pt1ErrVal + pt2ErrVal*pt2ErrVal ) );
00267             ddxyNormalizedResiduals_tracker_->Fill( dcurvVal/sqrt( pow(pt1ErrVal,2)/pow(pt1,4) + pow(pt2ErrVal,2)/pow(pt2,4) ) );
00268             
00269             // if do the same for split muons
00270             if (plotMuons_ && splitMuons.isValid()){
00271               
00272               int gmCtr = 0; 
00273               bool topGlobalMuonFlag = false;
00274               bool bottomGlobalMuonFlag = false;
00275               int topGlobalMuon = -1;
00276               int bottomGlobalMuon = -1;
00277               double topGlobalMuonNorchi2 = 1e10;
00278               double bottomGlobalMuonNorchi2 = 1e10;
00279               
00280               // check if usable split global muons
00281               for (std::vector<reco::Muon>::const_iterator gmI = splitMuons->begin(); gmI != splitMuons->end(); gmI++){
00282                 if ( gmI->isTrackerMuon() && gmI->isStandAloneMuon() && gmI->isGlobalMuon() ){
00283                   
00284                   reco::TrackRef trackerTrackRef1( splitTracks, 0 );
00285                   reco::TrackRef trackerTrackRef2( splitTracks, 1 );
00286                   
00287                   if (gmI->innerTrack() == trackerTrackRef1){
00288                     if (gmI->globalTrack()->normalizedChi2() < topGlobalMuonNorchi2){
00289                       topGlobalMuonFlag = true;
00290                       topGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
00291                       topGlobalMuon = gmCtr;
00292                     }
00293                   }
00294                   if (gmI->innerTrack() == trackerTrackRef2){
00295                     if (gmI->globalTrack()->normalizedChi2() < bottomGlobalMuonNorchi2){
00296                       bottomGlobalMuonFlag = true;
00297                       bottomGlobalMuonNorchi2 = gmI->globalTrack()->normalizedChi2();
00298                       bottomGlobalMuon = gmCtr;
00299                     }
00300                   }
00301                 }
00302                 gmCtr++;
00303               } 
00304               
00305               if (bottomGlobalMuonFlag && topGlobalMuonFlag) {
00306                 
00307                 reco::Muon muonTop = splitMuons->at( topGlobalMuon );
00308                 reco::Muon muonBottom = splitMuons->at( bottomGlobalMuon );                            
00309                 
00310                 reco::TrackRef glb1 = muonTop.globalTrack();
00311                 reco::TrackRef glb2 = muonBottom.globalTrack();
00312                 
00313                 double ddxyValGlb = glb1->d0() - glb2->d0();
00314                 double ddzValGlb = glb1->dz() - glb2->dz();                                             
00315                 double dphiValGlb = glb1->phi() - glb2->phi();
00316                 double dthetaValGlb = glb1->theta() - glb2->theta();                                            
00317                 double dptValGlb = glb1->pt() - glb2->pt();
00318                 double dcurvValGlb = (1/glb1->pt()) - (1/glb2->pt());
00319                 
00320                 double d01ErrValGlb = glb1->d0Error();
00321                 double d02ErrValGlb = glb2->d0Error();
00322                 double dz1ErrValGlb = glb1->dzError();
00323                 double dz2ErrValGlb = glb2->dzError();
00324                 double phi1ErrValGlb = glb1->phiError();
00325                 double phi2ErrValGlb = glb2->phiError();
00326                 double theta1ErrValGlb = glb1->thetaError();
00327                 double theta2ErrValGlb = glb2->thetaError();
00328                 double pt1ErrValGlb = glb1->ptError();
00329                 double pt2ErrValGlb = glb2->ptError();
00330                 
00331                 ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddxyValGlb/sqrt(2.0) );
00332                 ddxyAbsoluteResiduals_global_->Fill( 10000.0*ddzValGlb/sqrt(2.0) );
00333                 ddxyAbsoluteResiduals_global_->Fill( 1000.0*dphiValGlb/sqrt(2.0) );
00334                 ddxyAbsoluteResiduals_global_->Fill( 1000.0*dthetaValGlb/sqrt(2.0) );
00335                 ddxyAbsoluteResiduals_global_->Fill( dptValGlb/sqrt(2.0) );
00336                 ddxyAbsoluteResiduals_global_->Fill( dcurvValGlb/sqrt(2.0) );
00337                 
00338                 ddxyNormalizedResiduals_global_->Fill( ddxyValGlb/sqrt( d01ErrValGlb*d01ErrValGlb + d02ErrValGlb*d02ErrValGlb ) );
00339                 ddxyNormalizedResiduals_global_->Fill( ddzValGlb/sqrt( dz1ErrValGlb*dz1ErrValGlb + dz2ErrValGlb*dz2ErrValGlb ) );
00340                 ddxyNormalizedResiduals_global_->Fill( dphiValGlb/sqrt( phi1ErrValGlb*phi1ErrValGlb + phi2ErrValGlb*phi2ErrValGlb ) );
00341                 ddxyNormalizedResiduals_global_->Fill( dthetaValGlb/sqrt( theta1ErrValGlb*theta1ErrValGlb + theta2ErrValGlb*theta2ErrValGlb ) );
00342                 ddxyNormalizedResiduals_global_->Fill( dptValGlb/sqrt( pt1ErrValGlb*pt1ErrValGlb + pt2ErrValGlb*pt2ErrValGlb ) );
00343                 ddxyNormalizedResiduals_global_->Fill( dcurvValGlb/sqrt( pow(pt1ErrValGlb,2)/pow(pt1,4) + pow(pt2ErrValGlb,2)/pow(pt2,4) ) );
00344                 
00345               }
00346               
00347               
00348             } // end of split muons loop
00349           }
00350         }
00351       }
00352     }    
00353   }
00354 }
00355 
00356 
00357 void TrackSplittingMonitor::endJob(void) {
00358   bool outputMEsInRootFile = conf_.getParameter<bool>("OutputMEsInRootFile");
00359   std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
00360   if(outputMEsInRootFile){
00361     dqmStore_->showDirStructure();
00362     dqmStore_->save(outputFileName);
00363   }
00364 }
00365 
00366 DEFINE_FWK_MODULE(TrackSplittingMonitor);