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 TotemDQM and TOTEM offline software.
0004 * Authors:
0005 *   Jan Kašpar (jan.kaspar@gmail.com)
0006 *   Rafał Leszko (rafal.leszko@gmail.com)
0007 *
0008 ****************************************************************************/
0009 
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.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/DetSetVector.h"
0021 #include "DataFormats/CTPPSDetId/interface/TotemRPDetId.h"
0022 #include "DataFormats/CTPPSDigi/interface/TotemRPDigi.h"
0023 #include "DataFormats/CTPPSDigi/interface/TotemVFATStatus.h"
0024 #include "DataFormats/CTPPSReco/interface/TotemRPCluster.h"
0025 #include "DataFormats/CTPPSReco/interface/TotemRPRecHit.h"
0026 #include "DataFormats/CTPPSReco/interface/TotemRPUVPattern.h"
0027 #include "DataFormats/CTPPSReco/interface/TotemRPLocalTrack.h"
0028 
0029 #include "Geometry/Records/interface/VeryForwardRealGeometryRecord.h"
0030 #include "Geometry/VeryForwardGeometryBuilder/interface/CTPPSGeometry.h"
0031 #include "Geometry/VeryForwardRPTopology/interface/RPTopology.h"
0032 
0033 #include <string>
0034 
0035 //----------------------------------------------------------------------------------------------------
0036 
0037 class TotemRPDQMSource : public DQMEDAnalyzer {
0038 public:
0039   TotemRPDQMSource(const edm::ParameterSet &ps);
0040   ~TotemRPDQMSource() override;
0041 
0042 protected:
0043   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0044   void analyze(edm::Event const &e, edm::EventSetup const &eSetup) override;
0045 
0046 private:
0047   unsigned int verbosity;
0048 
0049   edm::EDGetTokenT<edm::DetSetVector<TotemVFATStatus>> tokenStatus;
0050   edm::EDGetTokenT<edm::DetSetVector<TotemRPDigi>> tokenDigi;
0051   edm::EDGetTokenT<edm::DetSetVector<TotemRPCluster>> tokenCluster;
0052   edm::EDGetTokenT<edm::DetSetVector<TotemRPRecHit>> tokenRecHit;
0053   edm::EDGetTokenT<edm::DetSetVector<TotemRPUVPattern>> tokenUVPattern;
0054   edm::EDGetTokenT<edm::DetSetVector<TotemRPLocalTrack>> tokenLocalTrack;
0055 
0056   edm::ESGetToken<CTPPSGeometry, VeryForwardRealGeometryRecord> geometryToken_;
0057 
0058   /// plots related to one RP
0059   struct PotPlots {
0060     MonitorElement *vfat_problem = nullptr, *vfat_missing = nullptr, *vfat_ec_bc_error = nullptr,
0061                    *vfat_corruption = nullptr;
0062 
0063     MonitorElement *activity = nullptr, *activity_u = nullptr, *activity_v = nullptr;
0064     MonitorElement *activity_per_bx = nullptr, *activity_per_bx_short = nullptr;
0065     MonitorElement *hit_plane_hist = nullptr;
0066     MonitorElement *patterns_u = nullptr, *patterns_v = nullptr;
0067     MonitorElement *h_planes_fit_u = nullptr, *h_planes_fit_v = nullptr;
0068     MonitorElement *event_category = nullptr;
0069     MonitorElement *trackHitsCumulativeHist = nullptr;
0070     MonitorElement *track_u_profile = nullptr, *track_v_profile = nullptr;
0071     MonitorElement *triggerSectorUVCorrelation_all = nullptr, *triggerSectorUVCorrelation_mult2 = nullptr,
0072                    *triggerSectorUVCorrelation_track = nullptr;
0073 
0074     PotPlots() {}
0075     PotPlots(DQMStore::IBooker &ibooker, unsigned int id);
0076   };
0077 
0078   std::map<unsigned int, PotPlots> potPlots;
0079 
0080   /// plots related to one RP plane
0081   struct PlanePlots {
0082     MonitorElement *digi_profile_cumulative = nullptr;
0083     MonitorElement *cluster_profile_cumulative = nullptr;
0084     MonitorElement *hit_multiplicity = nullptr;
0085     MonitorElement *cluster_size = nullptr;
0086     MonitorElement *efficiency_num = nullptr, *efficiency_den = nullptr;
0087 
0088     PlanePlots() {}
0089     PlanePlots(DQMStore::IBooker &ibooker, unsigned int id);
0090   };
0091 
0092   std::map<unsigned int, PlanePlots> planePlots;
0093 };
0094 
0095 //----------------------------------------------------------------------------------------------------
0096 //----------------------------------------------------------------------------------------------------
0097 
0098 using namespace std;
0099 using namespace edm;
0100 
0101 //----------------------------------------------------------------------------------------------------
0102 
0103 TotemRPDQMSource::PotPlots::PotPlots(DQMStore::IBooker &ibooker, unsigned int id) {
0104   string path;
0105   TotemRPDetId(id).rpName(path, TotemRPDetId::nPath);
0106   ibooker.setCurrentFolder(path);
0107 
0108   string title;
0109   TotemRPDetId(id).rpName(title, TotemRPDetId::nFull);
0110 
0111   vfat_problem = ibooker.book2D("vfats with any problem", title + ";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
0112   vfat_missing = ibooker.book2D("vfats missing", title + ";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
0113   vfat_ec_bc_error =
0114       ibooker.book2D("vfats with EC or BC error", title + ";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
0115   vfat_corruption =
0116       ibooker.book2D("vfats with data corruption", title + ";plane;vfat index", 10, -0.5, 9.5, 4, -0.5, 3.5);
0117 
0118   activity = ibooker.book1D("active planes", title + ";number of active planes", 11, -0.5, 10.5);
0119   activity_u = ibooker.book1D("active planes U", title + ";number of active U planes", 11, -0.5, 10.5);
0120   activity_v = ibooker.book1D("active planes V", title + ";number of active V planes", 11, -0.5, 10.5);
0121 
0122   activity_per_bx = ibooker.book1D("activity per BX", title + ";Event.BX", 4002, -1.5, 4000. + 0.5);
0123   activity_per_bx_short = ibooker.book1D("activity per BX (short)", title + ";Event.BX", 102, -1.5, 100. + 0.5);
0124 
0125   hit_plane_hist =
0126       ibooker.book2D("activity in planes (2D)", title + ";plane number;strip number", 10, -0.5, 9.5, 32, -0.5, 511.5);
0127 
0128   patterns_u = ibooker.book1D("recognized patterns U", title + ";number of recognized U patterns", 11, -0.5, 10.5);
0129   patterns_v = ibooker.book1D("recognized patterns V", title + ";number of recognized V patterns", 11, -0.5, 10.5);
0130 
0131   h_planes_fit_u =
0132       ibooker.book1D("planes contributing to fit U", title + ";number of planes contributing to U fit", 6, -0.5, 5.5);
0133   h_planes_fit_v =
0134       ibooker.book1D("planes contributing to fit V", title + ";number of planes contributing to V fit", 6, -0.5, 5.5);
0135 
0136   event_category = ibooker.book1D("event category", title + ";event category", 5, -0.5, 4.5);
0137   TH1F *event_category_h = event_category->getTH1F();
0138   event_category_h->GetXaxis()->SetBinLabel(1, "empty");
0139   event_category_h->GetXaxis()->SetBinLabel(2, "insufficient");
0140   event_category_h->GetXaxis()->SetBinLabel(3, "single-track");
0141   event_category_h->GetXaxis()->SetBinLabel(4, "multi-track");
0142   event_category_h->GetXaxis()->SetBinLabel(5, "shower");
0143 
0144   trackHitsCumulativeHist =
0145       ibooker.book2D("track XY profile", title + ";x   (mm);y   (mm)", 100, -18., +18., 100, -18., +18.);
0146 
0147   track_u_profile = ibooker.book1D("track profile U", title + "; U   (mm)", 512, -256 * 66E-3, +256 * 66E-3);
0148   track_v_profile = ibooker.book1D("track profile V", title + "; V   (mm)", 512, -256 * 66E-3, +256 * 66E-3);
0149 
0150   triggerSectorUVCorrelation_all = ibooker.book2D(
0151       "trigger sector UV correlation (no cond)", title + ";V sector;U sector", 16, -0.5, 15.5, 16, -0.5, 15.5);
0152   triggerSectorUVCorrelation_mult2 = ibooker.book2D(
0153       "trigger sector UV correlation (max mult 2)", title + ";V sector;U sector", 16, -0.5, 15.5, 16, -0.5, 15.5);
0154   triggerSectorUVCorrelation_track = ibooker.book2D(
0155       "trigger sector UV correlation (track)", title + ";V sector;U sector", 16, -0.5, 15.5, 16, -0.5, 15.5);
0156 }
0157 
0158 //----------------------------------------------------------------------------------------------------
0159 
0160 TotemRPDQMSource::PlanePlots::PlanePlots(DQMStore::IBooker &ibooker, unsigned int id) {
0161   string path;
0162   TotemRPDetId(id).planeName(path, TotemRPDetId::nPath);
0163   ibooker.setCurrentFolder(path);
0164 
0165   string title;
0166   TotemRPDetId(id).planeName(title, TotemRPDetId::nFull);
0167 
0168   digi_profile_cumulative = ibooker.book1D("digi profile", title + ";strip number", 512, -0.5, 511.5);
0169   cluster_profile_cumulative = ibooker.book1D("cluster profile", title + ";cluster center", 1024, -0.25, 511.75);
0170   hit_multiplicity = ibooker.book1D("hit multiplicity", title + ";hits/detector/event", 6, -0.5, 5.5);
0171   cluster_size = ibooker.book1D("cluster size", title + ";hits per cluster", 5, 0.5, 5.5);
0172 
0173   efficiency_num = ibooker.book1D("efficiency num", title + ";track position   (mm)", 30, -15., 0.);
0174   efficiency_den = ibooker.book1D("efficiency den", title + ";track position   (mm)", 30, -15., 0.);
0175 }
0176 
0177 //----------------------------------------------------------------------------------------------------
0178 //----------------------------------------------------------------------------------------------------
0179 
0180 TotemRPDQMSource::TotemRPDQMSource(const edm::ParameterSet &ps)
0181     : verbosity(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
0182       geometryToken_(esConsumes())
0183 
0184 {
0185   tokenStatus = consumes<DetSetVector<TotemVFATStatus>>(ps.getUntrackedParameter<edm::InputTag>("tagStatus"));
0186 
0187   tokenDigi = consumes<DetSetVector<TotemRPDigi>>(ps.getUntrackedParameter<edm::InputTag>("tagDigi"));
0188   tokenCluster = consumes<edm::DetSetVector<TotemRPCluster>>(ps.getUntrackedParameter<edm::InputTag>("tagCluster"));
0189   tokenRecHit = consumes<edm::DetSetVector<TotemRPRecHit>>(ps.getUntrackedParameter<edm::InputTag>("tagRecHit"));
0190   tokenUVPattern = consumes<DetSetVector<TotemRPUVPattern>>(ps.getUntrackedParameter<edm::InputTag>("tagUVPattern"));
0191   tokenLocalTrack = consumes<DetSetVector<TotemRPLocalTrack>>(ps.getUntrackedParameter<edm::InputTag>("tagLocalTrack"));
0192 }
0193 
0194 //----------------------------------------------------------------------------------------------------
0195 
0196 TotemRPDQMSource::~TotemRPDQMSource() {}
0197 
0198 //----------------------------------------------------------------------------------------------------
0199 
0200 void TotemRPDQMSource::bookHistograms(DQMStore::IBooker &ibooker, edm::Run const &, edm::EventSetup const &) {
0201   ibooker.cd();
0202   ibooker.setCurrentFolder("CTPPS");
0203 
0204   // loop over arms
0205   for (unsigned int arm : {0, 1}) {
0206     // loop over RPs
0207     for (unsigned int st_rp : {2, 3, 4, 5, 24, 25}) {
0208       const unsigned int st = st_rp / 10;
0209       const unsigned int rp = st_rp % 10;
0210 
0211       TotemRPDetId rpId(arm, st, rp);
0212       potPlots[rpId] = PotPlots(ibooker, rpId);
0213 
0214       // loop over planes
0215       for (unsigned int pl = 0; pl < 10; ++pl) {
0216         TotemRPDetId plId(arm, st, rp, pl);
0217         planePlots[plId] = PlanePlots(ibooker, plId);
0218       }
0219     }
0220   }
0221 }
0222 
0223 //----------------------------------------------------------------------------------------------------
0224 
0225 void TotemRPDQMSource::analyze(edm::Event const &event, edm::EventSetup const &eventSetup) {
0226   // get event setup data
0227   auto const &geometry = eventSetup.getData(geometryToken_);
0228 
0229   // get event data
0230   Handle<DetSetVector<TotemVFATStatus>> status;
0231   event.getByToken(tokenStatus, status);
0232 
0233   Handle<DetSetVector<TotemRPDigi>> digi;
0234   event.getByToken(tokenDigi, digi);
0235 
0236   Handle<DetSetVector<TotemRPCluster>> digCluster;
0237   event.getByToken(tokenCluster, digCluster);
0238 
0239   Handle<DetSetVector<TotemRPRecHit>> hits;
0240   event.getByToken(tokenRecHit, hits);
0241 
0242   Handle<DetSetVector<TotemRPUVPattern>> patterns;
0243   event.getByToken(tokenUVPattern, patterns);
0244 
0245   Handle<DetSetVector<TotemRPLocalTrack>> tracks;
0246   event.getByToken(tokenLocalTrack, tracks);
0247 
0248   // check validity
0249   bool valid = true;
0250   valid &= status.isValid();
0251   valid &= digi.isValid();
0252   valid &= digCluster.isValid();
0253   valid &= hits.isValid();
0254   valid &= patterns.isValid();
0255   valid &= tracks.isValid();
0256 
0257   if (!valid) {
0258     if (verbosity) {
0259       LogProblem("TotemRPDQMSource")
0260           << "ERROR in TotemDQMModuleRP::analyze > some of the required inputs are not valid. Skipping this event.\n"
0261           << "    status.isValid = " << status.isValid() << "\n"
0262           << "    digi.isValid = " << digi.isValid() << "\n"
0263           << "    digCluster.isValid = " << digCluster.isValid() << "\n"
0264           << "    hits.isValid = " << hits.isValid() << "\n"
0265           << "    patterns.isValid = " << patterns.isValid() << "\n"
0266           << "    tracks.isValid = " << tracks.isValid();
0267     }
0268 
0269     return;
0270   }
0271 
0272   //------------------------------
0273   // Status Plots
0274 
0275   for (auto &ds : *status) {
0276     TotemRPDetId detId(ds.detId());
0277     unsigned int plNum = detId.plane();
0278     CTPPSDetId rpId = detId.rpId();
0279 
0280     auto it = potPlots.find(rpId);
0281     if (it == potPlots.end())
0282       continue;
0283     auto &plots = it->second;
0284 
0285     for (auto &s : ds) {
0286       if (s.isMissing()) {
0287         plots.vfat_problem->Fill(plNum, s.chipPosition());
0288         plots.vfat_missing->Fill(plNum, s.chipPosition());
0289       }
0290 
0291       if (s.isECProgressError() || s.isBCProgressError()) {
0292         plots.vfat_problem->Fill(plNum, s.chipPosition());
0293         plots.vfat_ec_bc_error->Fill(plNum, s.chipPosition());
0294       }
0295 
0296       if (s.isIDMismatch() || s.isFootprintError() || s.isCRCError()) {
0297         plots.vfat_problem->Fill(plNum, s.chipPosition());
0298         plots.vfat_corruption->Fill(plNum, s.chipPosition());
0299       }
0300     }
0301   }
0302 
0303   //------------------------------
0304   // Plane Plots
0305 
0306   // digi profile cumulative
0307   for (DetSetVector<TotemRPDigi>::const_iterator it = digi->begin(); it != digi->end(); ++it) {
0308     TotemRPDetId detId(it->detId());
0309 
0310     auto plIt = planePlots.find(detId);
0311     if (plIt == planePlots.end())
0312       continue;
0313     auto &plots = plIt->second;
0314 
0315     for (DetSet<TotemRPDigi>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
0316       plots.digi_profile_cumulative->Fill(dit->stripNumber());
0317   }
0318 
0319   // cluster plots
0320   for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++) {
0321     TotemRPDetId detId(it->detId());
0322 
0323     auto plIt = planePlots.find(detId);
0324     if (plIt == planePlots.end())
0325       continue;
0326     auto &plots = plIt->second;
0327 
0328     // hit multiplicity
0329     plots.hit_multiplicity->Fill(it->size());
0330 
0331     for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit) {
0332       // profile cumulative
0333       plots.cluster_profile_cumulative->Fill(dit->centerStripPosition());
0334 
0335       // cluster size
0336       plots.cluster_size->Fill(dit->numberOfStrips());
0337     }
0338   }
0339 
0340   // plane efficiency plots
0341   for (auto &ds : *tracks) {
0342     CTPPSDetId rpId(ds.detId());
0343 
0344     for (auto &ft : ds) {
0345       if (!ft.isValid())
0346         continue;
0347 
0348       double rp_z = geometry.rpTranslation(rpId).z();
0349 
0350       for (unsigned int plNum = 0; plNum < 10; ++plNum) {
0351         TotemRPDetId plId = rpId;
0352         plId.setPlane(plNum);
0353 
0354         auto plIt = planePlots.find(plId);
0355         if (plIt == planePlots.end())
0356           continue;
0357         auto &plots = plIt->second;
0358 
0359         double ft_z = ft.z0();
0360         double ft_x = ft.x0() + ft.tx() * (ft_z - rp_z);
0361         double ft_y = ft.y0() + ft.ty() * (ft_z - rp_z);
0362 
0363         double ft_v = geometry.globalToLocal(plId, CTPPSGeometry::Vector(ft_x, ft_y, ft_z)).y();
0364 
0365         bool hasMatchingHit = false;
0366         const auto &hit_ds_it = hits->find(plId);
0367         if (hit_ds_it != hits->end()) {
0368           for (const auto &h : *hit_ds_it) {
0369             bool match = (fabs(ft_v - h.position()) < 2. * 0.066);
0370             if (match) {
0371               hasMatchingHit = true;
0372               break;
0373             }
0374           }
0375         }
0376 
0377         plots.efficiency_den->Fill(ft_v);
0378         if (hasMatchingHit)
0379           plots.efficiency_num->Fill(ft_v);
0380       }
0381     }
0382   }
0383 
0384   //------------------------------
0385   // Roman Pots Plots
0386 
0387   // determine active planes (from RecHits and VFATStatus)
0388   map<unsigned int, set<unsigned int>> planes;
0389   map<unsigned int, set<unsigned int>> planes_u;
0390   map<unsigned int, set<unsigned int>> planes_v;
0391   for (const auto &ds : *hits) {
0392     if (ds.empty())
0393       continue;
0394 
0395     TotemRPDetId detId(ds.detId());
0396     unsigned int planeNum = detId.plane();
0397     CTPPSDetId rpId = detId.rpId();
0398 
0399     planes[rpId].insert(planeNum);
0400     if (detId.isStripsCoordinateUDirection())
0401       planes_u[rpId].insert(planeNum);
0402     else
0403       planes_v[rpId].insert(planeNum);
0404   }
0405 
0406   for (const auto &ds : *status) {
0407     bool activity = false;
0408     for (const auto &s : ds) {
0409       if (s.isNumberOfClustersSpecified() && s.numberOfClusters() > 0) {
0410         activity = true;
0411         break;
0412       }
0413     }
0414 
0415     if (!activity)
0416       continue;
0417 
0418     TotemRPDetId detId(ds.detId());
0419     unsigned int planeNum = detId.plane();
0420     CTPPSDetId rpId = detId.rpId();
0421 
0422     planes[rpId].insert(planeNum);
0423     if (detId.isStripsCoordinateUDirection())
0424       planes_u[rpId].insert(planeNum);
0425     else
0426       planes_v[rpId].insert(planeNum);
0427   }
0428 
0429   // plane activity histogram
0430   for (std::map<unsigned int, PotPlots>::iterator it = potPlots.begin(); it != potPlots.end(); it++) {
0431     it->second.activity->Fill(planes[it->first].size());
0432     it->second.activity_u->Fill(planes_u[it->first].size());
0433     it->second.activity_v->Fill(planes_v[it->first].size());
0434 
0435     if (planes[it->first].size() >= 6) {
0436       it->second.activity_per_bx->Fill(event.bunchCrossing());
0437       it->second.activity_per_bx_short->Fill(event.bunchCrossing());
0438     }
0439   }
0440 
0441   for (DetSetVector<TotemRPCluster>::const_iterator it = digCluster->begin(); it != digCluster->end(); it++) {
0442     TotemRPDetId detId(it->detId());
0443     unsigned int planeNum = detId.plane();
0444     CTPPSDetId rpId = detId.rpId();
0445 
0446     auto plIt = potPlots.find(rpId);
0447     if (plIt == potPlots.end())
0448       continue;
0449     auto &plots = plIt->second;
0450 
0451     for (DetSet<TotemRPCluster>::const_iterator dit = it->begin(); dit != it->end(); ++dit)
0452       plots.hit_plane_hist->Fill(planeNum, dit->centerStripPosition());
0453   }
0454 
0455   // recognized pattern histograms
0456   for (auto &ds : *patterns) {
0457     CTPPSDetId rpId(ds.detId());
0458 
0459     auto plIt = potPlots.find(rpId);
0460     if (plIt == potPlots.end())
0461       continue;
0462     auto &plots = plIt->second;
0463 
0464     // count U and V patterns
0465     unsigned int u = 0, v = 0;
0466     for (auto &p : ds) {
0467       if (!p.fittable())
0468         continue;
0469 
0470       if (p.projection() == TotemRPUVPattern::projU)
0471         u++;
0472 
0473       if (p.projection() == TotemRPUVPattern::projV)
0474         v++;
0475     }
0476 
0477     plots.patterns_u->Fill(u);
0478     plots.patterns_v->Fill(v);
0479   }
0480 
0481   // event-category histogram
0482   for (auto &it : potPlots) {
0483     TotemRPDetId rpId(it.first);
0484     auto &pp = it.second;
0485 
0486     // process hit data for this plot
0487     unsigned int pl_u = planes_u[rpId].size();
0488     unsigned int pl_v = planes_v[rpId].size();
0489 
0490     // process pattern data for this pot
0491     const auto &rp_pat_it = patterns->find(rpId);
0492 
0493     unsigned int pat_u = 0, pat_v = 0;
0494     if (rp_pat_it != patterns->end()) {
0495       for (auto &p : *rp_pat_it) {
0496         if (!p.fittable())
0497           continue;
0498 
0499         if (p.projection() == TotemRPUVPattern::projU)
0500           pat_u++;
0501 
0502         if (p.projection() == TotemRPUVPattern::projV)
0503           pat_v++;
0504       }
0505     }
0506 
0507     // determine category
0508     signed int category = -1;
0509 
0510     if (pl_u == 0 && pl_v == 0)
0511       category = 0;  // empty
0512 
0513     if (category == -1 && pat_u + pat_v <= 1) {
0514       if (pl_u + pl_v < 6)
0515         category = 1;  // insuff
0516       else
0517         category = 4;  // shower
0518     }
0519 
0520     if (pat_u == 1 && pat_v == 1)
0521       category = 2;  // 1-track
0522 
0523     if (category == -1)
0524       category = 3;  // multi-track
0525 
0526     pp.event_category->Fill(category);
0527   }
0528 
0529   // RP track-fit plots
0530   set<unsigned int> rps_with_tracks;
0531 
0532   for (auto &ds : *tracks) {
0533     CTPPSDetId rpId(ds.detId());
0534 
0535     rps_with_tracks.insert(rpId);
0536 
0537     auto plIt = potPlots.find(rpId);
0538     if (plIt == potPlots.end())
0539       continue;
0540     auto &plots = plIt->second;
0541 
0542     for (auto &ft : ds) {
0543       if (!ft.isValid())
0544         continue;
0545 
0546       // number of planes contributing to (valid) fits
0547       unsigned int n_pl_in_fit_u = 0, n_pl_in_fit_v = 0;
0548       for (auto &hds : ft.hits()) {
0549         TotemRPDetId plId(hds.detId());
0550         bool uProj = plId.isStripsCoordinateUDirection();
0551 
0552         for (auto &h : hds) {
0553           h.position();  // just to keep compiler silent
0554           if (uProj)
0555             n_pl_in_fit_u++;
0556           else
0557             n_pl_in_fit_v++;
0558         }
0559       }
0560 
0561       plots.h_planes_fit_u->Fill(n_pl_in_fit_u);
0562       plots.h_planes_fit_v->Fill(n_pl_in_fit_v);
0563 
0564       // mean position of U and V planes
0565       TotemRPDetId plId_V(rpId);
0566       plId_V.setPlane(0);
0567       TotemRPDetId plId_U(rpId);
0568       plId_U.setPlane(1);
0569 
0570       double rp_x = (geometry.sensor(plId_V)->translation().x() + geometry.sensor(plId_U)->translation().x()) / 2.;
0571       double rp_y = (geometry.sensor(plId_V)->translation().y() + geometry.sensor(plId_U)->translation().y()) / 2.;
0572 
0573       // mean read-out direction of U and V planes
0574       const auto &rod_U = geometry.localToGlobalDirection(plId_U, CTPPSGeometry::Vector(0., 1., 0.));
0575       const auto &rod_V = geometry.localToGlobalDirection(plId_V, CTPPSGeometry::Vector(0., 1., 0.));
0576 
0577       double x = ft.x0() - rp_x;
0578       double y = ft.y0() - rp_y;
0579 
0580       plots.trackHitsCumulativeHist->Fill(x, y);
0581 
0582       double U = x * rod_U.x() + y * rod_U.y();
0583       double V = x * rod_V.x() + y * rod_V.y();
0584 
0585       plots.track_u_profile->Fill(U);
0586       plots.track_v_profile->Fill(V);
0587     }
0588   }
0589 
0590   // restore trigger-sector map from digis
0591   map<unsigned int, map<unsigned int, map<unsigned int, unsigned int>>>
0592       triggerSectorMap;  // [rpId, U/V flag, sector] --> number of planes
0593   for (const auto &dp : *digi) {
0594     TotemRPDetId plId(dp.detId());
0595     CTPPSDetId rpId = plId.rpId();
0596     unsigned int uvFlag = (plId.isStripsCoordinateUDirection()) ? 0 : 1;
0597 
0598     set<unsigned int> sectors;
0599     for (const auto &d : dp) {
0600       unsigned int sector = d.stripNumber() / 32;
0601       sectors.insert(sector);
0602     }
0603 
0604     for (const auto &sector : sectors)
0605       triggerSectorMap[rpId][uvFlag][sector]++;
0606   }
0607 
0608   for (auto &rpp : triggerSectorMap) {
0609     const unsigned int rpId = rpp.first;
0610 
0611     // trigger sector is counted as active if at least 3 planes report activity
0612 
0613     set<unsigned int> triggerSectorsU;
0614     for (const auto sp : rpp.second[0]) {
0615       if (sp.second >= 3)
0616         triggerSectorsU.insert(sp.first);
0617     }
0618 
0619     set<unsigned int> triggerSectorsV;
0620     for (const auto sp : rpp.second[1]) {
0621       if (sp.second >= 3)
0622         triggerSectorsV.insert(sp.first);
0623     }
0624 
0625     auto plIt = potPlots.find(rpId);
0626     if (plIt == potPlots.end())
0627       continue;
0628     auto &plots = plIt->second;
0629 
0630     const bool high_mult = (triggerSectorsU.size() > 2 && triggerSectorsV.size() > 2);
0631 
0632     const bool has_track = (rps_with_tracks.find(rpId) != rps_with_tracks.end());
0633 
0634     for (const auto &secU : triggerSectorsU) {
0635       for (const auto &secV : triggerSectorsV) {
0636         plots.triggerSectorUVCorrelation_all->Fill(secV, secU);
0637 
0638         if (!high_mult)
0639           plots.triggerSectorUVCorrelation_mult2->Fill(secV, secU);
0640 
0641         if (has_track)
0642           plots.triggerSectorUVCorrelation_track->Fill(secV, secU);
0643       }
0644     }
0645   }
0646 }
0647 
0648 //----------------------------------------------------------------------------------------------------
0649 
0650 DEFINE_FWK_MODULE(TotemRPDQMSource);