Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:01

0001 /****************************************************************************
0002  *
0003  * This is a part of TOTEM offline software.
0004  * Authors:
0005  *   Laurent Forthomme
0006  *   Arkadiusz Cwikla
0007  *   Fredrik Oljemark
0008  *
0009  ****************************************************************************/
0010 
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 
0018 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0019 #include "DQMServices/Core/interface/DQMStore.h"
0020 
0021 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0022 
0023 #include "DataFormats/CTPPSDetId/interface/TotemT2DetId.h"
0024 #include "DataFormats/TotemReco/interface/TotemT2Digi.h"
0025 
0026 #include "Geometry/Records/interface/TotemGeometryRcd.h"
0027 #include "DQM/CTPPS/interface/TotemT2Segmentation.h"
0028 
0029 #include <string>
0030 #include <bitset>
0031 
0032 class TotemT2DQMSource : public DQMEDAnalyzer {
0033 public:
0034   TotemT2DQMSource(const edm::ParameterSet&);
0035   ~TotemT2DQMSource() override = default;
0036 
0037 protected:
0038   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0039   void bookHistograms(DQMStore::IBooker&, const edm::Run&, const edm::EventSetup&) override;
0040   void analyze(const edm::Event&, const edm::EventSetup&) override;
0041 
0042 private:
0043   void fillActivePlanes(std::unordered_map<unsigned int, std::set<unsigned int>>&, const TotemT2DetId&);
0044   void fillTriggerBitset(const TotemT2DetId&);
0045   void clearTriggerBitset();
0046   bool areChannelsTriggered(const TotemT2DetId&);
0047   void bookErrorFlagsHistogram(DQMStore::IBooker&);
0048   void fillErrorFlagsHistogram(const TotemT2Digi&, const TotemT2DetId&);
0049   void fillEdges(const TotemT2Digi&, const TotemT2DetId&);
0050   void fillBX(const TotemT2Digi&, const TotemT2DetId&, const int);
0051   void fillToT(const TotemT2Digi&, const TotemT2DetId&);
0052   void fillFlags(const TotemT2Digi&, const TotemT2DetId&);
0053 
0054   const edm::EDGetTokenT<edmNew::DetSetVector<TotemT2Digi>> digiToken_;
0055 
0056   std::unique_ptr<TotemT2Segmentation> segm_;
0057 
0058   static constexpr double T2_BIN_WIDTH_NS_ = 25. / 4;
0059   MonitorElement* totemT2ErrorFlags_2D_ = nullptr;
0060 
0061   static constexpr int t2TE = 0, t2LE = 1, t2MT = 2, t2ML = 3;
0062 
0063   const unsigned int nbinsx_, nbinsy_;
0064   const unsigned int windowsNum_;
0065 
0066   struct SectorPlots {
0067     MonitorElement* activePlanes = nullptr;
0068     MonitorElement* activePlanesCount = nullptr;
0069     MonitorElement* activityPerBX = nullptr;
0070 
0071     MonitorElement* triggerEmulator = nullptr;
0072     std::bitset<(TotemT2DetId::maxPlane + 1) * (TotemT2DetId::maxChannel + 1)> hitTilesArray;
0073     static const unsigned int MINIMAL_TRIGGER = 3;
0074 
0075     MonitorElement *leadingEdge = nullptr, *trailingEdge = nullptr, *timeOverTreshold = nullptr, *eventFlags = nullptr;
0076 
0077     SectorPlots() = default;
0078     SectorPlots(
0079         DQMStore::IBooker& ibooker, unsigned int id, unsigned int nbinsx, unsigned int nbinsy, unsigned int windowsNum);
0080   };
0081 
0082   struct PlanePlots {
0083     MonitorElement* digisMultiplicity = nullptr;
0084     MonitorElement* rechitMultiplicity = nullptr;
0085     MonitorElement* eventFlagsPl = nullptr;
0086 
0087     PlanePlots() = default;
0088     PlanePlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int nbinsx, unsigned int nbinsy);
0089   };
0090   struct ChannelPlots {
0091     MonitorElement* leadingEdgeCh = nullptr;
0092     MonitorElement* trailingEdgeCh = nullptr;
0093     MonitorElement* timeOverTresholdCh = nullptr;
0094     MonitorElement* eventFlagsCh = nullptr;
0095     MonitorElement* activityPerBXCh = nullptr;
0096 
0097     ChannelPlots() = default;
0098     ChannelPlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int windowsNum);
0099   };
0100 
0101   std::unordered_map<unsigned int, SectorPlots> sectorPlots_;
0102   std::unordered_map<unsigned int, PlanePlots> planePlots_;
0103   std::unordered_map<unsigned int, ChannelPlots> channelPlots_;
0104 };
0105 
0106 TotemT2DQMSource::SectorPlots::SectorPlots(
0107     DQMStore::IBooker& ibooker, unsigned int id, unsigned int nbinsx, unsigned int nbinsy, unsigned int windowsNum) {
0108   std::string title, path;
0109 
0110   TotemT2DetId(id).armName(path, TotemT2DetId::nPath);
0111   ibooker.setCurrentFolder(path);
0112 
0113   TotemT2DetId(id).armName(title, TotemT2DetId::nFull);
0114 
0115   activePlanes =
0116       ibooker.book1D("active planes", title + " which planes are active (LE+TE digis);plane number", 8, -0.5, 7.5);
0117 
0118   activePlanesCount = ibooker.book1D("number of active planes",
0119                                      title + " how many planes are active (LE+TE digis);number of active planes",
0120                                      9,
0121                                      -0.5,
0122                                      8.5);
0123 
0124   activityPerBX =
0125       ibooker.book1D("activity per BX CMS", title + " Activity (=LE) per BX;Event.BX", 3600, -1.5, 3598. + 0.5);
0126 
0127   triggerEmulator = ibooker.book2DD("trigger emulator",
0128                                     title + " trigger emulator (LE+TE digis);arbitrary unit;arbitrary unit",
0129                                     nbinsx,
0130                                     -0.5,
0131                                     double(nbinsx) - 0.5,
0132                                     nbinsy,
0133                                     -0.5,
0134                                     double(nbinsy) - 0.5);
0135   leadingEdge = ibooker.book1D(
0136       "leading edge", title + " leading edge (all DIGIs); leading edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0137   trailingEdge = ibooker.book1D(
0138       "trailing edge", title + " trailing edge (all DIGIs); trailing edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0139 
0140   timeOverTreshold = ibooker.book1D("time over threshold",
0141                                     title + " time over threshold (digi, =0 if LE or TE=0);time over threshold (ns)",
0142                                     500,
0143                                     -50,
0144                                     200);
0145 
0146   eventFlags = ibooker.book1D(
0147       "event flags", title + " event flags (digi);Event flags (TE/LE valid, TE/LE multiple)", 4, -0.5, 3.5);
0148 
0149   for (unsigned short flag_index = 1; flag_index <= 4; ++flag_index)
0150     eventFlags->setBinLabel(flag_index, "Flag " + std::to_string(flag_index));
0151 }
0152 
0153 TotemT2DQMSource::PlanePlots::PlanePlots(DQMStore::IBooker& ibooker,
0154                                          unsigned int id,
0155                                          unsigned int nbinsx,
0156                                          unsigned int nbinsy) {
0157   std::string title, path;
0158   TotemT2DetId(id).planeName(title, TotemT2DetId::nFull);
0159   TotemT2DetId(id).planeName(path, TotemT2DetId::nPath);
0160   ibooker.setCurrentFolder(path);
0161 
0162   digisMultiplicity = ibooker.book2DD("digis multiplicity",
0163                                       title + " digis multiplicity (=LE);arbitrary unit;arbitrary unit",
0164                                       nbinsx,
0165                                       -0.5,
0166                                       double(nbinsx) - 0.5,
0167                                       nbinsy,
0168                                       -0.5,
0169                                       double(nbinsy) - 0.5);
0170   rechitMultiplicity = ibooker.book2DD("rechits multiplicity",
0171                                        title + " rechits multiplicity;x;y",
0172                                        nbinsx,
0173                                        -0.5,
0174                                        double(nbinsx) - 0.5,
0175                                        nbinsy,
0176                                        -0.5,
0177                                        double(nbinsy) - 0.5);
0178 
0179   eventFlagsPl = ibooker.book1D(
0180       "event flags", title + " event flags (digi);Event flags (TE/LE valid, TE/LE multiple)", 4, -0.5, 3.5);
0181 
0182   for (unsigned short flag_index = 1; flag_index <= 4; ++flag_index)
0183     eventFlagsPl->setBinLabel(flag_index, "Flag " + std::to_string(flag_index));
0184 }
0185 
0186 TotemT2DQMSource::ChannelPlots::ChannelPlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int windowsNum) {
0187   std::string title, path;
0188   TotemT2DetId(id).channelName(title, TotemT2DetId::nFull);
0189   TotemT2DetId(id).channelName(path, TotemT2DetId::nPath);
0190   ibooker.setCurrentFolder(path);
0191 
0192   leadingEdgeCh = ibooker.book1D(
0193       "leading edge", title + " leading edge (all DIGIs); leading edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0194   trailingEdgeCh = ibooker.book1D(
0195       "trailing edge", title + " trailing edge (all DIGIs); trailing edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0196 
0197   timeOverTresholdCh = ibooker.book1D("time over threshold",
0198                                       title + " time over threshold (digi, =0 if LE or TE=0);time over threshold (ns)",
0199                                       500,
0200                                       -50,
0201                                       200);
0202 
0203   eventFlagsCh = ibooker.book1D(
0204       "event flags", title + " event flags (digi);Event flags (TE/LE valid, TE/LE multiple)", 4, -0.5, 3.5);
0205 
0206   activityPerBXCh =
0207       ibooker.book1D("activity per BX", title + " Activity (=LE) per BX;Event.BX", 1000, -1.5, 998. + 0.5);
0208 
0209   for (unsigned short flag_index = 1; flag_index <= 4; ++flag_index)
0210     eventFlagsCh->setBinLabel(flag_index, "Flag " + std::to_string(flag_index));
0211 }
0212 
0213 TotemT2DQMSource::TotemT2DQMSource(const edm::ParameterSet& iConfig)
0214     : digiToken_(consumes<edmNew::DetSetVector<TotemT2Digi>>(iConfig.getParameter<edm::InputTag>("digisTag"))),
0215       nbinsx_(iConfig.getParameter<unsigned int>("nbinsx")),
0216       nbinsy_(iConfig.getParameter<unsigned int>("nbinsy")),
0217       windowsNum_(iConfig.getParameter<unsigned int>("windowsNum")) {}
0218 
0219 void TotemT2DQMSource::dqmBeginRun(const edm::Run&, const edm::EventSetup&) {}
0220 
0221 void TotemT2DQMSource::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup& iSetup) {
0222   ibooker.cd();
0223   ibooker.setCurrentFolder("TotemT2");
0224 
0225   bookErrorFlagsHistogram(ibooker);
0226 
0227   for (unsigned int arm = 0; arm <= CTPPSDetId::maxArm; ++arm) {
0228     for (unsigned int pl = 0; pl <= TotemT2DetId::maxPlane; ++pl) {
0229       const TotemT2DetId detid(arm, pl, 0);
0230       const TotemT2DetId planeId(detid.planeId());
0231       planePlots_[planeId] = PlanePlots(ibooker, planeId, nbinsx_, nbinsy_);
0232       for (unsigned int ch = 0; ch <= TotemT2DetId::maxChannel; ++ch) {
0233         const TotemT2DetId detidCh(arm, pl, ch);
0234         channelPlots_[detidCh] = ChannelPlots(ibooker, detidCh, windowsNum_);
0235       }
0236     }
0237     const TotemT2DetId detid(arm, 0, 0);
0238     const TotemT2DetId secId(detid.armId());
0239     sectorPlots_[secId] = SectorPlots(ibooker, secId, nbinsx_, nbinsy_, windowsNum_);
0240   }
0241 
0242   // build a segmentation helper for the size of histograms previously booked
0243   segm_ = std::make_unique<TotemT2Segmentation>(nbinsx_, nbinsy_);
0244 }
0245 
0246 void TotemT2DQMSource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0247   // fill digis information
0248   std::unordered_map<unsigned int, std::set<unsigned int>> planes;
0249 
0250   for (const auto& ds_digis : iEvent.get(digiToken_)) {
0251     if (!ds_digis.empty()) {
0252       const TotemT2DetId detid(ds_digis.detId());
0253       const TotemT2DetId planeId(detid.planeId());
0254       for (const auto& digi : ds_digis) {
0255         if (digi.hasLE()) {
0256           segm_->fill(planePlots_[planeId].digisMultiplicity->getTH2D(), detid);
0257           if (digi.hasTE())
0258             fillTriggerBitset(detid);
0259         }
0260         fillErrorFlagsHistogram(digi, detid);
0261         fillEdges(digi, detid);
0262         fillToT(digi, detid);
0263         fillFlags(digi, detid);
0264         int bx = iEvent.bunchCrossing();
0265         fillBX(digi, detid, bx);
0266         if (digi.hasLE() && digi.hasTE())  //good digi
0267           fillActivePlanes(planes, detid);
0268       }
0269     }
0270   }
0271 
0272   for (const auto& plt : sectorPlots_)
0273     plt.second.activePlanesCount->Fill(planes[plt.first].size());
0274 
0275   for (unsigned short arm = 0; arm <= CTPPSDetId::maxArm; ++arm)
0276     for (unsigned short plane = 0; plane <= 1; ++plane)
0277       for (unsigned short id = 0; id <= TotemT2DetId::maxChannel; ++id) {
0278         const TotemT2DetId detid(arm, plane, id);
0279         if (areChannelsTriggered(detid)) {
0280           const TotemT2DetId secId(detid.armId());
0281           segm_->fill(sectorPlots_[secId].triggerEmulator->getTH2D(), detid);
0282         }
0283       }
0284 
0285   // fill rechits information
0286 
0287   clearTriggerBitset();
0288 }
0289 
0290 void TotemT2DQMSource::fillActivePlanes(std::unordered_map<unsigned int, std::set<unsigned int>>& planes,
0291                                         const TotemT2DetId& detid) {
0292   const TotemT2DetId secId(detid.armId());
0293   unsigned short pl = detid.plane();
0294 
0295   planes[secId].insert(pl);
0296 
0297   sectorPlots_[secId].activePlanes->Fill(pl);
0298 }
0299 
0300 void TotemT2DQMSource::fillTriggerBitset(const TotemT2DetId& detid) {
0301   const TotemT2DetId secId(detid.armId());
0302   unsigned short pl = detid.plane();
0303   unsigned short ch = detid.channel();
0304   sectorPlots_[secId].hitTilesArray[4 * pl + ch] = true;
0305 }
0306 
0307 void TotemT2DQMSource::clearTriggerBitset() {
0308   for (auto& sectorPlot : sectorPlots_)
0309     sectorPlot.second.hitTilesArray.reset();
0310 }
0311 
0312 bool TotemT2DQMSource::areChannelsTriggered(const TotemT2DetId& detid) {
0313   unsigned int channel = detid.channel();
0314 
0315   // prepare mask
0316   std::bitset<(TotemT2DetId::maxPlane + 1) * (TotemT2DetId::maxChannel + 1)> mask;
0317   // check if plane is even or not
0318   unsigned int pl = detid.plane() % 2 == 0 ? 0 : 1;
0319   // set only even or only odd plane bits for this channel
0320   for (; pl <= TotemT2DetId::maxPlane; pl += 2)
0321     mask[4 * pl + channel] = true;
0322   const TotemT2DetId secId(detid.armId());
0323   // check how many masked channels were hit
0324   unsigned int triggeredChannelsNumber = (mask & sectorPlots_[secId].hitTilesArray).count();
0325 
0326   return triggeredChannelsNumber >= SectorPlots::MINIMAL_TRIGGER;
0327 }
0328 
0329 void TotemT2DQMSource::bookErrorFlagsHistogram(DQMStore::IBooker& ibooker) {
0330   totemT2ErrorFlags_2D_ = ibooker.book2D("nt2 readout flags", " nt2 readout flags", 4, -0.5, 3.5, 2, -0.5, 1.5);
0331   for (unsigned short error_index = 1; error_index <= 4; ++error_index)
0332     totemT2ErrorFlags_2D_->setBinLabel(error_index, "Flag " + std::to_string(error_index));
0333 
0334   int tmpIndex = 0;
0335   totemT2ErrorFlags_2D_->setBinLabel(++tmpIndex, "arm 4-5", /* axis */ 2);
0336   totemT2ErrorFlags_2D_->setBinLabel(++tmpIndex, "arm 5-6", /* axis */ 2);
0337 }
0338 
0339 void TotemT2DQMSource::fillErrorFlagsHistogram(const TotemT2Digi& digi, const TotemT2DetId& detid) {
0340   // readout flags histogram filling
0341   for (unsigned int i = 0; i < 4; i++) {
0342     if (digi.status() & (1 << i))
0343       totemT2ErrorFlags_2D_->Fill(i + 0.0, detid.arm() + 0.0);
0344   }
0345 }
0346 
0347 void TotemT2DQMSource::fillEdges(const TotemT2Digi& digi, const TotemT2DetId& detid) {
0348   const TotemT2DetId secId(detid.armId());
0349   sectorPlots_[secId].leadingEdge->Fill(T2_BIN_WIDTH_NS_ * digi.leadingEdge());
0350   sectorPlots_[secId].trailingEdge->Fill(T2_BIN_WIDTH_NS_ * digi.trailingEdge());
0351   channelPlots_[detid].leadingEdgeCh->Fill(T2_BIN_WIDTH_NS_ * digi.leadingEdge());
0352   channelPlots_[detid].trailingEdgeCh->Fill(T2_BIN_WIDTH_NS_ * digi.trailingEdge());
0353 }
0354 
0355 void TotemT2DQMSource::fillBX(const TotemT2Digi& digi, const TotemT2DetId& detid, const int bx) {
0356   const TotemT2DetId secId(detid.armId());
0357   if (digi.hasLE()) {
0358     sectorPlots_[secId].activityPerBX->Fill(bx);
0359     channelPlots_[detid].activityPerBXCh->Fill(bx);
0360   }
0361 }
0362 
0363 void TotemT2DQMSource::fillToT(const TotemT2Digi& digi, const TotemT2DetId& detid) {
0364   const TotemT2DetId secId(detid.armId());
0365 
0366   const int t_lead = digi.leadingEdge(), t_trail = digi.trailingEdge();
0367   // don't skip no-edge digis
0368   double toT = 0.;
0369   if (digi.hasLE() && digi.hasTE()) {
0370     toT = (t_trail - t_lead) * T2_BIN_WIDTH_NS_;  // in ns
0371   }
0372 
0373   sectorPlots_[secId].timeOverTreshold->Fill(toT);
0374   channelPlots_[detid].timeOverTresholdCh->Fill(toT);
0375 }
0376 
0377 void TotemT2DQMSource::fillFlags(const TotemT2Digi& digi, const TotemT2DetId& detid) {
0378   const TotemT2DetId secId(detid.armId());
0379   const TotemT2DetId planeId(detid.planeId());
0380   if (digi.hasTE()) {
0381     sectorPlots_[secId].eventFlags->Fill(t2TE + 0.0);
0382     planePlots_[planeId].eventFlagsPl->Fill(t2TE + 0.0);
0383     channelPlots_[detid].eventFlagsCh->Fill(t2TE + 0.0);
0384   }
0385 
0386   if (digi.hasLE()) {
0387     sectorPlots_[secId].eventFlags->Fill(t2LE + 0.0);
0388     planePlots_[planeId].eventFlagsPl->Fill(t2LE + 0.0);
0389     channelPlots_[detid].eventFlagsCh->Fill(t2LE + 0.0);
0390   }
0391 
0392   if (digi.hasManyTE()) {
0393     sectorPlots_[secId].eventFlags->Fill(t2MT + 0.0);
0394     planePlots_[planeId].eventFlagsPl->Fill(t2MT + 0.0);
0395     channelPlots_[detid].eventFlagsCh->Fill(t2MT + 0.0);
0396   }
0397 
0398   if (digi.hasManyLE()) {
0399     sectorPlots_[secId].eventFlags->Fill(t2ML + 0.0);
0400     planePlots_[planeId].eventFlagsPl->Fill(t2ML + 0.0);
0401     channelPlots_[detid].eventFlagsCh->Fill(t2ML + 0.0);
0402   }
0403 }
0404 
0405 DEFINE_FWK_MODULE(TotemT2DQMSource);