CMS 3D CMS Logo

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