File indexing completed on 2024-04-06 12:08:16
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "FWCore/Framework/interface/Frameworkfwd.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018 #include "FWCore/Utilities/interface/InputTag.h"
0019 #include "DataFormats/Common/interface/Handle.h"
0020
0021 #include "DQMServices/Core/interface/MonitorElement.h"
0022 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0023 #include "DQMServices/Core/interface/DQMStore.h"
0024 #include "CUDADataFormats/Vertex/interface/ZVertexSoAHeterogeneousHost.h"
0025 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0026
0027 class SiPixelMonitorVertexSoA : public DQMEDAnalyzer {
0028 public:
0029 using IndToEdm = std::vector<uint16_t>;
0030 explicit SiPixelMonitorVertexSoA(const edm::ParameterSet&);
0031 ~SiPixelMonitorVertexSoA() override = default;
0032 void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0033 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0034 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035
0036 private:
0037 edm::EDGetTokenT<ZVertexSoAHost> tokenSoAVertex_;
0038 edm::EDGetTokenT<reco::BeamSpot> tokenBeamSpot_;
0039 std::string topFolderName_;
0040 MonitorElement* hnVertex;
0041 MonitorElement* hx;
0042 MonitorElement* hy;
0043 MonitorElement* hz;
0044 MonitorElement* hchi2;
0045 MonitorElement* hchi2oNdof;
0046 MonitorElement* hptv2;
0047 MonitorElement* hntrks;
0048 };
0049
0050
0051
0052
0053
0054 SiPixelMonitorVertexSoA::SiPixelMonitorVertexSoA(const edm::ParameterSet& iConfig) {
0055 tokenSoAVertex_ = consumes(iConfig.getParameter<edm::InputTag>("pixelVertexSrc"));
0056 tokenBeamSpot_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"));
0057 topFolderName_ = iConfig.getParameter<std::string>("topFolderName");
0058 }
0059
0060
0061
0062
0063 void SiPixelMonitorVertexSoA::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0064 const auto& vsoaHandle = iEvent.getHandle(tokenSoAVertex_);
0065 if (!vsoaHandle.isValid()) {
0066 edm::LogWarning("SiPixelMonitorVertexSoA") << "No Vertex SoA found \n returning!" << std::endl;
0067 return;
0068 }
0069
0070 auto const& vsoa = *vsoaHandle;
0071 int nVertices = vsoa.view().nvFinal();
0072 auto bsHandle = iEvent.getHandle(tokenBeamSpot_);
0073 float x0 = 0., y0 = 0., z0 = 0., dxdz = 0., dydz = 0.;
0074 if (!bsHandle.isValid()) {
0075 edm::LogWarning("SiPixelMonitorVertexSoA") << "No beamspot found. returning vertexes with (0,0,Z) ";
0076 } else {
0077 const reco::BeamSpot& bs = *bsHandle;
0078 x0 = bs.x0();
0079 y0 = bs.y0();
0080 z0 = bs.z0();
0081 dxdz = bs.dxdz();
0082 dydz = bs.dydz();
0083 }
0084
0085 for (int iv = 0; iv < nVertices; iv++) {
0086 auto si = vsoa.view()[iv].sortInd();
0087 auto z = vsoa.view()[si].zv();
0088 auto x = x0 + dxdz * z;
0089 auto y = y0 + dydz * z;
0090
0091 z += z0;
0092 hx->Fill(x);
0093 hy->Fill(y);
0094 hz->Fill(z);
0095 auto ndof = vsoa.view()[si].ndof();
0096 hchi2->Fill(vsoa.view()[si].chi2());
0097 hchi2oNdof->Fill(vsoa.view()[si].chi2() / ndof);
0098 hptv2->Fill(vsoa.view()[si].ptv2());
0099 hntrks->Fill(ndof + 1);
0100 }
0101 hnVertex->Fill(nVertices);
0102 }
0103
0104
0105
0106
0107 void SiPixelMonitorVertexSoA::bookHistograms(DQMStore::IBooker& ibooker,
0108 edm::Run const& iRun,
0109 edm::EventSetup const& iSetup) {
0110
0111 ibooker.cd();
0112 ibooker.setCurrentFolder(topFolderName_);
0113 hnVertex = ibooker.book1D("nVertex", ";# of Vertices;#entries", 101, -0.5, 100.5);
0114 hx = ibooker.book1D("vx", ";Vertex x;#entries", 10, -5., 5.);
0115 hy = ibooker.book1D("vy", ";Vertex y;#entries", 10, -5., 5.);
0116 hz = ibooker.book1D("vz", ";Vertex z;#entries", 30, -30., 30);
0117 hchi2 = ibooker.book1D("chi2", ";Vertex chi-squared;#entries", 40, 0., 20.);
0118 hchi2oNdof = ibooker.book1D("chi2oNdof", ";Vertex chi-squared/Ndof;#entries", 40, 0., 20.);
0119 hptv2 = ibooker.book1D("ptsq", ";Vertex #sum (p_{T})^{2};#entries", 200, 0., 200.);
0120 hntrks = ibooker.book1D("ntrk", ";#tracks associated;#entries", 100, -0.5, 99.5);
0121 }
0122
0123 void SiPixelMonitorVertexSoA::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0124
0125 edm::ParameterSetDescription desc;
0126 desc.add<edm::InputTag>("pixelVertexSrc", edm::InputTag("pixelVerticesSoA"));
0127 desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
0128 desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelVertexSoA");
0129 descriptions.addWithDefaultLabel(desc);
0130 }
0131
0132 DEFINE_FWK_MODULE(SiPixelMonitorVertexSoA);