Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-30 02:16:25

0001 /****************************************************************************
0002 *
0003 * This is a part of CTPPSDQM software.
0004 * Authors:
0005 *   Jan Kašpar (jan.kaspar@gmail.com)
0006 *   Nicola Minafra
0007 *   Laurent Forthomme
0008 *   Christopher Misan
0009 *
0010 ****************************************************************************/
0011 
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/Utilities/interface/ESGetToken.h"
0017 #include "FWCore/Utilities/interface/InputTag.h"
0018 #include "FWCore/Utilities/interface/Transition.h"
0019 #include "FWCore/Framework/interface/Run.h"
0020 
0021 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0022 #include "DQMServices/Core/interface/DQMStore.h"
0023 
0024 #include "DataFormats/Provenance/interface/EventRange.h"
0025 #include "DataFormats/CTPPSDigi/interface/TotemVFATStatus.h"
0026 #include "DataFormats/CTPPSDigi/interface/TotemFEDInfo.h"
0027 #include "DataFormats/Common/interface/DetSetVector.h"
0028 
0029 #include "DataFormats/CTPPSDetId/interface/CTPPSDiamondDetId.h"
0030 #include "DataFormats/CTPPSDigi/interface/CTPPSDiamondDigi.h"
0031 
0032 #include "DataFormats/CTPPSReco/interface/CTPPSPixelLocalTrack.h"
0033 #include "DataFormats/CTPPSDetId/interface/CTPPSPixelDetId.h"
0034 
0035 #include "DataFormats/CTPPSReco/interface/CTPPSDiamondRecHit.h"
0036 #include "DataFormats/CTPPSReco/interface/CTPPSDiamondLocalTrack.h"
0037 #include "Geometry/VeryForwardGeometryBuilder/interface/CTPPSGeometry.h"
0038 #include "Geometry/Records/interface/VeryForwardRealGeometryRecord.h"
0039 
0040 // #include "CondFormats/RunInfo/interface/LHCInfo.h"
0041 // #include "CondFormats/DataRecord/interface/LHCInfoRcd.h"
0042 
0043 #include <string>
0044 
0045 //----------------------------------------------------------------------------------------------------
0046 
0047 // Utility for efficiency computations
0048 bool channelAlignedWithTrack(const CTPPSGeometry* geom,
0049                              const CTPPSDiamondDetId& detid,
0050                              const CTPPSDiamondLocalTrack& localTrack,
0051                              const float tolerance = 1) {
0052   const DetGeomDesc* det = geom->sensor(detid);
0053   const float x_pos = det->translation().x(),
0054               x_width = 2.0 * det->getDiamondDimensions().xHalfWidth;  // parameters stand for half the size
0055   return ((x_pos + 0.5 * x_width > localTrack.x0() - localTrack.x0Sigma() - tolerance &&
0056            x_pos + 0.5 * x_width < localTrack.x0() + localTrack.x0Sigma() + tolerance) ||
0057           (x_pos - 0.5 * x_width > localTrack.x0() - localTrack.x0Sigma() - tolerance &&
0058            x_pos - 0.5 * x_width < localTrack.x0() + localTrack.x0Sigma() + tolerance) ||
0059           (x_pos - 0.5 * x_width < localTrack.x0() - localTrack.x0Sigma() - tolerance &&
0060            x_pos + 0.5 * x_width > localTrack.x0() + localTrack.x0Sigma() + tolerance));
0061 }
0062 
0063 namespace dds {
0064   struct Cache {
0065     std::unordered_map<unsigned int, std::unique_ptr<TH2F>> hitDistribution2dMap;
0066 
0067     std::unordered_map<unsigned int, unsigned long> hitsCounterMap;
0068   };
0069 }  // namespace dds
0070 
0071 class CTPPSDiamondDQMSource : public DQMOneEDAnalyzer<edm::LuminosityBlockCache<dds::Cache>> {
0072 public:
0073   CTPPSDiamondDQMSource(const edm::ParameterSet&);
0074 
0075 protected:
0076   void dqmBeginRun(const edm::Run&, const edm::EventSetup&) override;
0077   void dqmEndRun(edm::Run const&, edm::EventSetup const&) override;
0078   void bookHistograms(DQMStore::IBooker&, const edm::Run&, const edm::EventSetup&) override;
0079   void analyze(const edm::Event&, const edm::EventSetup&) override;
0080   std::shared_ptr<dds::Cache> globalBeginLuminosityBlock(const edm::LuminosityBlock&,
0081                                                          const edm::EventSetup&) const override;
0082   void globalEndLuminosityBlock(const edm::LuminosityBlock&, const edm::EventSetup&) override;
0083 
0084 private:
0085   // Constants
0086   /// Number of seconds per lumisection: used to compute hit rates in Hz
0087   static constexpr double SEC_PER_LUMI_SECTION = 23.31;
0088   /// Channel ID of the VFAT that contains clock data
0089   static constexpr unsigned short CHANNEL_OF_VFAT_CLOCK = 30;
0090   /// Bin width of histograms showing hits and tracks (in mm)
0091   static constexpr double DISPLAY_RESOLUTION_FOR_HITS_MM = 0.1;
0092   static constexpr double INV_DISPLAY_RESOLUTION_FOR_HITS_MM = 1. / DISPLAY_RESOLUTION_FOR_HITS_MM;
0093   /// ns per HPTDC bin
0094   static constexpr double HPTDC_BIN_WIDTH_NS = 25. / 1024;
0095   static constexpr unsigned short CTPPS_PIXEL_STATION_ID = 2;
0096   static constexpr unsigned short CTPPS_FAR_RP_ID = 3;
0097   static constexpr unsigned short CTPPS_DIAMOND_NUM_OF_CHANNELS = 12;
0098   static constexpr unsigned short CTPPS_FED_ID_45 = 583;
0099   static constexpr unsigned short CTPPS_FED_ID_56 = 582;
0100   static constexpr unsigned short HPTDC_0_CHANNEL = 6;
0101   static constexpr unsigned short HPTDC_1_CHANNEL = 7;
0102   /// Number of OOT indices monitored
0103   static constexpr unsigned int FIRST_RUN_W_PIXELS = 300000;
0104 
0105   bool perLSsaving_;  //to avoid nanoDQMIO crashing, driven by  DQMServices/Core/python/DQMStore_cfi.py
0106 
0107   /// plots related to the whole system
0108   struct GlobalPlots {
0109     GlobalPlots() = default;
0110     GlobalPlots(DQMStore::IBooker& ibooker);
0111   };
0112   /// plots related to one sector
0113   struct SectorPlots {
0114     // Tracks
0115     MonitorElement* trackCorrelation = nullptr;
0116     MonitorElement* trackCorrelationLowMultiplicity = nullptr;
0117     MonitorElement* digiSentPercentage = nullptr;
0118     SectorPlots(){};
0119     SectorPlots(DQMStore::IBooker& ibooker, unsigned int id, bool plotOnline);
0120   };
0121   /// plots related to one Diamond detector package
0122   struct PotPlots {
0123     std::unordered_map<unsigned int, MonitorElement*> activity_per_bx;
0124 
0125     MonitorElement* hitDistribution2d = nullptr;
0126     MonitorElement* hitDistribution2d_lumisection = nullptr;
0127     MonitorElement* hitDistribution2dOOT = nullptr;
0128     MonitorElement* hitDistribution2dOOT_le = nullptr;
0129     MonitorElement *activePlanes = nullptr, *activePlanesInclusive = nullptr;
0130 
0131     MonitorElement* trackDistribution = nullptr;
0132     MonitorElement* trackDistributionOOT = nullptr;
0133 
0134     std::unordered_map<unsigned int, MonitorElement*> pixelTomographyAll;
0135 
0136     MonitorElement *leadingEdgeCumulative_both = nullptr, *leadingEdgeCumulative_all = nullptr,
0137                    *leadingEdgeCumulative_le = nullptr, *trailingEdgeCumulative_te = nullptr;
0138     MonitorElement* timeOverThresholdCumulativePot = nullptr;  //, *leadingTrailingCorrelationPot = nullptr;
0139     MonitorElement* leadingWithoutTrailingCumulativePot = nullptr;
0140 
0141     MonitorElement* ECCheck = nullptr;
0142 
0143     MonitorElement* HPTDCErrorFlags_2D = nullptr;
0144     MonitorElement* MHComprensive = nullptr;
0145 
0146     MonitorElement* recHitTime = nullptr;
0147 
0148     // MonitorElement* clock_Digi1_le = nullptr;
0149     // MonitorElement* clock_Digi1_te = nullptr;
0150     // MonitorElement* clock_Digi3_le = nullptr;
0151     // MonitorElement* clock_Digi3_te = nullptr;
0152 
0153     unsigned int HitCounter, MHCounter, LeadingOnlyCounter, TrailingOnlyCounter, CompleteCounter;
0154 
0155     std::map<int, int> effTriplecountingChMap;
0156     std::map<int, int> effDoublecountingChMap;
0157     MonitorElement* EfficiencyOfChannelsInPot = nullptr;
0158     TH2F pixelTracksMap;
0159 
0160     // MonitorElement* TOTVsLS = nullptr;
0161     // MonitorElement* trackTimeVsLS = nullptr;
0162     MonitorElement* trackTimeVsBX = nullptr;
0163     // MonitorElement* trackTimeVsXAngle = nullptr;
0164 
0165     // MonitorElement* TOTVsLSProfile = nullptr;
0166     // MonitorElement* trackTimeVsLSProfile = nullptr;
0167     MonitorElement* trackTimeVsBXProfile = nullptr;
0168     // MonitorElement* trackTimeVsXAngleProfile = nullptr;
0169 
0170     PotPlots() = default;
0171     PotPlots(DQMStore::IBooker& ibooker,
0172              unsigned int id,
0173              unsigned int windowsNum,
0174              bool plotOnline,
0175              bool plotOffline,
0176              bool perLSsaving);
0177   };
0178   /// plots related to one Diamond plane
0179   struct PlanePlots {
0180     MonitorElement* digiProfileCumulativePerPlane = nullptr;
0181     MonitorElement* hitProfile = nullptr;
0182     MonitorElement* hit_multiplicity = nullptr;
0183 
0184     MonitorElement* pixelTomography_far = nullptr;
0185     MonitorElement* EfficiencyWRTPixelsInPlane = nullptr;
0186 
0187     TH2F pixelTracksMapWithDiamonds;
0188 
0189     PlanePlots() = default;
0190     PlanePlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int windowsNum);
0191   };
0192   /// plots related to one Diamond channel
0193   struct ChannelPlots {
0194     std::unordered_map<unsigned int, MonitorElement*> activity_per_bx;
0195 
0196     MonitorElement* HPTDCErrorFlags = nullptr;
0197     MonitorElement *leadingEdgeCumulative_both = nullptr, *leadingEdgeCumulative_le = nullptr,
0198                    *trailingEdgeCumulative_te = nullptr;
0199     MonitorElement* TimeOverThresholdCumulativePerChannel = nullptr;
0200     //MonitorElement* LeadingTrailingCorrelationPerChannel = nullptr;
0201     MonitorElement* leadingWithoutTrailing = nullptr;
0202     MonitorElement* pixelTomography_far = nullptr;
0203     MonitorElement* hit_rate = nullptr;
0204     MonitorElement* recHitTime = nullptr;
0205 
0206     unsigned int HitCounter, MHCounter, LeadingOnlyCounter, TrailingOnlyCounter, CompleteCounter;
0207 
0208     ChannelPlots() = default;
0209     ChannelPlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int windowsNum);
0210   };
0211 
0212   void checkEventNumber(const CTPPSDiamondDetId&, const TotemFEDInfo&, const TotemVFATStatus&, PotPlots&, int&) const;
0213 
0214   edm::EDGetTokenT<edm::DetSetVector<TotemVFATStatus>> tokenStatus_;
0215   edm::EDGetTokenT<edm::DetSetVector<CTPPSPixelLocalTrack>> tokenPixelTrack_;
0216   edm::EDGetTokenT<edm::DetSetVector<CTPPSDiamondDigi>> tokenDigi_;
0217   edm::EDGetTokenT<edm::DetSetVector<CTPPSDiamondRecHit>> tokenDiamondHit_;
0218   edm::EDGetTokenT<edm::DetSetVector<CTPPSDiamondLocalTrack>> tokenDiamondTrack_;
0219   edm::EDGetTokenT<std::vector<TotemFEDInfo>> tokenFEDInfo_;
0220 
0221   edm::ESGetToken<CTPPSGeometry, VeryForwardRealGeometryRecord> ctppsGeometryRunToken_;
0222   edm::ESGetToken<CTPPSGeometry, VeryForwardRealGeometryRecord> ctppsGeometryEventToken_;
0223   // edm::ESGetToken<LHCInfo, LHCInfoRcd> ctppsLhcInfoToken_;
0224 
0225   bool excludeMultipleHits_;
0226   const bool extract_digi_info_;
0227   struct DiamondShifts {
0228     double global, withPixels;
0229   };
0230   std::unordered_map<CTPPSDetId, DiamondShifts> diamShifts_;
0231   std::vector<std::pair<edm::EventRange, int>> runParameters_;
0232   int centralOOT_;
0233   unsigned int verbosity_;
0234   const bool plotOnline_;
0235   const bool plotOffline_;
0236   unsigned int windowsNum_;
0237   unsigned int trackCorrelationThreshold_;
0238 
0239   GlobalPlots globalPlot_;
0240   std::unordered_map<unsigned int, PotPlots> potPlots_;
0241   std::unordered_map<unsigned int, SectorPlots> sectorPlots_;
0242   std::unordered_map<unsigned int, PlanePlots> planePlots_;
0243   std::unordered_map<unsigned int, ChannelPlots> channelPlots_;
0244 
0245   int EC_difference_56_, EC_difference_45_;
0246 };
0247 
0248 //----------------------------------------------------------------------------------------------------
0249 
0250 CTPPSDiamondDQMSource::GlobalPlots::GlobalPlots(DQMStore::IBooker& ibooker) { ibooker.setCurrentFolder("CTPPS"); }
0251 
0252 //----------------------------------------------------------------------------------------------------
0253 
0254 CTPPSDiamondDQMSource::SectorPlots::SectorPlots(DQMStore::IBooker& ibooker, unsigned int id, bool plotOnline) {
0255   std::string path, title;
0256   CTPPSDiamondDetId(id).armName(path, CTPPSDiamondDetId::nPath);
0257   ibooker.setCurrentFolder(path);
0258 
0259   CTPPSDiamondDetId(id).armName(title, CTPPSDiamondDetId::nFull);
0260 
0261   trackCorrelation = ibooker.book2D("tracks correlation near-far",
0262                                     title + " tracks correlation near-far;x (mm);x (mm)",
0263                                     19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0264                                     -1,
0265                                     18,
0266                                     19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0267                                     -1,
0268                                     18);
0269   trackCorrelationLowMultiplicity =
0270       ibooker.book2D("tracks correlation with low multiplicity near-far",
0271                      title + " tracks correlation with low multiplicity near-far;x (mm);x (mm)",
0272                      19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0273                      -1,
0274                      18,
0275                      19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0276                      -1,
0277                      18);
0278 }
0279 
0280 //----------------------------------------------------------------------------------------------------
0281 CTPPSDiamondDQMSource::PotPlots::PotPlots(DQMStore::IBooker& ibooker,
0282                                           unsigned int id,
0283                                           unsigned int windowsNum,
0284                                           bool plotOnline,
0285                                           bool plotOffline,
0286                                           bool perLSsaving)
0287     : HitCounter(0),
0288       MHCounter(0),
0289       LeadingOnlyCounter(0),
0290       TrailingOnlyCounter(0),
0291       CompleteCounter(0),
0292       pixelTracksMap("Pixel track maps for efficiency", "Pixel track maps for efficiency", 25, 0, 25, 12, -2, 10) {
0293   std::string path, title;
0294   CTPPSDiamondDetId(id).rpName(path, CTPPSDiamondDetId::nPath);
0295   ibooker.setCurrentFolder(path);
0296 
0297   CTPPSDiamondDetId(id).rpName(title, CTPPSDiamondDetId::nFull);
0298 
0299   if (plotOnline) {
0300     hitDistribution2d_lumisection =
0301         ibooker.book2D("hits in planes lumisection",
0302                        title + " hits in planes in the last lumisection;plane number;x (mm)",
0303                        10,
0304                        -0.5,
0305                        4.5,
0306                        19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0307                        -0.5,
0308                        18.5);
0309 
0310     hitDistribution2dOOT_le =
0311         ibooker.book2D("hits with OOT in planes (le only)",
0312                        title + " hits with OOT in planes (le only);plane number, OOT index;x (mm)",
0313                        1 + windowsNum * 4,
0314                        -1. / windowsNum,
0315                        4,
0316                        19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0317                        -0.5,
0318                        18.5);
0319 
0320     activePlanesInclusive =
0321         ibooker.book1D("active planes inclusive",
0322                        title + " active planes, MH and le only included (per event);number of active planes",
0323                        6,
0324                        -0.5,
0325                        5.5);
0326 
0327     ECCheck = ibooker.book1D("optorxEC(8bit) - vfatEC", title + " EC Error;optorxEC-vfatEC", 50, -25, 25);
0328 
0329     EfficiencyOfChannelsInPot =
0330         ibooker.book2D("Efficiency in channels",
0331                        title + " Efficiency (%) in channels (diamonds only);plane number;ch number",
0332                        10,
0333                        -0.5,
0334                        4.5,
0335                        14,
0336                        -1,
0337                        13);
0338   }
0339 
0340   if (plotOffline && !perLSsaving) {
0341     ibooker.setCurrentFolder(path + "/timing_profiles");
0342     // TOTVsLS=ibooker.book2D("ToT vs LS",title +" ToT vs LS;LS;ToT(ns)",4000,0,4000, 200,5,25);
0343     // trackTimeVsLS=ibooker.book2D("track time vs LS",title+" track time vs LS;LS;track_time(ns)",4000,0,4000, 500, -25, 25);
0344     trackTimeVsBX =
0345         ibooker.book2D("track time vs BX", title + " track time vs BX;BX;track_time(ns)", 4000, 0, 4000, 500, -25, 25);
0346     // trackTimeVsXAngle = ibooker.book2D(
0347     //     "track time vs xangle", title + " track time vs xangle;xangle;track_time(ns)", 60, 120, 180, 500, -25, 25);
0348 
0349     // TOTVsLSProfile=ibooker.bookProfile("ToT vs LS profile",title+" ToT vs LS profile;LS;track_time(ns)", 500, -25, 25,4000,0,4000);
0350     // trackTimeVsLSProfile=ibooker.bookProfile("track time vs LS profile",title+" track time vs LS profile;LS;track_time(ns)", 500, -25, 25,4000,0,4000);
0351     trackTimeVsBXProfile = ibooker.bookProfile(
0352         "track time vs BX profile", title + " track time vs BX profile;BX;track_time(ns)", 500, -25, 25, 4000, 0, 4000);
0353     // trackTimeVsXAngleProfile = ibooker.bookProfile("track time vs xangle profile",
0354     //                                                title + " track time vs xangle profile;xangle;track_time(ns)",
0355     //                                                500,
0356     //                                                -25,
0357     //                                                25,
0358     //                                                60,
0359     //                                                120,
0360     //                                                180);
0361     ibooker.setCurrentFolder(path);
0362   }
0363 
0364   for (unsigned int i = 0; i < windowsNum; i++) {
0365     std::string window = std::to_string(i * 25) + "-" + std::to_string((i + 1) * 25);
0366     activity_per_bx[i] = ibooker.book1D(
0367         "activity per BX " + window, title + " Activity per BX " + window + " ns;Event.BX", 3600, -1.5, 3598. + 0.5);
0368     pixelTomographyAll[i] =
0369         ibooker.book2D("tomography pixel " + window,
0370                        title + " tomography with pixel " + window + " ns (all planes);x + 25*plane(mm);y (mm)",
0371                        100,
0372                        0,
0373                        100,
0374                        10,
0375                        -5,
0376                        5);
0377   }
0378 
0379   hitDistribution2d = ibooker.book2D("hits in planes",
0380                                      title + " hits in planes;plane number;x (mm)",
0381                                      10,
0382                                      -0.5,
0383                                      4.5,
0384                                      19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0385                                      -0.5,
0386                                      18.5);
0387 
0388   hitDistribution2dOOT = ibooker.book2D("hits with OOT in planes",
0389                                         title + " hits with OOT in planes;plane number, OOT index;x (mm)",
0390                                         1 + windowsNum * 4,
0391                                         -1. / windowsNum,
0392                                         4,
0393                                         19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0394                                         -0.5,
0395                                         18.5);
0396 
0397   {               // bin labelling (for clarity)
0398     int idx = 2;  // start counting at 1, first bin is empty
0399     for (int pl = 0; pl < 4; ++pl)
0400       for (unsigned int oot = 0; oot < windowsNum; ++oot) {
0401         const std::string bin_label =
0402             (oot == 0 ? "Plane " + std::to_string(pl) + ", " : "") + "OOT" + std::to_string(oot);
0403         hitDistribution2dOOT->setBinLabel(idx, bin_label);
0404         if (plotOnline)
0405           hitDistribution2dOOT_le->setBinLabel(idx, bin_label);
0406         ++idx;
0407       }
0408   }
0409 
0410   recHitTime = ibooker.book1D("recHit time", title + " recHit time; t (ns)", 500, -25, 25);
0411 
0412   activePlanes =
0413       ibooker.book1D("active planes", title + " active planes (per event);number of active planes", 6, -0.5, 5.5);
0414 
0415   trackDistribution =
0416       ibooker.book1D("tracks", title + " tracks;x (mm)", 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM, -0.5, 18.5);
0417   trackDistributionOOT = ibooker.book2D("tracks with OOT",
0418                                         title + " tracks with OOT;plane number;x (mm)",
0419                                         9,
0420                                         -0.5,
0421                                         4,
0422                                         19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0423                                         -0.5,
0424                                         18.5);
0425 
0426   leadingEdgeCumulative_both = ibooker.book1D("leading edge (le and te)",
0427                                               title + " leading edge (le and te) (recHits); leading edge (ns)",
0428                                               25 * windowsNum,
0429                                               0,
0430                                               25 * windowsNum);
0431   leadingEdgeCumulative_all = ibooker.book1D("leading edge (all)",
0432                                              title + " leading edge (with or without te) (DIGIs); leading edge (ns)",
0433                                              25 * windowsNum,
0434                                              0,
0435                                              25 * windowsNum);
0436   leadingEdgeCumulative_le = ibooker.book1D("leading edge (le only)",
0437                                             title + " leading edge (le only) (DIGIs); leading edge (ns)",
0438                                             25 * windowsNum,
0439                                             0,
0440                                             25 * windowsNum);
0441   trailingEdgeCumulative_te = ibooker.book1D("trailing edge (te only)",
0442                                              title + " trailing edge (te only) (DIGIs); trailing edge (ns)",
0443                                              25 * windowsNum,
0444                                              0,
0445                                              25 * windowsNum);
0446   timeOverThresholdCumulativePot =
0447       ibooker.book1D("time over threshold", title + " time over threshold;time over threshold (ns)", 250, -25, 100);
0448   // leadingTrailingCorrelationPot =
0449   //     ibooker.book2D("leading trailing correlation",
0450   //                    title + " leading trailing correlation;leading edge (ns);trailing edge (ns)",
0451   //                    75,
0452   //                    0,
0453   //                    75,
0454   //                    75,
0455   //                    0,
0456   //                    75);
0457 
0458   leadingWithoutTrailingCumulativePot =
0459       ibooker.book1D("event category", title + " leading edges without trailing;;%", 3, 0.5, 3.5);
0460   leadingWithoutTrailingCumulativePot->setBinLabel(1, "Leading only");
0461   leadingWithoutTrailingCumulativePot->setBinLabel(2, "Trailing only");
0462   leadingWithoutTrailingCumulativePot->setBinLabel(3, "Both");
0463 
0464   HPTDCErrorFlags_2D = ibooker.book2D("HPTDC Errors", title + " HPTDC Errors", 16, -0.5, 16.5, 9, -0.5, 8.5);
0465   for (unsigned short error_index = 1; error_index < 16; ++error_index)
0466     HPTDCErrorFlags_2D->setBinLabel(error_index, HPTDCErrorFlags::hptdcErrorName(error_index - 1));
0467   HPTDCErrorFlags_2D->setBinLabel(16, "Wrong EC");
0468 
0469   int tmpIndex = 0;
0470   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 0 TDC 18", /* axis */ 2);
0471   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 0 TDC 17", /* axis */ 2);
0472   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 0 TDC 16", /* axis */ 2);
0473   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 0 TDC 15", /* axis */ 2);
0474   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 1 TDC 18", /* axis */ 2);
0475   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 1 TDC 17", /* axis */ 2);
0476   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 1 TDC 16", /* axis */ 2);
0477   HPTDCErrorFlags_2D->setBinLabel(++tmpIndex, "DB 1 TDC 15", /* axis */ 2);
0478 
0479   MHComprensive =
0480       ibooker.book2D("MH in channels", title + " MH (%) in channels;plane number;ch number", 10, -0.5, 4.5, 14, -1, 13);
0481 
0482   // ibooker.setCurrentFolder( path+"/clock/" );
0483   // clock_Digi1_le = ibooker.book1D( "clock1 leading edge", title+" clock1;leading edge (ns)", 250, 0, 25 );
0484   // clock_Digi1_te = ibooker.book1D( "clock1 trailing edge", title+" clock1;trailing edge (ns)", 75, 0, 75 );
0485   // clock_Digi3_le = ibooker.book1D( "clock3 leading edge", title+" clock3;leading edge (ns)", 250, 0, 25 );
0486   // clock_Digi3_te = ibooker.book1D( "clock3 trailing edge", title+" clock3;trailing edge (ns)", 75, 0, 75 );
0487 }
0488 
0489 //----------------------------------------------------------------------------------------------------
0490 
0491 CTPPSDiamondDQMSource::PlanePlots::PlanePlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int windowsNum)
0492     : pixelTracksMapWithDiamonds("Pixel track maps for efficiency with coincidence",
0493                                  "Pixel track maps for efficiency with coincidence",
0494                                  25,
0495                                  0,
0496                                  25,
0497                                  12,
0498                                  -2,
0499                                  10) {
0500   std::string path, title;
0501   CTPPSDiamondDetId(id).planeName(path, CTPPSDiamondDetId::nPath);
0502   ibooker.setCurrentFolder(path);
0503 
0504   CTPPSDiamondDetId(id).planeName(title, CTPPSDiamondDetId::nFull);
0505 
0506   digiProfileCumulativePerPlane = ibooker.book1D("digi profile",
0507                                                  title + " digi profile; ch number",
0508                                                  CTPPS_DIAMOND_NUM_OF_CHANNELS,
0509                                                  -0.5,
0510                                                  CTPPS_DIAMOND_NUM_OF_CHANNELS - 0.5);
0511   hitProfile = ibooker.book1D(
0512       "hit profile", title + " hit profile;x (mm)", 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM, -0.5, 18.5);
0513   hit_multiplicity = ibooker.book1D("channels per plane", title + " channels per plane; ch per plane", 13, -0.5, 12.5);
0514 
0515   pixelTomography_far = ibooker.book2D("tomography pixel",
0516                                        title + " tomography with pixel;x + 25 OOT (mm);y (mm)",
0517                                        25 * windowsNum,
0518                                        0,
0519                                        25 * windowsNum,
0520                                        8,
0521                                        0,
0522                                        8);
0523   EfficiencyWRTPixelsInPlane =
0524       ibooker.book2D("Efficiency wrt pixels", title + " Efficiency wrt pixels;x (mm);y (mm)", 25, 0, 25, 12, -2, 10);
0525 }
0526 
0527 //----------------------------------------------------------------------------------------------------
0528 
0529 CTPPSDiamondDQMSource::ChannelPlots::ChannelPlots(DQMStore::IBooker& ibooker, unsigned int id, unsigned int windowsNum)
0530     : HitCounter(0), MHCounter(0), LeadingOnlyCounter(0), TrailingOnlyCounter(0), CompleteCounter(0) {
0531   std::string path, title;
0532   CTPPSDiamondDetId(id).channelName(path, CTPPSDiamondDetId::nPath);
0533   ibooker.setCurrentFolder(path);
0534 
0535   CTPPSDiamondDetId(id).channelName(title, CTPPSDiamondDetId::nFull);
0536 
0537   leadingWithoutTrailing = ibooker.book1D("event category", title + " Event Category;;%", 3, 0.5, 3.5);
0538   leadingWithoutTrailing->setBinLabel(1, "Leading only");
0539   leadingWithoutTrailing->setBinLabel(2, "Trailing only");
0540   leadingWithoutTrailing->setBinLabel(3, "Full");
0541 
0542   for (unsigned int i = 0; i < windowsNum; i++) {
0543     std::string window = std::to_string(i * 25) + "-" + std::to_string((i + 1) * 25);
0544     activity_per_bx[i] = ibooker.book1D(
0545         "activity per BX " + window, title + " Activity per BX " + window + " ns;Event.BX", 3600, -1.5, 3598. + 0.5);
0546   }
0547 
0548   HPTDCErrorFlags = ibooker.book1D("hptdc_Errors", title + " HPTDC Errors", 16, -0.5, 16.5);
0549   for (unsigned short error_index = 1; error_index < 16; ++error_index)
0550     HPTDCErrorFlags->setBinLabel(error_index, HPTDCErrorFlags::hptdcErrorName(error_index - 1));
0551   HPTDCErrorFlags->setBinLabel(16, "MH  (%)");
0552 
0553   leadingEdgeCumulative_both = ibooker.book1D("leading edge (le and te)",
0554                                               title + " leading edge (recHits); leading edge (ns)",
0555                                               25 * windowsNum,
0556                                               0,
0557                                               25 * windowsNum);
0558   leadingEdgeCumulative_le = ibooker.book1D(
0559       "leading edge (le only)", title + " leading edge (DIGIs); leading edge (ns)", 25 * windowsNum, 0, 25 * windowsNum);
0560   trailingEdgeCumulative_te = ibooker.book1D("trailing edge (te only)",
0561                                              title + " trailing edge (te only) (DIGIs); trailing edge (ns)",
0562                                              25 * windowsNum,
0563                                              0,
0564                                              25 * windowsNum);
0565   TimeOverThresholdCumulativePerChannel =
0566       ibooker.book1D("time over threshold", title + " time over threshold;time over threshold (ns)", 75, -25, 50);
0567   // LeadingTrailingCorrelationPerChannel =
0568   //     ibooker.book2D("leading trailing correlation",
0569   //                    title + " leading trailing correlation;leading edge (ns);trailing edge (ns)",
0570   //                    75,
0571   //                    0,
0572   //                    75,
0573   //                    75,
0574   //                    0,
0575   //                    75);
0576 
0577   pixelTomography_far = ibooker.book2D(
0578       "tomography pixel", "tomography with pixel;x + 25 OOT (mm);y (mm)", 25 * windowsNum, 0, 25 * windowsNum, 8, 0, 8);
0579 
0580   hit_rate = ibooker.book1D("hit rate", title + "hit rate;rate (Hz)", 40, 0, 20);
0581 
0582   recHitTime = ibooker.book1D("recHit Time", title + " recHit Time; t (ns)", 500, -25, 25);
0583 }
0584 
0585 //----------------------------------------------------------------------------------------------------
0586 
0587 CTPPSDiamondDQMSource::CTPPSDiamondDQMSource(const edm::ParameterSet& ps)
0588     : perLSsaving_(ps.getUntrackedParameter<bool>("perLSsaving", false)),
0589       tokenPixelTrack_(consumes<edm::DetSetVector<CTPPSPixelLocalTrack>>(
0590           ps.getUntrackedParameter<edm::InputTag>("tagPixelLocalTracks"))),
0591       tokenDiamondHit_(consumes<edm::DetSetVector<CTPPSDiamondRecHit>>(
0592           ps.getUntrackedParameter<edm::InputTag>("tagDiamondRecHits"))),
0593       tokenDiamondTrack_(consumes<edm::DetSetVector<CTPPSDiamondLocalTrack>>(
0594           ps.getUntrackedParameter<edm::InputTag>("tagDiamondLocalTracks"))),
0595       ctppsGeometryRunToken_(esConsumes<CTPPSGeometry, VeryForwardRealGeometryRecord, edm::Transition::BeginRun>()),
0596       ctppsGeometryEventToken_(esConsumes<CTPPSGeometry, VeryForwardRealGeometryRecord>()),
0597       // ctppsLhcInfoToken_(esConsumes<LHCInfo, LHCInfoRcd>()),
0598       excludeMultipleHits_(ps.getParameter<bool>("excludeMultipleHits")),
0599       extract_digi_info_(ps.getParameter<bool>("extractDigiInfo")),
0600       centralOOT_(-999),
0601       verbosity_(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
0602       plotOnline_(ps.getUntrackedParameter<bool>("plotOnline", true)),
0603       plotOffline_(ps.getUntrackedParameter<bool>("plotOffline", false)),
0604       windowsNum_(ps.getUntrackedParameter<unsigned int>("windowsNum", 3)),
0605       trackCorrelationThreshold_(ps.getUntrackedParameter<unsigned int>("trackCorrelationThreshold", 3)),
0606       EC_difference_56_(-500),
0607       EC_difference_45_(-500) {
0608   if (extract_digi_info_) {
0609     tokenStatus_ = consumes<edm::DetSetVector<TotemVFATStatus>>(ps.getUntrackedParameter<edm::InputTag>("tagStatus"));
0610     tokenFEDInfo_ = consumes<std::vector<TotemFEDInfo>>(ps.getUntrackedParameter<edm::InputTag>("tagFEDInfo"));
0611     tokenDigi_ = consumes<edm::DetSetVector<CTPPSDiamondDigi>>(ps.getUntrackedParameter<edm::InputTag>("tagDigi"));
0612   }
0613   for (const auto& pset : ps.getParameter<std::vector<edm::ParameterSet>>("offsetsOOT")) {
0614     runParameters_.emplace_back(
0615         std::make_pair(pset.getParameter<edm::EventRange>("validityRange"), pset.getParameter<int>("centralOOT")));
0616   }
0617 }
0618 
0619 //----------------------------------------------------------------------------------------------------
0620 
0621 void CTPPSDiamondDQMSource::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0622   centralOOT_ = -999;
0623   for (const auto& oot : runParameters_) {
0624     if (edm::contains(oot.first, edm::EventID(iRun.run(), 0, 1))) {
0625       centralOOT_ = oot.second;
0626       break;
0627     }
0628   }
0629 
0630   // Get detector shifts from the geometry
0631   const CTPPSGeometry& geom = iSetup.getData(ctppsGeometryRunToken_);
0632   for (auto it = geom.beginRP(); it != geom.endRP(); ++it)
0633     if (CTPPSDiamondDetId::check(it->first)) {
0634       const CTPPSDiamondDetId diam_id(it->first);
0635       const auto diam = geom.sensor(it->first);
0636       diamShifts_[diam_id].global = diam->translation().x() - diam->getDiamondDimensions().xHalfWidth;
0637       if (iRun.run() > FIRST_RUN_W_PIXELS) {  // pixel installed
0638         const CTPPSPixelDetId pixid(diam_id.arm(), CTPPS_PIXEL_STATION_ID, CTPPS_FAR_RP_ID);
0639         auto pix = geom.sensor(pixid);
0640         // Rough alignement of pixel detector for diamond tomography
0641         diamShifts_[diam_id].withPixels =
0642             pix->translation().x() - pix->getDiamondDimensions().xHalfWidth - diamShifts_[diam_id].global - 1.;
0643       }
0644     }
0645 }
0646 
0647 //----------------------------------------------------------------------------------------------------
0648 
0649 void CTPPSDiamondDQMSource::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup& iSetup) {
0650   ibooker.cd();
0651   ibooker.setCurrentFolder("CTPPS");
0652 
0653   globalPlot_ = GlobalPlots(ibooker);
0654 
0655   // book plots from the geometry
0656   const CTPPSGeometry& geom = iSetup.getData(ctppsGeometryRunToken_);
0657   for (auto it = geom.beginSensor(); it != geom.endSensor(); ++it) {
0658     if (!CTPPSDiamondDetId::check(it->first))
0659       continue;
0660     // per-channel plots
0661     const CTPPSDiamondDetId chId(it->first);
0662     if (plotOnline_ && channelPlots_.count(chId) == 0)
0663       channelPlots_[chId] = ChannelPlots(ibooker, chId, windowsNum_);
0664 
0665     // per-plane plots
0666     const CTPPSDiamondDetId plId(chId.planeId());
0667     if (planePlots_.count(plId) == 0)
0668       planePlots_[plId] = PlanePlots(ibooker, plId, windowsNum_);
0669     // per-pot plots
0670     const CTPPSDiamondDetId rpId(chId.rpId());
0671     if (potPlots_.count(rpId) == 0)
0672       potPlots_[rpId] = PotPlots(ibooker, rpId, windowsNum_, plotOnline_, plotOffline_, perLSsaving_);
0673 
0674     // per-sector plots
0675     const CTPPSDiamondDetId secId(chId.armId());
0676     if (plotOffline_ && sectorPlots_.count(secId) == 0)
0677       sectorPlots_[secId] = SectorPlots(ibooker, secId, plotOnline_);
0678   }
0679 }
0680 
0681 //----------------------------------------------------------------------------------------------------
0682 
0683 std::shared_ptr<dds::Cache> CTPPSDiamondDQMSource::globalBeginLuminosityBlock(const edm::LuminosityBlock&,
0684                                                                               const edm::EventSetup&) const {
0685   auto d = std::make_shared<dds::Cache>();
0686   d->hitDistribution2dMap.reserve(potPlots_.size());
0687   if (!perLSsaving_ && plotOnline_) {
0688     for (auto& plot : potPlots_) {
0689       d->hitDistribution2dMap[plot.first] = std::make_unique<TH2F>(
0690           "hits in planes lumisection",
0691           (std::string(plot.second.hitDistribution2d_lumisection->getTH2F()->GetTitle()) + ";plane number;x (mm)")
0692               .c_str(),
0693           10,
0694           -0.5,
0695           4.5,
0696           19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0697           -0.5,
0698           18.5);
0699     }
0700   }
0701   return d;
0702 }
0703 
0704 //----------------------------------------------------------------------------------------------------
0705 
0706 void CTPPSDiamondDQMSource::analyze(const edm::Event& event, const edm::EventSetup& iSetup) {
0707   // get event data
0708 
0709   edm::Handle<edm::DetSetVector<TotemVFATStatus>> diamondVFATStatus;
0710   edm::Handle<edm::DetSetVector<CTPPSDiamondDigi>> diamondDigis;
0711   edm::Handle<std::vector<TotemFEDInfo>> fedInfo;
0712   if (extract_digi_info_) {
0713     event.getByToken(tokenStatus_, diamondVFATStatus);
0714     event.getByToken(tokenDigi_, diamondDigis);
0715     event.getByToken(tokenFEDInfo_, fedInfo);
0716   }
0717 
0718   edm::Handle<edm::DetSetVector<CTPPSPixelLocalTrack>> pixelTracks;
0719   event.getByToken(tokenPixelTrack_, pixelTracks);
0720 
0721   edm::Handle<edm::DetSetVector<CTPPSDiamondRecHit>> diamondRecHits;
0722   event.getByToken(tokenDiamondHit_, diamondRecHits);
0723 
0724   edm::Handle<edm::DetSetVector<CTPPSDiamondLocalTrack>> diamondLocalTracks;
0725   event.getByToken(tokenDiamondTrack_, diamondLocalTracks);
0726 
0727   const CTPPSGeometry* ctppsGeometry = &iSetup.getData(ctppsGeometryEventToken_);
0728   // const LHCInfo* hLhcInfo = &iSetup.getData(ctppsLhcInfoToken_);
0729 
0730   // check validity
0731   bool valid = true;
0732   if (extract_digi_info_) {  // drop DIGI-level validity checks if not monitored
0733     valid &= diamondVFATStatus.isValid();
0734     valid &= diamondDigis.isValid();
0735     valid &= fedInfo.isValid();
0736   }
0737   valid &= pixelTracks.isValid();
0738   valid &= diamondRecHits.isValid();
0739   valid &= diamondLocalTracks.isValid();
0740 
0741   if (!valid) {
0742     if (verbosity_)
0743       edm::LogProblem("CTPPSDiamondDQMSource")
0744           << "ERROR in CTPPSDiamondDQMSource::analyze > some of the required inputs are not valid. Skipping this "
0745              "event.\n"
0746           << "  DIGI-level: (checked? " << std::boolalpha << extract_digi_info_ << ")\n"
0747           << "    diamondVFATStatus.isValid = " << diamondVFATStatus.isValid() << "\n"
0748           << "    diamondDigis.isValid = " << diamondDigis.isValid() << "\n"
0749           << "    fedInfo.isValid = " << fedInfo.isValid() << "\n"
0750           << "  RECO-level:\n"
0751           << "    pixelTracks.isValid = " << pixelTracks.isValid() << "\n"
0752           << "    diamondRecHits.isValid = " << diamondRecHits.isValid() << "\n"
0753           << "    diamondLocalTracks.isValid = " << diamondLocalTracks.isValid();
0754     return;
0755   }
0756 
0757   //------------------------------
0758   // Sector Plots
0759   //------------------------------
0760   // Using CTPPSDiamondLocalTrack
0761   if (plotOffline_)
0762     for (const auto& tracks : *diamondLocalTracks) {
0763       const CTPPSDiamondDetId detId_near(tracks.detId());
0764 
0765       for (const auto& track : tracks) {
0766         if (!track.isValid())
0767           continue;
0768         if (potPlots_.count(detId_near.rpId()) == 0)
0769           continue;
0770         TH1F* trackHistoInTimeTmp = potPlots_[detId_near.rpId()].trackDistribution->getTH1F();
0771         int startBin_near =
0772             trackHistoInTimeTmp->FindBin(track.x0() - diamShifts_[detId_near.rpId()].global - track.x0Sigma());
0773         int numOfBins_near = 2 * track.x0Sigma() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0774 
0775         for (const auto& tracks_far : *diamondLocalTracks) {
0776           CTPPSDiamondDetId detId_far(tracks_far.detId());
0777           if (detId_near.arm() != detId_far.arm() || detId_near.station() == detId_far.station())
0778             continue;
0779           for (const auto& track_far : tracks_far) {
0780             if (!track.isValid())
0781               continue;
0782             if (sectorPlots_.count(detId_far.armId()) == 0)
0783               continue;
0784             TH2F* trackHistoTmp = sectorPlots_[detId_far.armId()].trackCorrelation->getTH2F();
0785             TAxis* trackHistoTmpXAxis = trackHistoTmp->GetXaxis();
0786             TAxis* trackHistoTmpYAxis = trackHistoTmp->GetYaxis();
0787             int startBin_far = trackHistoTmpYAxis->FindBin(track_far.x0() - diamShifts_[detId_far.rpId()].global -
0788                                                            track_far.x0Sigma());
0789             int numOfBins_far = 2 * track_far.x0Sigma() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0790             for (int i = 0; i < numOfBins_near; ++i)
0791               for (int y = 0; y < numOfBins_far; ++y) {
0792                 trackHistoTmp->Fill(trackHistoTmpXAxis->GetBinCenter(startBin_near + i),
0793                                     trackHistoTmpYAxis->GetBinCenter(startBin_far + y));
0794                 if (tracks.size() < 3 && tracks_far.size() < trackCorrelationThreshold_)
0795                   sectorPlots_[detId_far.armId()].trackCorrelationLowMultiplicity->Fill(
0796                       trackHistoTmpXAxis->GetBinCenter(startBin_near + i),
0797                       trackHistoTmpYAxis->GetBinCenter(startBin_far + y));
0798               }
0799           }
0800         }
0801       }
0802     }
0803 
0804   //------------------------------
0805   // RP Plots
0806   //------------------------------
0807 
0808   //------------------------------
0809   // Correlation Plots
0810   //------------------------------
0811 
0812   if (extract_digi_info_) {
0813     // Using CTPPSDiamondDigi
0814     for (const auto& digis : *diamondDigis) {
0815       const CTPPSDiamondDetId detId(digis.detId()), detId_pot(detId.rpId());
0816       for (const auto& digi : digis) {
0817         if (detId.channel() == CHANNEL_OF_VFAT_CLOCK)
0818           continue;
0819         if (potPlots_.count(detId_pot) == 0)
0820           continue;
0821         //Leading without trailing investigation
0822         if (digi.leadingEdge() != 0 || digi.trailingEdge() != 0) {
0823           ++potPlots_[detId_pot].HitCounter;
0824           if (digi.leadingEdge() != 0) {
0825             potPlots_[detId_pot].leadingEdgeCumulative_all->Fill(HPTDC_BIN_WIDTH_NS * digi.leadingEdge());
0826           }
0827           if (digi.leadingEdge() != 0 && digi.trailingEdge() == 0) {
0828             ++potPlots_[detId_pot].LeadingOnlyCounter;
0829             potPlots_[detId_pot].leadingEdgeCumulative_le->Fill(HPTDC_BIN_WIDTH_NS * digi.leadingEdge());
0830           }
0831           if (digi.leadingEdge() == 0 && digi.trailingEdge() != 0) {
0832             ++potPlots_[detId_pot].TrailingOnlyCounter;
0833             potPlots_[detId_pot].trailingEdgeCumulative_te->Fill(HPTDC_BIN_WIDTH_NS * digi.trailingEdge());
0834           }
0835           if (digi.leadingEdge() != 0 && digi.trailingEdge() != 0) {
0836             ++potPlots_[detId_pot].CompleteCounter;
0837             // potPlots_[detId_pot].leadingTrailingCorrelationPot->Fill(HPTDC_BIN_WIDTH_NS * digi.leadingEdge(),
0838             //                                                          HPTDC_BIN_WIDTH_NS * digi.trailingEdge());
0839           }
0840         }
0841 
0842         // HPTDC Errors
0843         const HPTDCErrorFlags hptdcErrors = digi.hptdcErrorFlags();
0844         if (detId.channel() == HPTDC_0_CHANNEL ||
0845             detId.channel() == HPTDC_1_CHANNEL) {  // ch6 for HPTDC 0 and ch7 for HPTDC 1
0846           int verticalIndex = 2 * detId.plane() + (detId.channel() - HPTDC_0_CHANNEL);
0847           for (unsigned short hptdcErrorIndex = 1; hptdcErrorIndex < 16; ++hptdcErrorIndex)
0848             if (hptdcErrors.errorId(hptdcErrorIndex - 1))
0849               potPlots_[detId_pot].HPTDCErrorFlags_2D->Fill(hptdcErrorIndex, verticalIndex);
0850         }
0851         if (digi.multipleHit())
0852           ++potPlots_[detId_pot].MHCounter;
0853       }
0854     }
0855   }
0856 
0857   // EC Errors
0858   if (extract_digi_info_) {
0859     for (const auto& vfat_status : *diamondVFATStatus) {
0860       const CTPPSDiamondDetId detId(vfat_status.detId());
0861       for (const auto& status : vfat_status) {
0862         if (!status.isOK())
0863           continue;
0864         if (potPlots_.count(detId.rpId()) == 0)
0865           continue;
0866         if (channelPlots_.count(detId) == 0)
0867           continue;
0868 
0869         // Check Event Number
0870         for (const auto& optorx : *fedInfo) {
0871           if (detId.arm() == 1 && optorx.fedId() == CTPPS_FED_ID_56)
0872             checkEventNumber(detId, optorx, status, potPlots_[detId.rpId()], EC_difference_56_);
0873           else if (detId.arm() == 0 && optorx.fedId() == CTPPS_FED_ID_45)
0874             checkEventNumber(detId, optorx, status, potPlots_[detId.rpId()], EC_difference_45_);
0875         }
0876       }
0877     }
0878   }
0879 
0880   // Using CTPPSDiamondRecHit
0881   std::unordered_map<unsigned int, std::set<unsigned int>> planes;
0882   std::unordered_map<unsigned int, std::set<unsigned int>> planes_inclusive;
0883 
0884   auto lumiCache = luminosityBlockCache(event.getLuminosityBlock().index());
0885   for (const auto& rechits : *diamondRecHits) {
0886     const CTPPSDiamondDetId detId(rechits.detId()), detId_pot(detId.rpId());
0887     const auto& x_shift = diamShifts_.at(detId_pot);
0888 
0889     for (const auto& rechit : rechits) {
0890       planes_inclusive[detId_pot].insert(detId.plane());
0891       if (excludeMultipleHits_ && rechit.multipleHits() > 0)
0892         continue;
0893       if (rechit.toT() != 0 && centralOOT_ != -999 && rechit.ootIndex() == centralOOT_)
0894         planes[detId_pot].insert(detId.plane());
0895 
0896       if (potPlots_.count(detId_pot) == 0)
0897         continue;
0898 
0899       potPlots_[detId_pot].recHitTime->Fill(rechit.time());
0900 
0901       float UFSDShift = 0.0;
0902       if (rechit.yWidth() < 3)
0903         UFSDShift = 0.5;  // Display trick for UFSD that have 2 pixels with same X
0904 
0905       if (rechit.toT() != 0 && centralOOT_ != -999 && rechit.ootIndex() == centralOOT_) {
0906         TH2F* hitHistoTmp = potPlots_[detId_pot].hitDistribution2d->getTH2F();
0907         TAxis* hitHistoTmpYAxis = hitHistoTmp->GetYaxis();
0908         int startBin = hitHistoTmpYAxis->FindBin(rechit.x() - x_shift.global - 0.5 * rechit.xWidth());
0909         int numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0910         for (int i = 0; i < numOfBins; ++i)
0911           hitHistoTmp->Fill(detId.plane() + UFSDShift, hitHistoTmpYAxis->GetBinCenter(startBin + i));
0912 
0913         if (!perLSsaving_ && plotOnline_) {
0914           hitHistoTmp = lumiCache->hitDistribution2dMap[detId_pot].get();
0915           hitHistoTmpYAxis = hitHistoTmp->GetYaxis();
0916           startBin = hitHistoTmpYAxis->FindBin(rechit.x() - x_shift.global - 0.5 * rechit.xWidth());
0917           numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0918           for (int i = 0; i < numOfBins; ++i)
0919             hitHistoTmp->Fill(detId.plane() + UFSDShift, hitHistoTmpYAxis->GetBinCenter(startBin + i));
0920         }
0921       }
0922 
0923       if (rechit.toT() > 0) {
0924         // Both
0925         potPlots_[detId_pot].leadingEdgeCumulative_both->Fill(rechit.time() + 25 * rechit.ootIndex());
0926         potPlots_[detId_pot].timeOverThresholdCumulativePot->Fill(rechit.toT());
0927 
0928         TH2F* hitHistoOOTTmp = potPlots_[detId_pot].hitDistribution2dOOT->getTH2F();
0929         TAxis* hitHistoOOTTmpYAxis = hitHistoOOTTmp->GetYaxis();
0930         int startBin = hitHistoOOTTmpYAxis->FindBin(rechit.x() - x_shift.global - 0.5 * rechit.xWidth());
0931         int numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0932         for (int i = 0; i < numOfBins; ++i)
0933           hitHistoOOTTmp->Fill(detId.plane() + 1. / windowsNum_ * rechit.ootIndex(),
0934                                hitHistoOOTTmpYAxis->GetBinCenter(startBin + i));
0935 
0936       } else if (rechit.ootIndex() != CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING && plotOnline_) {
0937         // Only leading
0938         TH2F* hitHistoOOTTmp = potPlots_[detId_pot].hitDistribution2dOOT_le->getTH2F();
0939         TAxis* hitHistoOOTTmpYAxis = hitHistoOOTTmp->GetYaxis();
0940         int startBin = hitHistoOOTTmpYAxis->FindBin(rechit.x() - x_shift.global - 0.5 * rechit.xWidth());
0941         int numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0942         for (int i = 0; i < numOfBins; ++i)
0943           hitHistoOOTTmp->Fill(detId.plane() + 1. / windowsNum_ * rechit.ootIndex(),
0944                                hitHistoOOTTmpYAxis->GetBinCenter(startBin + i));
0945       }
0946       if (rechit.ootIndex() != CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING &&
0947           potPlots_[detId_pot].activity_per_bx.count(rechit.ootIndex()) > 0)
0948         potPlots_[detId_pot].activity_per_bx.at(rechit.ootIndex())->Fill(event.bunchCrossing());
0949 
0950       // if(plotOffline_)
0951       //   potPlots_[detId_pot].TOTVsLS->Fill(event.luminosityBlock(),rechit.toT());
0952     }
0953   }
0954 
0955   for (const auto& plt : potPlots_) {
0956     plt.second.activePlanes->Fill(planes[plt.first].size());
0957     if (plotOnline_)
0958       plt.second.activePlanesInclusive->Fill(planes_inclusive[plt.first].size());
0959   }
0960 
0961   // Using CTPPSDiamondLocalTrack
0962   for (const auto& tracks : *diamondLocalTracks) {
0963     const CTPPSDiamondDetId detId(tracks.detId()), detId_pot(detId.rpId());
0964     const auto& x_shift = diamShifts_.at(detId_pot);
0965 
0966     for (const auto& track : tracks) {
0967       if (!track.isValid())
0968         continue;
0969       if (track.ootIndex() == CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING)
0970         continue;
0971       if (excludeMultipleHits_ && track.multipleHits() > 0)
0972         continue;
0973       if (potPlots_.count(detId_pot) == 0)
0974         continue;
0975 
0976       TH2F* trackHistoOOTTmp = potPlots_[detId_pot].trackDistributionOOT->getTH2F();
0977       TAxis* trackHistoOOTTmpYAxis = trackHistoOOTTmp->GetYaxis();
0978       int startBin = trackHistoOOTTmpYAxis->FindBin(track.x0() - x_shift.global - track.x0Sigma());
0979       int numOfBins = 2 * track.x0Sigma() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0980       for (int i = 0; i < numOfBins; ++i)
0981         trackHistoOOTTmp->Fill(track.ootIndex(), trackHistoOOTTmpYAxis->GetBinCenter(startBin + i));
0982 
0983       if (centralOOT_ != -999 && track.ootIndex() == centralOOT_) {
0984         TH1F* trackHistoInTimeTmp = potPlots_[detId_pot].trackDistribution->getTH1F();
0985         int startBin = trackHistoInTimeTmp->FindBin(track.x0() - x_shift.global - track.x0Sigma());
0986         int numOfBins = 2 * track.x0Sigma() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0987         for (int i = 0; i < numOfBins; ++i)
0988           trackHistoInTimeTmp->Fill(trackHistoInTimeTmp->GetBinCenter(startBin + i));
0989       }
0990       if (plotOffline_ && !perLSsaving_) {
0991         // potPlots_[detId_pot].trackTimeVsLS->Fill(event.luminosityBlock(),track.time());
0992         potPlots_[detId_pot].trackTimeVsBX->Fill(event.bunchCrossing(), track.time());
0993         //potPlots_[detId_pot].trackTimeVsXAngle->Fill(hLhcInfo->crossingAngle(), track.time());
0994       }
0995     }
0996   }
0997 
0998   // Channel efficiency using CTPPSDiamondLocalTrack
0999   for (const auto& tracks : *diamondLocalTracks) {
1000     const CTPPSDiamondDetId detId(tracks.detId()), detId_pot(detId.rpId());
1001     for (const auto& track : tracks) {
1002       // Find hits and planes in the track
1003       int numOfHits = 0;
1004       std::set<int> planesInTrackSet;
1005       for (const auto& vec : *diamondRecHits) {
1006         const CTPPSDiamondDetId detid(vec.detId());
1007         if (detid.arm() != detId_pot.arm())
1008           continue;
1009 
1010         for (const auto& hit : vec) {
1011           // first check if the hit contributes to the track
1012           if (track.containsHit(hit)) {
1013             ++numOfHits;
1014             planesInTrackSet.insert(detid.plane());
1015           }
1016         }
1017       }
1018 
1019       if (numOfHits > 0 && numOfHits <= 10 && planesInTrackSet.size() > 2) {
1020         for (const auto& plt_vs_ch : channelPlots_) {  // loop over all channels registered in geometry
1021           const CTPPSDiamondDetId detId(plt_vs_ch.first);
1022           if (detId.rpId() != detId_pot)
1023             continue;
1024           const unsigned short map_index = detId.plane() * 100 + detId.channel();
1025           if (potPlots_[detId_pot].effDoublecountingChMap.count(map_index) == 0) {
1026             potPlots_[detId_pot].effTriplecountingChMap[map_index] = 0;
1027             potPlots_[detId_pot].effDoublecountingChMap[map_index] = 0;
1028           }
1029           if (channelAlignedWithTrack(ctppsGeometry, detId, track, 0.2)) {
1030             // Channel should fire
1031             ++potPlots_[detId_pot].effDoublecountingChMap[map_index];
1032             for (const auto& rechits : *diamondRecHits) {
1033               const CTPPSDiamondDetId detId_hit(rechits.detId());
1034               if (detId_hit == detId)
1035                 for (const auto& rechit : rechits)
1036                   if (track.containsHit(rechit, 1))  // Channel fired
1037                     ++potPlots_[detId_pot].effTriplecountingChMap[map_index];
1038             }
1039           }
1040         }
1041       }
1042     }
1043   }
1044 
1045   // Tomography of diamonds using pixel
1046   for (const auto& rechits : *diamondRecHits) {
1047     const CTPPSDiamondDetId detId(rechits.detId()), detId_pot(detId.rpId());
1048     const auto pix_shift = diamShifts_.at(detId_pot).withPixels;
1049     for (const auto& rechit : rechits) {
1050       if (excludeMultipleHits_ && rechit.multipleHits() > 0)
1051         continue;
1052       if (rechit.toT() == 0)
1053         continue;
1054       if (!pixelTracks.isValid())
1055         continue;
1056       if (potPlots_.count(detId_pot) == 0)
1057         continue;
1058 
1059       for (const auto& ds : *pixelTracks) {
1060         if (ds.size() > 1)
1061           continue;
1062         const CTPPSPixelDetId pixId(ds.detId());
1063         if (pixId.station() != CTPPS_PIXEL_STATION_ID || pixId.rp() != CTPPS_FAR_RP_ID)
1064           continue;
1065         for (const auto& lt : ds) {
1066           if (lt.isValid() && pixId.arm() == detId_pot.arm()) {
1067             if (rechit.ootIndex() != CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING && rechit.ootIndex() >= 0 &&
1068                 potPlots_[detId_pot].pixelTomographyAll.count(rechit.ootIndex()) > 0 && lt.x0() - pix_shift < 24)
1069               potPlots_[detId_pot]
1070                   .pixelTomographyAll.at(rechit.ootIndex())
1071                   ->Fill(lt.x0() - pix_shift + 25 * detId.plane(), lt.y0());
1072           }
1073         }
1074       }
1075     }
1076   }
1077 
1078   //------------------------------
1079   // Clock Plots
1080   //------------------------------
1081   // Commented out to save space in the DQM files, but code should be kept
1082   // for ( const auto& digis : *diamondDigis ) {
1083   //   const CTPPSDiamondDetId detId(digis.detId()), detId_pot(detId.rpId());
1084   //   if ( detId.channel() == CHANNEL_OF_VFAT_CLOCK ) {
1085   //     for ( const auto& digi : digis ) {
1086   //       if ( digi.leadingEdge() != 0 )  {
1087   //         if ( detId.plane() == 1 ) {
1088   //           potPlots_[detId_pot].clock_Digi1_le->Fill( HPTDC_BIN_WIDTH_NS * digi.leadingEdge() );
1089   //           potPlots_[detId_pot].clock_Digi1_te->Fill( HPTDC_BIN_WIDTH_NS * digi.trailingEdge() );
1090   //         }
1091   //         if ( detId.plane() == 3 ) {
1092   //           potPlots_[detId_pot].clock_Digi3_le->Fill( HPTDC_BIN_WIDTH_NS * digi.leadingEdge() );
1093   //           potPlots_[detId_pot].clock_Digi3_te->Fill( HPTDC_BIN_WIDTH_NS * digi.trailingEdge() );
1094   //         }
1095   //       }
1096   //     }
1097   //   }
1098   // }
1099 
1100   //------------------------------
1101   // Plane Plots
1102   //------------------------------
1103 
1104   // Using CTPPSDiamondDigi
1105   std::unordered_map<unsigned int, unsigned int> channelsPerPlane;
1106   if (extract_digi_info_) {
1107     for (const auto& digis : *diamondDigis) {
1108       const CTPPSDiamondDetId detId(digis.detId()), detId_plane(detId.planeId());
1109       for (const auto& digi : digis) {
1110         if (detId.channel() == CHANNEL_OF_VFAT_CLOCK)
1111           continue;
1112         if (planePlots_.count(detId_plane) == 0)
1113           continue;
1114 
1115         if (digi.leadingEdge() != 0) {
1116           planePlots_[detId_plane].digiProfileCumulativePerPlane->Fill(detId.channel());
1117           channelsPerPlane[detId_plane]++;
1118         }
1119       }
1120     }
1121   }
1122 
1123   for (const auto& plt : channelsPerPlane)
1124     planePlots_[plt.first].hit_multiplicity->Fill(plt.second);
1125 
1126   // Using CTPPSDiamondRecHit
1127   for (const auto& rechits : *diamondRecHits) {
1128     const CTPPSDiamondDetId detId(rechits.detId()), detId_plane(detId.planeId()), detId_pot(detId.rpId());
1129     for (const auto& rechit : rechits) {
1130       if (excludeMultipleHits_ && rechit.multipleHits() > 0)
1131         continue;
1132       if (rechit.toT() == 0)
1133         continue;
1134       if (planePlots_.count(detId_plane) != 0) {
1135         if (centralOOT_ != -999 && rechit.ootIndex() == centralOOT_) {
1136           TH1F* hitHistoTmp = planePlots_[detId_plane].hitProfile->getTH1F();
1137           int startBin = hitHistoTmp->FindBin(rechit.x() - diamShifts_.at(detId_pot).global - 0.5 * rechit.xWidth());
1138           int numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
1139           for (int i = 0; i < numOfBins; ++i)
1140             hitHistoTmp->Fill(hitHistoTmp->GetBinCenter(startBin + i));
1141         }
1142       }
1143     }
1144   }
1145 
1146   //Tomography of diamonds using pixel and Efficiency WRT Pixels
1147   for (const auto& ds : *pixelTracks) {
1148     const CTPPSPixelDetId pixId(ds.detId());
1149     if (pixId.station() != CTPPS_PIXEL_STATION_ID || pixId.rp() != CTPPS_FAR_RP_ID)
1150       continue;
1151     if (ds.size() > 1)
1152       continue;
1153     for (const auto& lt : ds) {
1154       if (lt.isValid()) {
1155         for (const auto& sh_vs_id : diamShifts_) {
1156           const CTPPSDiamondDetId& detId_pot = sh_vs_id.first;
1157           const auto& shift = sh_vs_id.second;
1158           if (detId_pot.arm() == pixId.arm())
1159             // For efficieny
1160             potPlots_[detId_pot].pixelTracksMap.Fill(lt.x0() - shift.withPixels, lt.y0());
1161         }
1162 
1163         std::set<CTPPSDiamondDetId> planesWitHits_set;
1164         for (const auto& rechits : *diamondRecHits) {
1165           const CTPPSDiamondDetId detId(rechits.detId()), detId_plane(detId.planeId()), detId_pot(detId.rpId());
1166           const auto pix_shift = diamShifts_.at(detId_pot).withPixels;
1167           for (const auto& rechit : rechits) {
1168             if (excludeMultipleHits_ && rechit.multipleHits() > 0)
1169               continue;
1170             if (rechit.ootIndex() == CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING || rechit.toT() == 0)
1171               continue;
1172             if (planePlots_.count(detId_plane) == 0)
1173               continue;
1174             if (pixId.arm() == detId_plane.arm() && lt.x0() - pix_shift < 24) {
1175               planePlots_[detId_plane].pixelTomography_far->Fill(lt.x0() - pix_shift + 25 * rechit.ootIndex(), lt.y0());
1176               if (centralOOT_ != -999 && rechit.ootIndex() == centralOOT_)
1177                 planesWitHits_set.insert(detId_plane);
1178             }
1179           }
1180         }
1181 
1182         for (auto& planeId : planesWitHits_set)
1183           planePlots_[planeId].pixelTracksMapWithDiamonds.Fill(lt.x0() - diamShifts_.at(planeId.rpId()).withPixels,
1184                                                                lt.y0());
1185       }
1186     }
1187   }
1188 
1189   //------------------------------
1190   // Channel Plots
1191   //------------------------------
1192 
1193   // digi profile cumulative
1194   if (extract_digi_info_) {
1195     for (const auto& digis : *diamondDigis) {
1196       const CTPPSDiamondDetId detId(digis.detId());
1197       for (const auto& digi : digis) {
1198         if (detId.channel() == CHANNEL_OF_VFAT_CLOCK)
1199           continue;
1200         if (channelPlots_.count(detId) != 0) {
1201           // HPTDC Errors
1202           const HPTDCErrorFlags hptdcErrors = digi.hptdcErrorFlags();
1203           for (unsigned short hptdcErrorIndex = 1; hptdcErrorIndex < 16; ++hptdcErrorIndex)
1204             if (hptdcErrors.errorId(hptdcErrorIndex - 1))
1205               channelPlots_[detId].HPTDCErrorFlags->Fill(hptdcErrorIndex);
1206           if (digi.multipleHit())
1207             ++channelPlots_[detId].MHCounter;
1208 
1209           // Check dropped trailing edges
1210           if (digi.leadingEdge() != 0 || digi.trailingEdge() != 0) {
1211             ++channelPlots_[detId].HitCounter;
1212             if (digi.trailingEdge() == 0) {
1213               ++channelPlots_[detId].LeadingOnlyCounter;
1214               channelPlots_[detId].leadingEdgeCumulative_le->Fill(HPTDC_BIN_WIDTH_NS * digi.leadingEdge());
1215             } else if (digi.leadingEdge() == 0) {
1216               ++channelPlots_[detId].TrailingOnlyCounter;
1217               channelPlots_[detId].trailingEdgeCumulative_te->Fill(HPTDC_BIN_WIDTH_NS * digi.trailingEdge());
1218             } else {
1219               ++channelPlots_[detId].CompleteCounter;
1220               // channelPlots_[detId].LeadingTrailingCorrelationPerChannel->Fill(HPTDC_BIN_WIDTH_NS * digi.leadingEdge(),
1221               //                                                                 HPTDC_BIN_WIDTH_NS * digi.trailingEdge());
1222             }
1223           }
1224         }
1225       }
1226     }
1227   }
1228 
1229   // Using CTPPSDiamondRecHit
1230 
1231   for (const auto& rechits : *diamondRecHits) {
1232     CTPPSDiamondDetId detId(rechits.detId());
1233     for (const auto& rechit : rechits) {
1234       if (excludeMultipleHits_ && rechit.multipleHits() > 0)
1235         continue;
1236       if (channelPlots_.count(detId) != 0) {
1237         channelPlots_[detId].recHitTime->Fill(rechit.time());
1238         if (rechit.ootIndex() != CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING && rechit.toT() != 0) {
1239           channelPlots_[detId].leadingEdgeCumulative_both->Fill(rechit.time() + 25 * rechit.ootIndex());
1240           channelPlots_[detId].TimeOverThresholdCumulativePerChannel->Fill(rechit.toT());
1241         }
1242         ++lumiCache->hitsCounterMap[detId];
1243 
1244         if (rechit.ootIndex() != CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING &&
1245             channelPlots_[detId].activity_per_bx.count(rechit.ootIndex()) > 0)
1246           channelPlots_[detId].activity_per_bx.at(rechit.ootIndex())->Fill(event.bunchCrossing());
1247       }
1248     }
1249   }
1250 
1251   // Tomography of diamonds using pixel
1252   for (const auto& rechits : *diamondRecHits) {
1253     const CTPPSDiamondDetId detId(rechits.detId()), detId_pot(detId.rpId());
1254     const auto shift_pix = diamShifts_.at(detId_pot).withPixels;
1255     for (const auto& rechit : rechits) {
1256       if (excludeMultipleHits_ && rechit.multipleHits() > 0)
1257         continue;
1258       if (rechit.ootIndex() == CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING || rechit.toT() == 0)
1259         continue;
1260       if (!pixelTracks.isValid())
1261         continue;
1262       if (channelPlots_.count(detId) == 0)
1263         continue;
1264 
1265       for (const auto& ds : *pixelTracks) {
1266         const CTPPSPixelDetId pixId(ds.detId());
1267         if (pixId.station() != CTPPS_PIXEL_STATION_ID || pixId.rp() != CTPPS_FAR_RP_ID)
1268           continue;
1269         if (ds.size() > 1)
1270           continue;
1271         for (const auto& lt : ds) {
1272           if (lt.isValid() && pixId.arm() == detId.arm() && lt.x0() - shift_pix < 24)
1273             channelPlots_[detId].pixelTomography_far->Fill(lt.x0() - shift_pix + 25 * rechit.ootIndex(), lt.y0());
1274         }
1275       }
1276     }
1277   }
1278 }
1279 
1280 //----------------------------------------------------------------------------------------------------
1281 
1282 void CTPPSDiamondDQMSource::globalEndLuminosityBlock(const edm::LuminosityBlock& iLumi, const edm::EventSetup&) {
1283   auto lumiCache = luminosityBlockCache(iLumi.index());
1284   if (!perLSsaving_) {
1285     if (plotOnline_)
1286       for (auto& plot : potPlots_)
1287         *(plot.second.hitDistribution2d_lumisection->getTH2F()) = *(lumiCache->hitDistribution2dMap[plot.first]);
1288 
1289     for (auto& plot : channelPlots_) {
1290       auto hitsCounterPerLumisection = lumiCache->hitsCounterMap[plot.first];
1291       if (hitsCounterPerLumisection != 0) {
1292         plot.second.hit_rate->Fill((double)hitsCounterPerLumisection / SEC_PER_LUMI_SECTION);
1293       }
1294 
1295       double HundredOverHitCounter = .0;
1296       if (plot.second.HitCounter != 0)
1297         HundredOverHitCounter = 100. / plot.second.HitCounter;
1298       plot.second.HPTDCErrorFlags->setBinContent(16, HundredOverHitCounter * plot.second.MHCounter);
1299       plot.second.leadingWithoutTrailing->setBinContent(1, HundredOverHitCounter * plot.second.LeadingOnlyCounter);
1300       plot.second.leadingWithoutTrailing->setBinContent(2, HundredOverHitCounter * plot.second.TrailingOnlyCounter);
1301       plot.second.leadingWithoutTrailing->setBinContent(3, HundredOverHitCounter * plot.second.CompleteCounter);
1302     }
1303 
1304     for (auto& plot : potPlots_) {
1305       double HundredOverHitCounterPot = 0.;
1306       if (plot.second.HitCounter != 0)
1307         HundredOverHitCounterPot = 100. / plot.second.HitCounter;
1308       plot.second.leadingWithoutTrailingCumulativePot->setBinContent(
1309           1, HundredOverHitCounterPot * plot.second.LeadingOnlyCounter);
1310       plot.second.leadingWithoutTrailingCumulativePot->setBinContent(
1311           2, HundredOverHitCounterPot * plot.second.TrailingOnlyCounter);
1312       plot.second.leadingWithoutTrailingCumulativePot->setBinContent(
1313           3, HundredOverHitCounterPot * plot.second.CompleteCounter);
1314 
1315       plot.second.MHComprensive->Reset();
1316       const CTPPSDiamondDetId rpId(plot.first);
1317       for (auto& chPlot : channelPlots_) {
1318         const CTPPSDiamondDetId chId(chPlot.first);
1319         if (chId.arm() == rpId.arm() && chId.rp() == rpId.rp()) {
1320           plot.second.MHComprensive->Fill(
1321               chId.plane(), chId.channel(), chPlot.second.HPTDCErrorFlags->getBinContent(16));
1322         }
1323       }
1324     }
1325     // Efficiencies of single channels
1326     if (plotOnline_)
1327       for (auto& plot : potPlots_) {
1328         plot.second.EfficiencyOfChannelsInPot->Reset();
1329         for (auto& element : plot.second.effTriplecountingChMap) {
1330           if (plot.second.effDoublecountingChMap[element.first] > 0) {
1331             int plane = element.first / 100;
1332             int channel = element.first % 100;
1333             double counted = element.second;
1334             double total = plot.second.effDoublecountingChMap[element.first];
1335             double efficiency = counted / total;
1336             //         double error = std::sqrt( efficiency * ( 1 - efficiency ) / total );
1337 
1338             plot.second.EfficiencyOfChannelsInPot->Fill(plane, channel, 100 * efficiency);
1339           }
1340         }
1341       }
1342 
1343     // Efficeincy wrt pixels  //TODO
1344     for (auto& plot : planePlots_) {
1345       TH2F* hitHistoTmp = plot.second.EfficiencyWRTPixelsInPlane->getTH2F();
1346 
1347       const CTPPSDiamondDetId detId(plot.first), detId_pot(detId.rpId());
1348       hitHistoTmp->Divide(&(plot.second.pixelTracksMapWithDiamonds), &(potPlots_[detId_pot].pixelTracksMap));
1349     }
1350   }  //perLSsaving
1351 }
1352 
1353 //----------------------------------------------------------------------------------------------------
1354 
1355 void CTPPSDiamondDQMSource::checkEventNumber(const CTPPSDiamondDetId& detId,
1356                                              const TotemFEDInfo& optorx,
1357                                              const TotemVFATStatus& status,
1358                                              CTPPSDiamondDQMSource::PotPlots& plots,
1359                                              int& EC_difference) const {
1360   const CTPPSDiamondDetId detId_pot(detId.rpId());
1361   if (plotOnline_)
1362     plots.ECCheck->Fill((int)((optorx.lv1() & 0xFF) - ((unsigned int)status.ec() & 0xFF)) & 0xFF);
1363   if ((static_cast<int>((optorx.lv1() & 0xFF) - status.ec()) != EC_difference) &&
1364       (static_cast<uint8_t>((optorx.lv1() & 0xFF) - status.ec()) < 128))
1365     EC_difference = static_cast<int>(optorx.lv1() & 0xFF) - (static_cast<unsigned int>(status.ec()) & 0xFF);
1366   if (EC_difference != 1 && EC_difference != -500 && std::abs(EC_difference) < 127) {
1367     if (detId.channel() == HPTDC_0_CHANNEL || detId.channel() == HPTDC_1_CHANNEL)
1368       plots.HPTDCErrorFlags_2D->Fill(16, 2 * detId.plane() + (detId.channel() - HPTDC_0_CHANNEL));
1369     if (verbosity_)
1370       edm::LogProblem("CTPPSDiamondDQMSource")
1371           << "FED " << optorx.fedId() << ": ECError at EV: 0x" << std::hex << optorx.lv1() << "\t\tVFAT EC: 0x"
1372           << static_cast<unsigned int>(status.ec()) << "\twith ID: " << std::dec << detId
1373           << "\tdiff: " << EC_difference;
1374   }
1375 }
1376 
1377 //----------------------------------------------------------------------------------------------------
1378 
1379 void CTPPSDiamondDQMSource::dqmEndRun(edm::Run const&, edm::EventSetup const&) {
1380   if (plotOffline_ && !perLSsaving_)
1381     for (const auto& rpPlots : potPlots_) {
1382       auto plots = rpPlots.second;
1383       // *(plots.TOTVsLSProfile->getTProfile())=*plots.TOTVsLS->getTH2F()->ProfileX();
1384       // *(plots.trackTimeVsLSProfile->getTProfile())=*plots.trackTimeVsLS->getTH2F()->ProfileX();
1385       *(plots.trackTimeVsBXProfile->getTProfile()) = *plots.trackTimeVsBX->getTH2F()->ProfileX();
1386       // *(plots.trackTimeVsXAngleProfile->getTProfile()) = *plots.trackTimeVsXAngle->getTH2F()->ProfileX();
1387     }
1388 }
1389 
1390 //----------------------------------------------------------------------------------------------------
1391 DEFINE_FWK_MODULE(CTPPSDiamondDQMSource);