CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/DQMOffline/RecoB/plugins/PrimaryVertexMonitor.cc

Go to the documentation of this file.
00001 #include "DQMOffline/RecoB/plugins/PrimaryVertexMonitor.h"
00002 #include "FWCore/ServiceRegistry/interface/Service.h"
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 
00005 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00006 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00007 
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"     
00009 
00010 using namespace reco;
00011 using namespace edm;
00012 
00013 PrimaryVertexMonitor::PrimaryVertexMonitor(const edm::ParameterSet& pSet)
00014 {
00015   moduleLabel = pSet.getParameter<InputTag>("vertexLabel");
00016   beamSpotLabel = pSet.getParameter<InputTag>("beamSpotLabel");
00017 
00018   //
00019   // Book all histograms.
00020   //
00021 
00022   //  get the store
00023   dqmStore_ = edm::Service<DQMStore>().operator->();
00024   dqmLabel = "OfflinePV/"+moduleLabel.label();
00025   dqmStore_->setCurrentFolder(dqmLabel);
00026 
00027 //   xPos = dqmStore_->book1D ("xPos","x Coordinate" ,100, -0.1, 0.1);
00028 
00029   nbvtx      = dqmStore_->book1D("vtxNbr","Reconstructed Vertices in Event",20,-0.5,19.5);
00030 
00031   nbtksinvtx[0] = dqmStore_->book1D("otherVtxTrksNbr","Reconstructed Tracks in Vertex (other Vtx)",40,-0.5,99.5); 
00032   trksWeight[0] = dqmStore_->book1D("otherVtxTrksWeight","Total weight of Tracks in Vertex (other Vtx)",40,0,100.); 
00033   vtxchi2[0]    = dqmStore_->book1D("otherVtxChi2","#chi^{2} (other Vtx)",100,0.,200.);
00034   vtxndf[0]     = dqmStore_->book1D("otherVtxNdf","ndof (other Vtx)",100,0.,200.);
00035   vtxprob[0]    = dqmStore_->book1D("otherVtxProb","#chi^{2} probability (other Vtx)",100,0.,1.);
00036   nans[0]       = dqmStore_->book1D("otherVtxNans","Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (other Vtx)",9,0.5,9.5);
00037 
00038   nbtksinvtx[1] = dqmStore_->book1D("tagVtxTrksNbr","Reconstructed Tracks in Vertex (tagged Vtx)",100,-0.5,99.5); 
00039   trksWeight[1] = dqmStore_->book1D("tagVtxTrksWeight","Total weight of Tracks in Vertex (tagged Vtx)",100,0,100.); 
00040   vtxchi2[1]    = dqmStore_->book1D("tagVtxChi2","#chi^{2} (tagged Vtx)",100,0.,200.);
00041   vtxndf[1]     = dqmStore_->book1D("tagVtxNdf","ndof (tagged Vtx)",100,0.,200.);
00042   vtxprob[1]    = dqmStore_->book1D("tagVtxProb","#chi^{2} probability (tagged Vtx)",100,0.,1.);
00043   nans[1]       = dqmStore_->book1D("tagVtxNans","Illegal values for x,y,z,xx,xy,xz,yy,yz,zz (tagged Vtx)",9,0.5,9.5);
00044 
00045   xrec[0]        = dqmStore_->book1D("otherPosX","Position x Coordinate (other Vtx)",100,-0.1,0.1);
00046   yrec[0]        = dqmStore_->book1D("otherPosY","Position y Coordinate (other Vtx)",100,-0.1,0.1);
00047   zrec[0]        = dqmStore_->book1D("otherPosZ","Position z Coordinate (other Vtx)",100,-20.,20.);
00048   xDiff[0]       = dqmStore_->book1D("otherDiffX","X distance from BeamSpot (other Vtx)",100,-500,500);
00049   yDiff[0]       = dqmStore_->book1D("otherDiffY","Y distance from BeamSpot (other Vtx)",100,-500,500);
00050   xerr[0]        = dqmStore_->book1D("otherErrX","Uncertainty x Coordinate (other Vtx)",100,-0.1,0.1);
00051   yerr[0]        = dqmStore_->book1D("otherErrY","Uncertainty y Coordinate (other Vtx)",100,-0.1,0.1);
00052   zerr[0]        = dqmStore_->book1D("otherErrZ","Uncertainty z Coordinate (other Vtx)",100,-20.,20.);
00053   xerrVsTrks[0]  = dqmStore_->book2D("otherErrVsWeightX","Uncertainty x Coordinate vs. track weight (other Vtx)",100,0,100.,100,-0.1,0.1);
00054   yerrVsTrks[0]  = dqmStore_->book2D("otherErrVsWeightY","Uncertainty y Coordinate vs. track weight (other Vtx)",100,0,100.,100,-0.1,0.1);
00055   zerrVsTrks[0]  = dqmStore_->book2D("otherErrVsWeightZ","Uncertainty z Coordinate vs. track weight (other Vtx)",100,0,100.,100,-0.1,0.1);
00056 
00057 
00058   xrec[1]     = dqmStore_->book1D("tagPosX","Position x Coordinate (tagged Vtx)",100,-0.1,0.1);
00059   yrec[1]     = dqmStore_->book1D("tagPosY","Position y Coordinate (tagged Vtx)",100,-0.1,0.1);
00060   zrec[1]     = dqmStore_->book1D("tagPosZ","Position z Coordinate (tagged Vtx)",100,-20.,20.);
00061   xDiff[1]    = dqmStore_->book1D("tagDiffX","X distance from BeamSpot (tagged Vtx)",100,-500, 500);
00062   yDiff[1]    = dqmStore_->book1D("tagDiffY","Y distance from BeamSpot (tagged Vtx)",100,-500, 500);
00063   xerr[1]     = dqmStore_->book1D("tagErrX","Uncertainty x Coordinate (tagged Vtx)",100,0.,100);
00064   yerr[1]     = dqmStore_->book1D("tagErrY","Uncertainty y Coordinate (tagged Vtx)",100,0.,100);
00065   zerr[1]     = dqmStore_->book1D("tagErrZ","Uncertainty z Coordinate (tagged Vtx)",100,0.,100);
00066   xerrVsTrks[1]  = dqmStore_->book2D("tagErrVsWeightX","Uncertainty x Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
00067   yerrVsTrks[1]  = dqmStore_->book2D("tagErrVsWeightY","Uncertainty y Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
00068   zerrVsTrks[1]  = dqmStore_->book2D("tagErrVsWeightZ","Uncertainty z Coordinate vs. track weight (tagged Vtx)",100,0,100.,100,0.,100);
00069 
00070   type[0] = dqmStore_->book1D("otherType","Vertex type (other Vtx)",3,-0.5,2.5);
00071   type[1] = dqmStore_->book1D("tagType","Vertex type (tagged Vtx)",3,-0.5,2.5);
00072   for (int i=0;i<2;++i){
00073     type[i]->getTH1F()->GetXaxis()->SetBinLabel(1,"Valid, real");
00074     type[i]->getTH1F()->GetXaxis()->SetBinLabel(2,"Valid, fake");
00075     type[i]->getTH1F()->GetXaxis()->SetBinLabel(3,"Invalid");
00076   }
00077 
00078   bsX           = dqmStore_->book1D("bsX", "BeamSpot x0", 100,-0.1,0.1);
00079   bsY           = dqmStore_->book1D("bsY", "BeamSpot y0", 100,-0.1,0.1);
00080   bsZ           = dqmStore_->book1D("bsZ", "BeamSpot z0", 100,-2.,2.);
00081   bsSigmaZ      = dqmStore_->book1D("bsSigmaZ", "BeamSpot sigmaZ", 100, 0., 10. );
00082   bsDxdz        = dqmStore_->book1D("bsDxdz", "BeamSpot dxdz", 100, -0.0003, 0.0003);
00083   bsDydz        = dqmStore_->book1D("bsDydz", "BeamSpot dydz", 100, -0.0003, 0.0003);
00084   bsBeamWidthX  = dqmStore_->book1D("bsBeamWidthX", "BeamSpot BeamWidthX", 100, 0., 100.);
00085   bsBeamWidthY  = dqmStore_->book1D("bsBeamWidthY", "BeamSpot BeamWidthY", 100, 0., 100.);
00086   bsType        = dqmStore_->book1D("bsType", "BeamSpot type", 4, -1.5, 2.5);
00087   bsType->getTH1F()->GetXaxis()->SetBinLabel(1,"Unknown");
00088   bsType->getTH1F()->GetXaxis()->SetBinLabel(2,"Fake");
00089   bsType->getTH1F()->GetXaxis()->SetBinLabel(3,"LHC");
00090   bsType->getTH1F()->GetXaxis()->SetBinLabel(4,"Tracker");
00091 }
00092 
00093 
00094 PrimaryVertexMonitor::~PrimaryVertexMonitor()
00095 {}
00096 
00097 void PrimaryVertexMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00098 {
00099   Handle<reco::VertexCollection> recVtxs;
00100   iEvent.getByLabel(moduleLabel, recVtxs);
00101 
00102   edm::Handle<reco::BeamSpot> beamSpotHandle;
00103   iEvent.getByLabel(beamSpotLabel,beamSpotHandle);
00104 
00105   //
00106   // check for absent products and simply "return" in that case
00107   //
00108   if (recVtxs.isValid() == false || beamSpotHandle.isValid()== false){
00109     edm::LogWarning("PrimaryVertexMonitor")
00110       <<" Some products not available in the event: VertexCollection "
00111       <<moduleLabel<<" " 
00112       <<recVtxs.isValid() <<" BeamSpot "
00113       <<beamSpotLabel<<" "
00114       <<beamSpotHandle.isValid()<<". Skipping plots for this event";
00115     return;
00116   }
00117 
00118   BeamSpot beamSpot = *beamSpotHandle;
00119 
00120   nbvtx->Fill(recVtxs->size()*1.);
00121 
00122   vertexPlots(recVtxs->front(), beamSpot, 1);
00123 
00124   for(reco::VertexCollection::const_iterator v=recVtxs->begin()+1; 
00125       v!=recVtxs->end(); ++v){
00126     vertexPlots(*v, beamSpot, 0);
00127   }
00128   // Beamline plots:
00129   bsX->Fill(beamSpot.x0());
00130   bsY->Fill(beamSpot.y0());
00131   bsZ->Fill(beamSpot.z0());
00132   bsSigmaZ->Fill(beamSpot.sigmaZ());
00133   bsDxdz->Fill(beamSpot.dxdz());
00134   bsDydz->Fill(beamSpot.dydz());
00135   bsBeamWidthX->Fill(beamSpot.BeamWidthX()*10000);
00136   bsBeamWidthY->Fill(beamSpot.BeamWidthY()*10000);
00137   // bsType->Fill(beamSpot.type());
00138 
00139 }
00140 
00141 void PrimaryVertexMonitor::vertexPlots(const Vertex & v, const BeamSpot& beamSpot, int i)
00142 {
00143     float weight = 0;
00144     for(reco::Vertex::trackRef_iterator t = v.tracks_begin(); 
00145         t!=v.tracks_end(); t++) weight+= v.trackWeight(*t);
00146     trksWeight[i]->Fill(weight);
00147     nbtksinvtx[i]->Fill(v.tracksSize());
00148 
00149     vtxchi2[i]->Fill(v.chi2());
00150     vtxndf[i]->Fill(v.ndof());
00151     vtxprob[i]->Fill(ChiSquaredProbability(v.chi2() ,v.ndof()));
00152 
00153     if (!v.isValid()) type[i]->Fill(2.);
00154     else if (v.isFake()) type[i]->Fill(1.);
00155     else type[i]->Fill(0.);
00156 
00157     xrec[i]->Fill(v.position().x());
00158     yrec[i]->Fill(v.position().y());
00159     zrec[i]->Fill(v.position().z());
00160 
00161     float xb = beamSpot.x0() + beamSpot.dxdz() * (v.position().z() - beamSpot.z0());
00162     float yb = beamSpot.y0() + beamSpot.dydz() * (v.position().z() - beamSpot.z0());
00163     xDiff[i]->Fill((v.position().x() - xb)*10000);
00164     yDiff[i]->Fill((v.position().y() - yb)*10000);
00165 
00166     xerr[i]->Fill(v.xError()*10000);
00167     yerr[i]->Fill(v.yError()*10000);
00168     zerr[i]->Fill(v.zError()*10000);
00169     xerrVsTrks[i]->Fill(weight, v.xError()*10000);
00170     yerrVsTrks[i]->Fill(weight, v.yError()*10000);
00171     zerrVsTrks[i]->Fill(weight, v.zError()*10000);
00172 
00173     nans[i]->Fill(1.,std::isnan(v.position().x())*1.);
00174     nans[i]->Fill(2.,std::isnan(v.position().y())*1.);
00175     nans[i]->Fill(3.,std::isnan(v.position().z())*1.);
00176 
00177     int index = 3;
00178     for (int i = 0; i != 3; i++) {
00179       for (int j = i; j != 3; j++) {
00180         index++;
00181         nans[i]->Fill(index*1., std::isnan(v.covariance(i, j))*1.);
00182         // in addition, diagonal element must be positive
00183         if (j == i && v.covariance(i, j) < 0) {
00184           nans[i]->Fill(index*1., 1.);
00185         }
00186       }
00187     }
00188 }
00189 
00190 
00191 void PrimaryVertexMonitor::endJob()
00192 {
00193 }
00194 
00195 
00196 //define this as a plug-in
00197 DEFINE_FWK_MODULE(PrimaryVertexMonitor);