Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:16

0001 // -*- C++ -*-
0002 // Package:    SiPixelMonitorTrackSoA
0003 // Class:      SiPixelMonitorTrackSoA
0004 //
0005 /**\class SiPixelMonitorTrackSoA SiPixelMonitorTrackSoA.cc
0006 */
0007 //
0008 // Author: Suvankar Roy Chowdhury
0009 //
0010 #include "DataFormats/Common/interface/Handle.h"
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/Frameworkfwd.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 // DQM Histograming
0020 #include "DQMServices/Core/interface/MonitorElement.h"
0021 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0022 #include "DQMServices/Core/interface/DQMStore.h"
0023 #include "CUDADataFormats/Track/interface/PixelTrackUtilities.h"
0024 #include "CUDADataFormats/Track/interface/TrackSoAHeterogeneousHost.h"
0025 // for string manipulations
0026 #include <fmt/printf.h>
0027 
0028 template <typename T>
0029 class SiPixelMonitorTrackSoA : public DQMEDAnalyzer {
0030 public:
0031   using PixelTrackHeterogeneous = TrackSoAHeterogeneousHost<T>;
0032   explicit SiPixelMonitorTrackSoA(const edm::ParameterSet&);
0033   ~SiPixelMonitorTrackSoA() override = default;
0034   void bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) override;
0035   void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0036   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0037 
0038 private:
0039   edm::EDGetTokenT<PixelTrackHeterogeneous> tokenSoATrack_;
0040   std::string topFolderName_;
0041   bool useQualityCut_;
0042   pixelTrack::Quality minQuality_;
0043   MonitorElement* hnTracks;
0044   MonitorElement* hnLooseAndAboveTracks;
0045   MonitorElement* hnHits;
0046   MonitorElement* hnHitsVsPhi;
0047   MonitorElement* hnHitsVsEta;
0048   MonitorElement* hnLayers;
0049   MonitorElement* hnLayersVsPhi;
0050   MonitorElement* hnLayersVsEta;
0051   MonitorElement* hchi2;
0052   MonitorElement* hChi2VsPhi;
0053   MonitorElement* hChi2VsEta;
0054   MonitorElement* hpt;
0055   MonitorElement* heta;
0056   MonitorElement* hphi;
0057   MonitorElement* hz;
0058   MonitorElement* htip;
0059   MonitorElement* hquality;
0060 };
0061 
0062 //
0063 // constructors
0064 //
0065 
0066 template <typename T>
0067 SiPixelMonitorTrackSoA<T>::SiPixelMonitorTrackSoA(const edm::ParameterSet& iConfig) {
0068   tokenSoATrack_ = consumes<PixelTrackHeterogeneous>(iConfig.getParameter<edm::InputTag>("pixelTrackSrc"));
0069   topFolderName_ = iConfig.getParameter<std::string>("topFolderName");  //"SiPixelHeterogeneous/PixelTrackSoA";
0070   useQualityCut_ = iConfig.getParameter<bool>("useQualityCut");
0071   minQuality_ = pixelTrack::qualityByName(iConfig.getParameter<std::string>("minQuality"));
0072 }
0073 
0074 //
0075 // -- Analyze
0076 //
0077 template <typename T>
0078 void SiPixelMonitorTrackSoA<T>::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0079   const auto& tsoaHandle = iEvent.getHandle(tokenSoATrack_);
0080   if (!tsoaHandle.isValid()) {
0081     edm::LogWarning("SiPixelMonitorTrackSoA") << "No Track SoA found \n returning!" << std::endl;
0082     return;
0083   }
0084 
0085   using helper = TracksUtilities<T>;
0086   auto const& tsoa = *tsoaHandle.product();
0087   auto maxTracks = tsoa.view().metadata().size();
0088   auto const* quality = tsoa.view().quality();
0089   int32_t nTracks = 0;
0090   int32_t nLooseAndAboveTracks = 0;
0091 
0092   for (int32_t it = 0; it < maxTracks; ++it) {
0093     auto nHits = helper::nHits(tsoa.const_view(), it);
0094     auto nLayers = tsoa.view()[it].nLayers();
0095     if (nHits == 0)
0096       break;  // this is a guard
0097     float pt = tsoa.view()[it].pt();
0098     if (!(pt > 0.))
0099       continue;
0100 
0101     // fill the quality for all tracks
0102     pixelTrack::Quality qual = quality[it];
0103     hquality->Fill(int(qual));
0104     nTracks++;
0105 
0106     if (useQualityCut_ && quality[it] < minQuality_)
0107       continue;
0108 
0109     // fill parameters only for quality >= loose
0110     float chi2 = tsoa.view()[it].chi2();
0111     float phi = helper::phi(tsoa.const_view(), it);
0112     float zip = helper::zip(tsoa.const_view(), it);
0113     float eta = tsoa.view()[it].eta();
0114     float tip = helper::tip(tsoa.const_view(), it);
0115 
0116     hchi2->Fill(chi2);
0117     hChi2VsPhi->Fill(phi, chi2);
0118     hChi2VsEta->Fill(eta, chi2);
0119     hnHits->Fill(nHits);
0120     hnLayers->Fill(nLayers);
0121     hnHitsVsPhi->Fill(phi, nHits);
0122     hnHitsVsEta->Fill(eta, nHits);
0123     hnLayersVsPhi->Fill(phi, nLayers);
0124     hnLayersVsEta->Fill(eta, nLayers);
0125     hpt->Fill(pt);
0126     heta->Fill(eta);
0127     hphi->Fill(phi);
0128     hz->Fill(zip);
0129     htip->Fill(tip);
0130     nLooseAndAboveTracks++;
0131   }
0132   hnTracks->Fill(nTracks);
0133   hnLooseAndAboveTracks->Fill(nLooseAndAboveTracks);
0134 }
0135 
0136 //
0137 // -- Book Histograms
0138 //
0139 template <typename T>
0140 void SiPixelMonitorTrackSoA<T>::bookHistograms(DQMStore::IBooker& iBook,
0141                                                edm::Run const& iRun,
0142                                                edm::EventSetup const& iSetup) {
0143   iBook.cd();
0144   iBook.setCurrentFolder(topFolderName_);
0145 
0146   // clang-format off
0147   std::string toRep = "Number of tracks";
0148   hnTracks = iBook.book1D("nTracks", fmt::sprintf(";%s per event;#events",toRep), 1001, -0.5, 1000.5);
0149   hnLooseAndAboveTracks = iBook.book1D("nLooseAndAboveTracks", fmt::sprintf(";%s (quality #geq loose) per event;#events",toRep), 1001, -0.5, 1000.5);
0150 
0151   toRep = "Number of all RecHits per track (quality #geq loose)";
0152   hnHits = iBook.book1D("nRecHits", fmt::sprintf(";%s;#tracks",toRep), 15, -0.5, 14.5);
0153   hnHitsVsPhi = iBook.bookProfile("nHitsPerTrackVsPhi", fmt::sprintf("%s vs track #phi;Track #phi;%s",toRep,toRep), 30, -M_PI, M_PI,0., 15.);
0154   hnHitsVsEta = iBook.bookProfile("nHitsPerTrackVsEta", fmt::sprintf("%s vs track #eta;Track #eta;%s",toRep,toRep), 30, -3., 3., 0., 15.);
0155 
0156   toRep = "Number of all layers per track (quality #geq loose)";
0157   hnLayers = iBook.book1D("nLayers", fmt::sprintf(";%s;#tracks",toRep), 15, -0.5, 14.5);
0158   hnLayersVsPhi = iBook.bookProfile("nLayersPerTrackVsPhi", fmt::sprintf("%s vs track #phi;Track #phi;%s",toRep,toRep), 30, -M_PI, M_PI,0., 15.);
0159   hnLayersVsEta = iBook.bookProfile("nLayersPerTrackVsEta", fmt::sprintf("%s vs track #eta;Track #eta;%s",toRep,toRep), 30, -3., 3., 0., 15.);
0160 
0161   toRep = "Track (quality #geq loose) #chi^{2}/ndof";
0162   hchi2 = iBook.book1D("nChi2ndof", fmt::sprintf(";%s;#tracks",toRep), 40, 0., 20.);
0163   hChi2VsPhi = iBook.bookProfile("nChi2ndofVsPhi", fmt::sprintf("%s vs track #phi;Track #phi;%s",toRep,toRep), 30, -M_PI, M_PI, 0., 20.);
0164   hChi2VsEta = iBook.bookProfile("nChi2ndofVsEta", fmt::sprintf("%s vs track #eta;Track #eta;%s",toRep,toRep), 30, -3., 3., 0., 20.);
0165   // clang-format on
0166 
0167   hpt = iBook.book1D("pt", ";Track (quality #geq loose) p_{T} [GeV];#tracks", 200, 0., 200.);
0168   heta = iBook.book1D("eta", ";Track (quality #geq loose) #eta;#tracks", 30, -3., 3.);
0169   hphi = iBook.book1D("phi", ";Track (quality #geq loose) #phi;#tracks", 30, -M_PI, M_PI);
0170   hz = iBook.book1D("z", ";Track (quality #geq loose) z [cm];#tracks", 30, -30., 30.);
0171   htip = iBook.book1D("tip", ";Track (quality #geq loose) TIP [cm];#tracks", 100, -0.5, 0.5);
0172   hquality = iBook.book1D("quality", ";Track Quality;#tracks", 7, -0.5, 6.5);
0173   uint i = 1;
0174   for (const auto& q : pixelTrack::qualityName) {
0175     hquality->setBinLabel(i, q, 1);
0176     i++;
0177   }
0178 }
0179 
0180 template <typename T>
0181 void SiPixelMonitorTrackSoA<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0182   // monitorpixelTrackSoA
0183   edm::ParameterSetDescription desc;
0184   desc.add<edm::InputTag>("pixelTrackSrc", edm::InputTag("pixelTracksSoA"));
0185   desc.add<std::string>("topFolderName", "SiPixelHeterogeneous/PixelTrackSoA");
0186   desc.add<bool>("useQualityCut", true);
0187   desc.add<std::string>("minQuality", "loose");
0188   descriptions.addWithDefaultLabel(desc);
0189 }
0190 
0191 using SiPixelPhase1MonitorTrackSoA = SiPixelMonitorTrackSoA<pixelTopology::Phase1>;
0192 using SiPixelPhase2MonitorTrackSoA = SiPixelMonitorTrackSoA<pixelTopology::Phase2>;
0193 using SiPixelHIonPhase1MonitorTrackSoA = SiPixelMonitorTrackSoA<pixelTopology::HIonPhase1>;
0194 
0195 DEFINE_FWK_MODULE(SiPixelPhase1MonitorTrackSoA);
0196 DEFINE_FWK_MODULE(SiPixelPhase2MonitorTrackSoA);
0197 DEFINE_FWK_MODULE(SiPixelHIonPhase1MonitorTrackSoA);