Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:54:05

0001 /****************************************************************************
0002  *
0003  * This is a part of TOTEM offline software.
0004  * Author:
0005  *   Laurent Forthomme
0006  *   Arkadiusz Cwikla
0007  *
0008  ****************************************************************************/
0009 
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 
0017 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0018 #include "DQMServices/Core/interface/DQMStore.h"
0019 
0020 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0021 
0022 #include "DataFormats/CTPPSDetId/interface/TotemT2DetId.h"
0023 #include "DataFormats/TotemReco/interface/TotemT2Digi.h"
0024 #include "DataFormats/TotemReco/interface/TotemT2RecHit.h"
0025 
0026 #include "Geometry/Records/interface/TotemGeometryRcd.h"
0027 #include "DQM/CTPPS/interface/TotemT2Segmentation.h"
0028 
0029 #include <string>
0030 
0031 class TotemT2DQMSource : public DQMEDAnalyzer {
0032 public:
0033   TotemT2DQMSource(const edm::ParameterSet&);
0034   ~TotemT2DQMSource() override = default;
0035 
0036 protected:
0037   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0038   void bookHistograms(DQMStore::IBooker&, const edm::Run&, const edm::EventSetup&) override;
0039   void analyze(const edm::Event&, const edm::EventSetup&) override;
0040 
0041 private:
0042   void fillActivePlanes(std::unordered_map<unsigned int, std::set<unsigned int>>&, const TotemT2DetId&);
0043   void fillTriggerBitset(const TotemT2DetId&);
0044   void clearTriggerBitset();
0045   bool areChannelsTriggered(const TotemT2DetId&);
0046   void bookErrorFlagsHistogram(DQMStore::IBooker&);
0047   void fillErrorFlagsHistogram(const TotemT2Digi&);
0048   void fillEdges(const TotemT2Digi&, const TotemT2DetId&);
0049   void fillToT(const TotemT2RecHit&, const TotemT2DetId&);
0050 
0051   const edm::ESGetToken<TotemGeometry, TotemGeometryRcd> geometryToken_;
0052   const edm::EDGetTokenT<edmNew::DetSetVector<TotemT2Digi>> digiToken_;
0053   const edm::EDGetTokenT<edmNew::DetSetVector<TotemT2RecHit>> rechitToken_;
0054 
0055   std::unique_ptr<TotemT2Segmentation> segm_;
0056 
0057   static constexpr double HPTDC_BIN_WIDTH_NS_ = 25. / 1024;
0058   MonitorElement* HPTDCErrorFlags_2D_ = nullptr;
0059 
0060   const unsigned int nbinsx_, nbinsy_;
0061   const unsigned int windowsNum_;
0062 
0063   struct SectorPlots {
0064     MonitorElement* activePlanes = nullptr;
0065     MonitorElement* activePlanesCount = nullptr;
0066 
0067     MonitorElement* triggerEmulator = nullptr;
0068     std::bitset<(TotemT2DetId::maxPlane + 1) * (TotemT2DetId::maxChannel + 1)> hitTilesArray;
0069     static const unsigned int MINIMAL_TRIGGER = 3;
0070 
0071     MonitorElement *leadingEdge = nullptr, *trailingEdge = nullptr, *timeOverTreshold = nullptr;
0072 
0073     SectorPlots() = default;
0074     SectorPlots(
0075         DQMStore::IBooker& ibooker, unsigned int id, unsigned int nbinsx, unsigned int nbinsy, unsigned int windowsNum);
0076   };
0077 
0078   struct PlanePlots {
0079     MonitorElement* digisMultiplicity = nullptr;
0080     MonitorElement* rechitMultiplicity = nullptr;
0081 
0082     PlanePlots() = default;
0083     PlanePlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int nbinsx, unsigned int nbinsy);
0084   };
0085 
0086   std::unordered_map<unsigned int, SectorPlots> sectorPlots_;
0087   std::unordered_map<unsigned int, PlanePlots> planePlots_;
0088 };
0089 
0090 TotemT2DQMSource::SectorPlots::SectorPlots(
0091     DQMStore::IBooker& ibooker, unsigned int id, unsigned int nbinsx, unsigned int nbinsy, unsigned int windowsNum) {
0092   std::string title, path;
0093 
0094   TotemT2DetId(id).armName(path, TotemT2DetId::nPath);
0095   ibooker.setCurrentFolder(path);
0096 
0097   TotemT2DetId(id).armName(title, TotemT2DetId::nFull);
0098 
0099   activePlanes = ibooker.book1D("active planes", title + " which planes are active;plane number", 8, -0.5, 7.5);
0100 
0101   activePlanesCount = ibooker.book1D(
0102       "number of active planes", title + " how many planes are active;number of active planes", 9, -0.5, 8.5);
0103 
0104   triggerEmulator = ibooker.book2DD("trigger emulator",
0105                                     title + " trigger emulator;arbitrary unit;arbitrary unit",
0106                                     nbinsx,
0107                                     -0.5,
0108                                     double(nbinsx) - 0.5,
0109                                     nbinsy,
0110                                     -0.5,
0111                                     double(nbinsy) - 0.5);
0112   leadingEdge = ibooker.book1D(
0113       "leading edge", title + " leading edge (DIGIs); leading edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0114   trailingEdge = ibooker.book1D(
0115       "trailing edge", title + " trailing edge (DIGIs); trailing edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0116 
0117   timeOverTreshold = ibooker.book1D(
0118       "time over threshold", title + " time over threshold (rechit);time over threshold (ns)", 250, -25, 100);
0119 }
0120 
0121 TotemT2DQMSource::PlanePlots::PlanePlots(DQMStore::IBooker& ibooker,
0122                                          unsigned int id,
0123                                          unsigned int nbinsx,
0124                                          unsigned int nbinsy) {
0125   std::string title, path;
0126   TotemT2DetId(id).planeName(title, TotemT2DetId::nFull);
0127   TotemT2DetId(id).planeName(path, TotemT2DetId::nPath);
0128   ibooker.setCurrentFolder(path);
0129 
0130   digisMultiplicity = ibooker.book2DD("digis multiplicity",
0131                                       title + " digis multiplicity;arbitrary unit;arbitrary unit",
0132                                       nbinsx,
0133                                       -0.5,
0134                                       double(nbinsx) - 0.5,
0135                                       nbinsy,
0136                                       -0.5,
0137                                       double(nbinsy) - 0.5);
0138   rechitMultiplicity = ibooker.book2DD("rechits multiplicity",
0139                                        title + " rechits multiplicity;x;y",
0140                                        nbinsx,
0141                                        -0.5,
0142                                        double(nbinsx) - 0.5,
0143                                        nbinsy,
0144                                        -0.5,
0145                                        double(nbinsy) - 0.5);
0146 }
0147 
0148 TotemT2DQMSource::TotemT2DQMSource(const edm::ParameterSet& iConfig)
0149     : geometryToken_(esConsumes<TotemGeometry, TotemGeometryRcd, edm::Transition::BeginRun>()),
0150       digiToken_(consumes<edmNew::DetSetVector<TotemT2Digi>>(iConfig.getParameter<edm::InputTag>("digisTag"))),
0151       rechitToken_(consumes<edmNew::DetSetVector<TotemT2RecHit>>(iConfig.getParameter<edm::InputTag>("rechitsTag"))),
0152       nbinsx_(iConfig.getParameter<unsigned int>("nbinsx")),
0153       nbinsy_(iConfig.getParameter<unsigned int>("nbinsy")),
0154       windowsNum_(iConfig.getParameter<unsigned int>("windowsNum")) {}
0155 
0156 void TotemT2DQMSource::dqmBeginRun(const edm::Run&, const edm::EventSetup&) {}
0157 
0158 void TotemT2DQMSource::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup& iSetup) {
0159   ibooker.cd();
0160   ibooker.setCurrentFolder("TotemT2");
0161 
0162   bookErrorFlagsHistogram(ibooker);
0163 
0164   for (unsigned int arm = 0; arm <= CTPPSDetId::maxArm; ++arm) {
0165     for (unsigned int pl = 0; pl <= TotemT2DetId::maxPlane; ++pl) {
0166       const TotemT2DetId detid(arm, pl, 0);
0167       const TotemT2DetId planeId(detid.planeId());
0168       planePlots_[planeId] = PlanePlots(ibooker, planeId, nbinsx_, nbinsy_);
0169     }
0170     const TotemT2DetId detid(arm, 0, 0);
0171     const TotemT2DetId secId(detid.armId());
0172     sectorPlots_[secId] = SectorPlots(ibooker, secId, nbinsx_, nbinsy_, windowsNum_);
0173   }
0174 
0175   // build a segmentation helper for the size of histograms previously booked
0176   segm_ = std::make_unique<TotemT2Segmentation>(iSetup.getData(geometryToken_), nbinsx_, nbinsy_);
0177 }
0178 
0179 void TotemT2DQMSource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0180   // fill digis information
0181   for (const auto& ds_digis : iEvent.get(digiToken_)) {
0182     const TotemT2DetId detid(ds_digis.detId());
0183     const TotemT2DetId planeId(detid.planeId());
0184     for (const auto& digi : ds_digis) {
0185       segm_->fill(planePlots_[planeId].digisMultiplicity->getTH2D(), detid);
0186       fillTriggerBitset(detid);
0187       fillErrorFlagsHistogram(digi);
0188       fillEdges(digi, detid);
0189     }
0190   }
0191 
0192   // fill rechits information
0193   std::unordered_map<unsigned int, std::set<unsigned int>> planes;
0194   for (const auto& ds_rechits : iEvent.get(rechitToken_)) {
0195     const TotemT2DetId detid(ds_rechits.detId());
0196     const TotemT2DetId planeId(detid.planeId());
0197     for (const auto& rechit : ds_rechits) {
0198       segm_->fill(planePlots_[planeId].rechitMultiplicity->getTH2D(), detid);
0199       fillToT(rechit, detid);
0200       fillActivePlanes(planes, detid);
0201     }
0202   }
0203 
0204   for (const auto& plt : sectorPlots_)
0205     plt.second.activePlanesCount->Fill(planes[plt.first].size());
0206 
0207   for (unsigned short arm = 0; arm <= CTPPSDetId::maxArm; ++arm)
0208     for (unsigned short plane = 0; plane <= 1; ++plane)
0209       for (unsigned short id = 0; id <= TotemT2DetId::maxChannel; ++id) {
0210         const TotemT2DetId detid(arm, plane, id);
0211         if (areChannelsTriggered(detid)) {
0212           const TotemT2DetId secId(detid.armId());
0213           segm_->fill(sectorPlots_[secId].triggerEmulator->getTH2D(), detid);
0214         }
0215       }
0216 
0217   clearTriggerBitset();
0218 }
0219 
0220 void TotemT2DQMSource::fillActivePlanes(std::unordered_map<unsigned int, std::set<unsigned int>>& planes,
0221                                         const TotemT2DetId& detid) {
0222   const TotemT2DetId secId(detid.armId());
0223   unsigned short pl = detid.plane();
0224 
0225   planes[secId].insert(pl);
0226 
0227   sectorPlots_[secId].activePlanes->Fill(pl);
0228 }
0229 
0230 void TotemT2DQMSource::fillTriggerBitset(const TotemT2DetId& detid) {
0231   const TotemT2DetId secId(detid.armId());
0232   unsigned short pl = detid.plane();
0233   unsigned short ch = detid.channel();
0234   sectorPlots_[secId].hitTilesArray[4 * pl + ch] = true;
0235 }
0236 
0237 void TotemT2DQMSource::clearTriggerBitset() {
0238   for (auto& sectorPlot : sectorPlots_)
0239     sectorPlot.second.hitTilesArray.reset();
0240 }
0241 
0242 bool TotemT2DQMSource::areChannelsTriggered(const TotemT2DetId& detid) {
0243   unsigned int channel = detid.channel();
0244 
0245   // prepare mask
0246   std::bitset<(TotemT2DetId::maxPlane + 1) * (TotemT2DetId::maxChannel + 1)> mask;
0247   // check if plane is even or not
0248   unsigned int pl = detid.plane() % 2 == 0 ? 0 : 1;
0249   // set only even or only odd plane bits for this channel
0250   for (; pl <= TotemT2DetId::maxPlane; pl += 2)
0251     mask[4 * pl + channel] = true;
0252   const TotemT2DetId secId(detid.armId());
0253   // check how many masked channels were hit
0254   unsigned int triggeredChannelsNumber = (mask & sectorPlots_[secId].hitTilesArray).count();
0255 
0256   return triggeredChannelsNumber >= SectorPlots::MINIMAL_TRIGGER;
0257 }
0258 
0259 void TotemT2DQMSource::bookErrorFlagsHistogram(DQMStore::IBooker& ibooker) {
0260   HPTDCErrorFlags_2D_ = ibooker.book2D("HPTDC Errors", " HPTDC Errors?", 8, -0.5, 7.5, 2, -0.5, 1.5);
0261   for (unsigned short error_index = 1; error_index <= 8; ++error_index)
0262     HPTDCErrorFlags_2D_->setBinLabel(error_index, "Flag " + std::to_string(error_index));
0263 
0264   int tmpIndex = 0;
0265   HPTDCErrorFlags_2D_->setBinLabel(++tmpIndex, "some id 0", /* axis */ 2);
0266   HPTDCErrorFlags_2D_->setBinLabel(++tmpIndex, "some id 1", /* axis */ 2);
0267 }
0268 
0269 void TotemT2DQMSource::fillErrorFlagsHistogram(const TotemT2Digi& digi) {
0270   // placeholder for error hitogram filling
0271   (void)digi;
0272 }
0273 
0274 void TotemT2DQMSource::fillEdges(const TotemT2Digi& digi, const TotemT2DetId& detid) {
0275   const TotemT2DetId secId(detid.armId());
0276   sectorPlots_[secId].leadingEdge->Fill(HPTDC_BIN_WIDTH_NS_ * digi.leadingEdge());
0277   sectorPlots_[secId].trailingEdge->Fill(HPTDC_BIN_WIDTH_NS_ * digi.trailingEdge());
0278 }
0279 
0280 void TotemT2DQMSource::fillToT(const TotemT2RecHit& rechit, const TotemT2DetId& detid) {
0281   const TotemT2DetId secId(detid.armId());
0282   sectorPlots_[secId].timeOverTreshold->Fill(rechit.toT());
0283 }
0284 
0285 DEFINE_FWK_MODULE(TotemT2DQMSource);