File indexing completed on 2023-03-17 10:54:04
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/MakerMacros.h"
0014 #include "FWCore/Framework/interface/Run.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017
0018 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0019 #include "DQMServices/Core/interface/DQMStore.h"
0020
0021 #include "DataFormats/CTPPSDigi/interface/TotemFEDInfo.h"
0022 #include "DataFormats/CTPPSDigi/interface/TotemVFATStatus.h"
0023 #include "DataFormats/CTPPSReco/interface/TotemRPLocalTrack.h"
0024 #include "DataFormats/Common/interface/DetSetVector.h"
0025 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0026 #include "DataFormats/Provenance/interface/EventRange.h"
0027
0028 #include "DataFormats/CTPPSDetId/interface/CTPPSDiamondDetId.h"
0029 #include "DataFormats/CTPPSDigi/interface/TotemTimingDigi.h"
0030 #include "DataFormats/CTPPSReco/interface/TotemTimingRecHit.h"
0031
0032 #include "DataFormats/CTPPSDetId/interface/TotemRPDetId.h"
0033 #include "DataFormats/CTPPSReco/interface/TotemTimingLocalTrack.h"
0034 #include "Geometry/Records/interface/VeryForwardRealGeometryRecord.h"
0035 #include "Geometry/VeryForwardGeometryBuilder/interface/CTPPSGeometry.h"
0036
0037 #include <string>
0038
0039
0040
0041 namespace totemds {
0042 struct Cache {
0043 std::unordered_map<unsigned int, std::unique_ptr<TH2F>> hitDistribution2dMap;
0044
0045 std::unordered_map<unsigned int, unsigned long> hitsCounterMap;
0046 };
0047 }
0048
0049 class DiamondSampicDQMSource : public DQMOneEDAnalyzer<edm::LuminosityBlockCache<totemds::Cache>> {
0050 public:
0051 DiamondSampicDQMSource(const edm::ParameterSet &);
0052 ~DiamondSampicDQMSource() override;
0053
0054 protected:
0055 void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override;
0056 void bookHistograms(DQMStore::IBooker &, const edm::Run &, const edm::EventSetup &) override;
0057 void analyze(const edm::Event &, const edm::EventSetup &) override;
0058 std::shared_ptr<totemds::Cache> globalBeginLuminosityBlock(const edm::LuminosityBlock &,
0059 const edm::EventSetup &) const override;
0060 void globalEndLuminosityBlock(const edm::LuminosityBlock &, const edm::EventSetup &) override;
0061
0062 private:
0063
0064 static const double SEC_PER_LUMI_SECTION;
0065
0066
0067 static const double LHC_CLOCK_PERIOD_NS;
0068 static const double DQM_FRACTION_OF_EVENTS;
0069
0070 static const double HIT_RATE_FACTOR;
0071 static const double DISPLAY_RESOLUTION_FOR_HITS_MM;
0072
0073
0074 static const double INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0075 static const double SAMPIC_ADC_V;
0076
0077 edm::EDGetTokenT<edm::DetSetVector<TotemRPLocalTrack>> tokenLocalTrack_;
0078 edm::EDGetTokenT<edm::DetSetVector<TotemTimingDigi>> tokenDigi_;
0079 edm::EDGetTokenT<edm::DetSetVector<TotemTimingRecHit>> tokenRecHit_;
0080 edm::EDGetTokenT<edm::DetSetVector<TotemTimingLocalTrack>> tokenTrack_;
0081 edm::EDGetTokenT<std::vector<TotemFEDInfo>> tokenFEDInfo_;
0082
0083 edm::ESGetToken<CTPPSGeometry, VeryForwardRealGeometryRecord> ctppsGeometryRunToken_;
0084
0085 unsigned int samplesForNoise_;
0086 unsigned int verbosity_;
0087 bool plotOnline_;
0088 bool perLSsaving_;
0089 unsigned int trackCorrelationThreshold_;
0090 edm::TimeValue_t timeOfPreviousEvent_;
0091
0092 std::unordered_map<unsigned int, double> horizontalShiftOfDiamond_;
0093
0094
0095 struct GlobalPlots {
0096 MonitorElement *digiSentPercentage = nullptr;
0097
0098 GlobalPlots() {}
0099 GlobalPlots(DQMStore::IBooker &ibooker);
0100 };
0101
0102 GlobalPlots globalPlot_;
0103
0104 struct SectorPlots {
0105
0106 MonitorElement *trackCorrelation = nullptr;
0107 MonitorElement *trackCorrelationLowMultiplicity = nullptr;
0108 MonitorElement *digiSentPercentage = nullptr;
0109 SectorPlots(){};
0110 SectorPlots(DQMStore::IBooker &ibooker, unsigned int id, bool plotOnline);
0111 };
0112 std::unordered_map<unsigned int, SectorPlots> sectorPlots_;
0113
0114
0115 struct PotPlots {
0116
0117 MonitorElement *activityPerBX = nullptr;
0118 MonitorElement *digiDistribution = nullptr;
0119 MonitorElement *dataSamplesRaw = nullptr;
0120 MonitorElement *baseline = nullptr;
0121 MonitorElement *noiseRMS = nullptr;
0122
0123 MonitorElement *digiSent = nullptr;
0124 MonitorElement *digiAll = nullptr;
0125 MonitorElement *digiSentPercentage = nullptr;
0126
0127
0128 MonitorElement *hitDistribution2d = nullptr;
0129 MonitorElement *hitDistribution2dWithTime = nullptr;
0130 MonitorElement *hitDistribution2d_lumisection = nullptr;
0131
0132 MonitorElement *recHitTime = nullptr;
0133 MonitorElement *amplitude = nullptr;
0134 MonitorElement *baselineRMS = nullptr;
0135 MonitorElement *triggerCellTime = nullptr;
0136 MonitorElement *meanAmplitude = nullptr;
0137 MonitorElement *cellOfMax = nullptr;
0138
0139 MonitorElement *hitRate = nullptr;
0140
0141 MonitorElement *planesWithDigis = nullptr;
0142 MonitorElement *planesWithTime = nullptr;
0143
0144 MonitorElement *trackDistribution = nullptr;
0145
0146 std::set<unsigned int> planesWithDigisSet;
0147 std::set<unsigned int> planesWithTimeSet;
0148
0149 PotPlots(){};
0150 PotPlots(DQMStore::IBooker &ibooker, unsigned int id, bool plotOnline);
0151 };
0152
0153 std::unordered_map<unsigned int, PotPlots> potPlots_;
0154
0155
0156 struct PlanePlots {
0157 MonitorElement *digiDistribution = nullptr;
0158
0159 MonitorElement *hitProfile = nullptr;
0160 MonitorElement *hitMultiplicity = nullptr;
0161 MonitorElement *hitMultiplicityWithTime = nullptr;
0162
0163 PlanePlots() {}
0164 PlanePlots(DQMStore::IBooker &ibooker, unsigned int id);
0165 };
0166
0167 std::unordered_map<unsigned int, PlanePlots> planePlots_;
0168
0169
0170 struct ChannelPlots {
0171
0172 MonitorElement *activityPerBX = nullptr;
0173 MonitorElement *dataSamplesRaw = nullptr;
0174 MonitorElement *cellOfMax = nullptr;
0175 MonitorElement *maxTimeAfterTrigger = nullptr;
0176
0177
0178 MonitorElement *triggerCellTime = nullptr;
0179 MonitorElement *recHitTime = nullptr;
0180 MonitorElement *amplitude = nullptr;
0181 MonitorElement *noiseSamples = nullptr;
0182
0183
0184 MonitorElement *hitRate = nullptr;
0185
0186 ChannelPlots() {}
0187 ChannelPlots(DQMStore::IBooker &ibooker, unsigned int id);
0188 };
0189
0190 std::unordered_map<unsigned int, ChannelPlots> channelPlots_;
0191 static std::string changePathToSampic(std::string path);
0192 };
0193
0194 std::string DiamondSampicDQMSource::changePathToSampic(std::string path) {
0195 std::string toReplace = "TimingDiamond";
0196 path = path.substr(path.find(toReplace) + toReplace.length());
0197 path = "CTPPS/DiamondSampic/" + path;
0198 return path;
0199 }
0200
0201
0202
0203 const double DiamondSampicDQMSource::SEC_PER_LUMI_SECTION = 23.31;
0204 const double DiamondSampicDQMSource::LHC_CLOCK_PERIOD_NS = 24.95;
0205 const double DiamondSampicDQMSource::DQM_FRACTION_OF_EVENTS = 1.;
0206 const double DiamondSampicDQMSource::HIT_RATE_FACTOR = DQM_FRACTION_OF_EVENTS / SEC_PER_LUMI_SECTION;
0207 const double DiamondSampicDQMSource::DISPLAY_RESOLUTION_FOR_HITS_MM = 0.05;
0208 const double DiamondSampicDQMSource::INV_DISPLAY_RESOLUTION_FOR_HITS_MM = 1. / DISPLAY_RESOLUTION_FOR_HITS_MM;
0209 const double DiamondSampicDQMSource::SAMPIC_ADC_V = 1. / 256;
0210
0211
0212 DiamondSampicDQMSource::GlobalPlots::GlobalPlots(DQMStore::IBooker &ibooker) {
0213 ibooker.setCurrentFolder("CTPPS/DiamondSampic");
0214 digiSentPercentage = ibooker.book2D(
0215 "sent digis percentage", "sent digis percentage (sampic);board + 0.5 sampic;channel", 14, -0.5, 6.5, 16, 0, 16);
0216 }
0217
0218
0219
0220 DiamondSampicDQMSource::SectorPlots::SectorPlots(DQMStore::IBooker &ibooker, unsigned int id, bool plotOnline) {
0221 std::string path, title;
0222 CTPPSDiamondDetId(id).armName(path, CTPPSDiamondDetId::nPath);
0223 ibooker.setCurrentFolder(DiamondSampicDQMSource::changePathToSampic(path));
0224
0225 CTPPSDiamondDetId(id).armName(title, CTPPSDiamondDetId::nFull);
0226
0227 trackCorrelation = ibooker.book2D("tracks correlation near-far",
0228 title + " tracks correlation near-far;x (mm);x (mm)",
0229 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0230 -1,
0231 18,
0232 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0233 -1,
0234 18);
0235 trackCorrelationLowMultiplicity =
0236 ibooker.book2D("tracks correlation with low multiplicity near-far",
0237 title + " tracks correlation with low multiplicity near-far;x (mm);x (mm)",
0238 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0239 -1,
0240 18,
0241 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0242 -1,
0243 18);
0244 if (plotOnline)
0245 digiSentPercentage = ibooker.book2D("sent digis percentage",
0246 title + " sent digis percentage (sampic);board + 0.5 sampic;channel",
0247 14,
0248 -0.5,
0249 6.5,
0250 16,
0251 0,
0252 16);
0253 }
0254
0255
0256
0257 DiamondSampicDQMSource::PotPlots::PotPlots(DQMStore::IBooker &ibooker, unsigned int id, bool plotOnline) {
0258 std::string path, title;
0259 CTPPSDiamondDetId(id).rpName(path, CTPPSDiamondDetId::nPath);
0260 ibooker.setCurrentFolder(DiamondSampicDQMSource::changePathToSampic(path));
0261
0262 CTPPSDiamondDetId(id).rpName(title, CTPPSDiamondDetId::nFull);
0263
0264 digiDistribution =
0265 ibooker.book2D("digi distribution", title + " digi distribution;plane;channel", 10, -0.5, 4.5, 12, 0, 12);
0266
0267 hitDistribution2d = ibooker.book2D("hits in planes",
0268 title + " hits in planes;plane number;x (mm)",
0269 10,
0270 -0.5,
0271 4.5,
0272 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0273 -0.5,
0274 18.5);
0275 hitDistribution2dWithTime = ibooker.book2D("hits in planes with time",
0276 title + " hits in planes with time;plane number;x (mm)",
0277 10,
0278 -0.5,
0279 4.5,
0280 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0281 -0.5,
0282 18.5);
0283
0284 recHitTime = ibooker.book1D("recHit time", title + " recHit time; t (ns)", 500, -25, 25);
0285 trackDistribution = ibooker.book1D(
0286 "tracks", title + " tracks;x (mm)", 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM, -1, 18);
0287
0288 if (plotOnline) {
0289 hitDistribution2d_lumisection =
0290 ibooker.book2D("hits in planes lumisection",
0291 title + " hits in planes in the last lumisection;plane number;x (mm)",
0292 18,
0293 -0.5,
0294 4,
0295 15. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0296 0,
0297 15);
0298 triggerCellTime = ibooker.book1D("trigger cell time", title + " Trigger Cell Time; t (ns)", 390, -25, 25);
0299 activityPerBX = ibooker.book1D("activity per BX CMS", title + " Activity per BX;Event.BX", 3600, -1.5, 3598. + 0.5);
0300 amplitude = ibooker.book1D("amplitude", title + " amplitude above baseline; amplitude (V)", 50, 0, 1);
0301 baselineRMS = ibooker.book2D("noise RMS", title + " noise RMS (V);plane;channel", 10, -0.5, 4.5, 12, 0, 12);
0302 meanAmplitude =
0303 ibooker.book2D("mean amplitude", title + " Mean Amplitude (V);plane;channel", 10, -0.5, 4.5, 12, 0, 12);
0304 cellOfMax = ibooker.book2D("cell of max", title + " cell of max (0-23);plane;channel", 10, -0.5, 4.5, 12, 0, 12);
0305
0306
0307
0308 planesWithDigis = ibooker.book1D("active planes digis",
0309 title + " active planes with digis sent (per event);number of active planes",
0310 6,
0311 -0.5,
0312 5.5);
0313 planesWithTime = ibooker.book1D(
0314 "active planes with time", title + " active planes with time (per event);number of active planes", 6, -0.5, 5.5);
0315
0316 dataSamplesRaw = ibooker.book1D("raw Samples", title + " Raw Samples; ADC", 256, 0, 256);
0317
0318 baseline = ibooker.book2D("baseline", title + " baseline (V);plane;channel", 10, -0.5, 4.5, 12, 0, 12);
0319 noiseRMS = ibooker.book2D("noise RMS", title + " noise RMS (V);plane;channel", 10, -0.5, 4.5, 12, 0, 12);
0320
0321 digiSent = ibooker.book2D(
0322 "digis sent", title + " digi sent (sampic);board + 0.5 sampic;channel", 14, -0.5, 6.5, 16, 0, 16);
0323 digiAll =
0324 ibooker.book2D("all digis", title + " all digis(sampic);board + 0.5 sampic;channel", 14, -0.5, 6.5, 16, 0, 16);
0325 digiSentPercentage = ibooker.book2D("sent digis percentage",
0326 title + " sent digis percentage (sampic);board + 0.5 sampic;channel",
0327 14,
0328 -0.5,
0329 6.5,
0330 16,
0331 0,
0332 16);
0333 }
0334 }
0335
0336
0337
0338 DiamondSampicDQMSource::PlanePlots::PlanePlots(DQMStore::IBooker &ibooker, unsigned int id) {
0339 std::string path, title;
0340 CTPPSDiamondDetId(id).planeName(path, CTPPSDiamondDetId::nPath);
0341 ibooker.setCurrentFolder(DiamondSampicDQMSource::changePathToSampic(path));
0342
0343 CTPPSDiamondDetId(id).planeName(title, CTPPSDiamondDetId::nFull);
0344
0345 digiDistribution = ibooker.book1D("digi distribution", title + " digi distribution;channel", 12, 0, 12);
0346
0347 hitProfile = ibooker.book1D("hit distribution with time",
0348 title + " hit distribution (with time);x(mm)",
0349 30. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0350 0,
0351 30);
0352
0353 hitMultiplicity = ibooker.book1D("channels per plane", title + " channels per plane; ch per plane", 13, -0.5, 12.5);
0354
0355 hitMultiplicityWithTime = ibooker.book1D(
0356 "channels per plane with time", title + " channels per plane with time; ch per plane", 13, -0.5, 12.5);
0357 }
0358
0359
0360
0361 DiamondSampicDQMSource::ChannelPlots::ChannelPlots(DQMStore::IBooker &ibooker, unsigned int id) {
0362 std::string path, title;
0363 CTPPSDiamondDetId(id).channelName(path, CTPPSDiamondDetId::nPath);
0364 ibooker.setCurrentFolder(DiamondSampicDQMSource::changePathToSampic(path));
0365
0366 CTPPSDiamondDetId(id).channelName(title, CTPPSDiamondDetId::nFull);
0367
0368 activityPerBX = ibooker.book1D("activity per BX", title + " Activity per BX;Event.BX", 1000, -1.5, 998. + 0.5);
0369 dataSamplesRaw = ibooker.book1D("raw samples", title + " Raw Samples; ADC", 256, 0, 256);
0370 cellOfMax = ibooker.book1D("cell of max", title + " cell of max; cell", 24, 0, 24);
0371
0372 triggerCellTime = ibooker.book1D("sampic trigger time", title + " Sampic Trigger Time; t (ns)", 100, -25, 25);
0373 recHitTime = ibooker.book1D("recHit Time", title + " recHit Time; t (ns)", 500, -25, 25);
0374 amplitude = ibooker.book1D("amplitude", title + " amplitude above baseline; amplitude (V)", 50, 0, 1);
0375 noiseSamples = ibooker.book1D("noise samples", title + " noise samples; V", 50, 0, 1);
0376
0377
0378
0379 }
0380
0381
0382
0383 DiamondSampicDQMSource::DiamondSampicDQMSource(const edm::ParameterSet &ps)
0384 : tokenLocalTrack_(
0385 consumes<edm::DetSetVector<TotemRPLocalTrack>>(ps.getUntrackedParameter<edm::InputTag>("tagLocalTrack"))),
0386 tokenDigi_(consumes<edm::DetSetVector<TotemTimingDigi>>(ps.getUntrackedParameter<edm::InputTag>("tagDigi"))),
0387 tokenRecHit_(
0388 consumes<edm::DetSetVector<TotemTimingRecHit>>(ps.getUntrackedParameter<edm::InputTag>("tagRecHits"))),
0389 tokenTrack_(
0390 consumes<edm::DetSetVector<TotemTimingLocalTrack>>(ps.getUntrackedParameter<edm::InputTag>("tagTracks"))),
0391 tokenFEDInfo_(consumes<std::vector<TotemFEDInfo>>(ps.getUntrackedParameter<edm::InputTag>("tagFEDInfo"))),
0392 ctppsGeometryRunToken_(esConsumes<CTPPSGeometry, VeryForwardRealGeometryRecord, edm::Transition::BeginRun>()),
0393 samplesForNoise_(ps.getUntrackedParameter<unsigned int>("samplesForNoise", 5)),
0394 verbosity_(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
0395 plotOnline_(ps.getUntrackedParameter<bool>("plotOnline", true)),
0396 perLSsaving_(ps.getUntrackedParameter<bool>("perLSsaving", false)),
0397 trackCorrelationThreshold_(ps.getUntrackedParameter<unsigned int>("trackCorrelationThreshold", 3)),
0398 timeOfPreviousEvent_(0) {}
0399
0400
0401
0402 DiamondSampicDQMSource::~DiamondSampicDQMSource() {}
0403
0404
0405
0406 void DiamondSampicDQMSource::dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) {
0407
0408 const CTPPSGeometry *geom = &iSetup.getData(ctppsGeometryRunToken_);
0409 for (auto it = geom->beginSensor(); it != geom->endSensor(); it++) {
0410 if (!CTPPSDiamondDetId::check(it->first))
0411 continue;
0412 const CTPPSDiamondDetId detid(it->first);
0413
0414 const DetGeomDesc *det = geom->sensorNoThrow(detid);
0415 if (det)
0416 horizontalShiftOfDiamond_[detid.rpId()] = det->translation().x() - det->getDiamondDimensions().xHalfWidth;
0417 else
0418 edm::LogProblem("DiamondSampicCalibrationDQMSource") << "ERROR: no descriptor for detId";
0419 }
0420
0421 }
0422
0423
0424
0425 void DiamondSampicDQMSource::bookHistograms(DQMStore::IBooker &ibooker,
0426 const edm::Run &,
0427 const edm::EventSetup &iSetup) {
0428 ibooker.cd();
0429 ibooker.setCurrentFolder("CTPPS/DiamondSampic");
0430
0431 const CTPPSGeometry *geom = &iSetup.getData(ctppsGeometryRunToken_);
0432 for (auto it = geom->beginSensor(); it != geom->endSensor(); ++it) {
0433 if (!CTPPSDiamondDetId::check(it->first))
0434 continue;
0435 const CTPPSDiamondDetId detid(it->first);
0436
0437 sectorPlots_[detid.armId()] = SectorPlots(ibooker, detid.armId(), plotOnline_);
0438
0439 const CTPPSDiamondDetId rpId(detid.arm(), detid.station(), detid.rp());
0440 potPlots_[rpId] = PotPlots(ibooker, rpId, plotOnline_);
0441
0442 if (plotOnline_) {
0443 globalPlot_ = GlobalPlots(ibooker);
0444 const CTPPSDiamondDetId plId(detid.arm(), detid.station(), detid.rp(), detid.plane());
0445 planePlots_[plId] = PlanePlots(ibooker, plId);
0446
0447 const CTPPSDiamondDetId chId(detid.arm(), detid.station(), detid.rp(), detid.plane(), detid.channel());
0448 channelPlots_[chId] = ChannelPlots(ibooker, chId);
0449 }
0450 }
0451 }
0452
0453
0454
0455 std::shared_ptr<totemds::Cache> DiamondSampicDQMSource::globalBeginLuminosityBlock(const edm::LuminosityBlock &,
0456 const edm::EventSetup &) const {
0457 auto d = std::make_shared<totemds::Cache>();
0458 d->hitDistribution2dMap.reserve(potPlots_.size());
0459 if (!perLSsaving_ && plotOnline_)
0460 for (auto &plot : potPlots_)
0461 d->hitDistribution2dMap[plot.first] =
0462 std::unique_ptr<TH2F>(static_cast<TH2F *>(plot.second.hitDistribution2d_lumisection->getTH2F()->Clone()));
0463 return d;
0464 }
0465
0466
0467
0468 void DiamondSampicDQMSource::analyze(const edm::Event &event, const edm::EventSetup &eventSetup) {
0469
0470 edm::Handle<edm::DetSetVector<TotemTimingDigi>> timingDigis;
0471 event.getByToken(tokenDigi_, timingDigis);
0472
0473 edm::Handle<std::vector<TotemFEDInfo>> fedInfo;
0474 event.getByToken(tokenFEDInfo_, fedInfo);
0475
0476 edm::Handle<edm::DetSetVector<TotemTimingRecHit>> timingRecHits;
0477 event.getByToken(tokenRecHit_, timingRecHits);
0478
0479 edm::Handle<edm::DetSetVector<TotemTimingLocalTrack>> timingLocalTracks;
0480 event.getByToken(tokenTrack_, timingLocalTracks);
0481
0482
0483 bool valid = true;
0484 valid &= timingDigis.isValid();
0485
0486
0487 if (!valid) {
0488 if (verbosity_) {
0489 edm::LogProblem("DiamondSampicDQMSource")
0490 << "ERROR in DiamondSampicDQMSource::analyze > some of the required inputs "
0491 "are not valid. Skipping this event.\n"
0492 << " timingDigis.isValid = " << timingDigis.isValid() << "\n"
0493 << " fedInfo.isValid = " << fedInfo.isValid();
0494 }
0495
0496 return;
0497 }
0498
0499
0500 std::set<uint8_t> boardSet;
0501 std::unordered_map<unsigned int, unsigned int> channelsPerPlane;
0502 std::unordered_map<unsigned int, unsigned int> channelsPerPlaneWithTime;
0503
0504 auto lumiCache = luminosityBlockCache(event.getLuminosityBlock().index());
0505
0506 for (const auto &digis : *timingDigis) {
0507 const CTPPSDiamondDetId detId(digis.detId());
0508 CTPPSDiamondDetId detId_pot(digis.detId());
0509 detId_pot.setPlane(0);
0510 detId_pot.setChannel(0);
0511 CTPPSDiamondDetId detId_plane(digis.detId());
0512 detId_plane.setChannel(0);
0513
0514 for (const auto &digi : digis) {
0515
0516 if (potPlots_.find(detId_pot) != potPlots_.end()) {
0517 potPlots_[detId_pot].digiDistribution->Fill(detId.plane(), detId.channel());
0518
0519 if (plotOnline_) {
0520 potPlots_[detId_pot].activityPerBX->Fill(event.bunchCrossing());
0521
0522 for (auto it = digi.samplesBegin(); it != digi.samplesEnd(); ++it) {
0523 potPlots_[detId_pot].dataSamplesRaw->Fill(*it);
0524 }
0525
0526 float boardId = digi.eventInfo().hardwareBoardId() + 0.5 * digi.eventInfo().hardwareSampicId();
0527 potPlots_[detId_pot].digiSent->Fill(boardId, digi.hardwareChannelId());
0528
0529 if (boardSet.find(digi.eventInfo().hardwareId()) == boardSet.end()) {
0530
0531 boardSet.insert(digi.eventInfo().hardwareId());
0532 std::bitset<16> chMap(digi.eventInfo().channelMap());
0533 for (int i = 0; i < 16; ++i) {
0534 if (chMap.test(i)) {
0535 potPlots_[detId_pot].digiAll->Fill(boardId, i);
0536 }
0537 }
0538 }
0539
0540 potPlots_[detId_pot].planesWithDigisSet.insert(detId.plane());
0541 }
0542 }
0543
0544 if (plotOnline_) {
0545
0546 if (planePlots_.find(detId_plane) != planePlots_.end()) {
0547 planePlots_[detId_plane].digiDistribution->Fill(detId.channel());
0548
0549 if (channelsPerPlane.find(detId_plane) != channelsPerPlane.end())
0550 channelsPerPlane[detId_plane]++;
0551 else
0552
0553 channelsPerPlane[detId_plane] = 1;
0554 }
0555
0556
0557 if (channelPlots_.find(detId) != channelPlots_.end()) {
0558 channelPlots_[detId].activityPerBX->Fill(event.bunchCrossing());
0559
0560 for (auto it = digi.samplesBegin(); it != digi.samplesEnd(); ++it)
0561 channelPlots_[detId].dataSamplesRaw->Fill(*it);
0562 for (unsigned short i = 0; i < samplesForNoise_; ++i)
0563 channelPlots_[detId].noiseSamples->Fill(SAMPIC_ADC_V * digi.sampleAt(i));
0564
0565 unsigned int cellOfMax = std::max_element(digi.samplesBegin(), digi.samplesEnd()) - digi.samplesBegin();
0566 channelPlots_[detId].cellOfMax->Fill((int)cellOfMax);
0567
0568
0569
0570
0571 ++(lumiCache->hitsCounterMap[detId]);
0572 }
0573 }
0574 }
0575 }
0576
0577
0578 for (const auto &rechits : *timingRecHits) {
0579 const CTPPSDiamondDetId detId(rechits.detId());
0580 CTPPSDiamondDetId detId_pot(rechits.detId());
0581 detId_pot.setPlane(0);
0582 detId_pot.setChannel(0);
0583 CTPPSDiamondDetId detId_plane(rechits.detId());
0584 detId_plane.setChannel(0);
0585
0586 for (const auto &rechit : rechits) {
0587 if (potPlots_.find(detId_pot) != potPlots_.end()) {
0588 float UFSDShift = 0.0;
0589 if (rechit.yWidth() < 3)
0590 UFSDShift = 0.5;
0591
0592 TH2F *hitHistoTmp = potPlots_[detId_pot].hitDistribution2d->getTH2F();
0593 TAxis *hitHistoTmpYAxis = hitHistoTmp->GetYaxis();
0594 int startBin =
0595 hitHistoTmpYAxis->FindBin(rechit.x() - horizontalShiftOfDiamond_[detId_pot] - 0.5 * rechit.xWidth());
0596 const int numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0597 for (int i = 0; i < numOfBins; ++i) {
0598 potPlots_[detId_pot].hitDistribution2d->Fill(detId.plane() + UFSDShift,
0599 hitHistoTmpYAxis->GetBinCenter(startBin + i));
0600 if (!perLSsaving_ && plotOnline_)
0601 potPlots_[detId_pot].hitDistribution2d_lumisection->Fill(detId.plane() + UFSDShift,
0602 hitHistoTmpYAxis->GetBinCenter(startBin + i));
0603 }
0604
0605
0606 if (rechit.time() != TotemTimingRecHit::NO_T_AVAILABLE) {
0607 for (int i = 0; i < numOfBins; ++i)
0608 potPlots_[detId_pot].hitDistribution2dWithTime->Fill(detId.plane() + UFSDShift,
0609 hitHistoTmpYAxis->GetBinCenter(startBin + i));
0610
0611 potPlots_[detId_pot].recHitTime->Fill(rechit.time());
0612 if (plotOnline_) {
0613 potPlots_[detId_pot].amplitude->Fill(rechit.amplitude());
0614 potPlots_[detId_pot].planesWithTimeSet.insert(detId.plane());
0615
0616
0617 if (planePlots_.find(detId_plane) != planePlots_.end()) {
0618 TH1F *hitProfileHistoTmp = planePlots_[detId_plane].hitProfile->getTH1F();
0619 const int startBin = hitProfileHistoTmp->FindBin(rechit.x() - horizontalShiftOfDiamond_[detId_pot] -
0620 0.5 * rechit.xWidth());
0621 for (int i = 0; i < numOfBins; ++i)
0622 hitProfileHistoTmp->Fill(hitProfileHistoTmp->GetBinCenter(startBin + i));
0623
0624 if (channelsPerPlaneWithTime.find(detId_plane) != channelsPerPlaneWithTime.end())
0625 channelsPerPlaneWithTime[detId_plane]++;
0626 else
0627
0628 channelsPerPlaneWithTime[detId_plane] = 1;
0629 }
0630
0631 if (channelPlots_.find(detId) != channelPlots_.end()) {
0632 potPlots_[detId_pot].triggerCellTime->Fill(rechit.sampicThresholdTime());
0633 channelPlots_[detId].triggerCellTime->Fill(rechit.sampicThresholdTime());
0634 channelPlots_[detId].recHitTime->Fill(rechit.time());
0635 channelPlots_[detId].amplitude->Fill(rechit.amplitude());
0636 }
0637 }
0638 }
0639 }
0640 }
0641 }
0642
0643
0644
0645 for (const auto &tracks : *timingLocalTracks) {
0646 CTPPSDiamondDetId detId_pot(tracks.detId());
0647 detId_pot.setPlane(0);
0648 detId_pot.setChannel(0);
0649 const CTPPSDiamondDetId detId_near(tracks.detId());
0650
0651 for (const auto &track : tracks) {
0652 if (!track.isValid())
0653 continue;
0654 if (potPlots_.find(detId_pot) == potPlots_.end())
0655 continue;
0656
0657 TH1F *trackHistoInTimeTmp = potPlots_[detId_pot].trackDistribution->getTH1F();
0658 const int startBin =
0659 trackHistoInTimeTmp->FindBin(track.x0() - horizontalShiftOfDiamond_[detId_pot] - track.x0Sigma());
0660 const int numOfBins = 2 * track.x0Sigma() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0661 for (int i = 0; i < numOfBins; ++i) {
0662 trackHistoInTimeTmp->Fill(trackHistoInTimeTmp->GetBinCenter(startBin + i));
0663 }
0664
0665
0666 for (const auto &tracks_far : *timingLocalTracks) {
0667 CTPPSDiamondDetId detId_far(tracks_far.detId());
0668 detId_far.setPlane(0);
0669 detId_far.setChannel(0);
0670 if (detId_near.arm() != detId_far.arm() || detId_near.station() == detId_far.station())
0671 continue;
0672 for (const auto &track_far : tracks_far) {
0673 if (!track.isValid())
0674 continue;
0675 if (sectorPlots_.find(detId_far.armId()) == sectorPlots_.end())
0676 continue;
0677 TH2F *trackHistoTmp = sectorPlots_[detId_far.armId()].trackCorrelation->getTH2F();
0678 TAxis *trackHistoTmpXAxis = trackHistoTmp->GetXaxis();
0679 TAxis *trackHistoTmpYAxis = trackHistoTmp->GetYaxis();
0680 const int startBin_far =
0681 trackHistoTmpYAxis->FindBin(track_far.x0() - horizontalShiftOfDiamond_[detId_far] - track_far.x0Sigma());
0682 const int numOfBins_far = 2 * track_far.x0Sigma() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0683 for (int i = 0; i < numOfBins; ++i) {
0684 for (int y = 0; y < numOfBins_far; ++y) {
0685 trackHistoTmp->Fill(trackHistoTmpXAxis->GetBinCenter(startBin + i),
0686 trackHistoTmpYAxis->GetBinCenter(startBin_far + y));
0687 if (tracks.size() < 3 && tracks_far.size() < trackCorrelationThreshold_)
0688 sectorPlots_[detId_far.armId()].trackCorrelationLowMultiplicity->Fill(
0689 trackHistoTmpXAxis->GetBinCenter(startBin + i), trackHistoTmpYAxis->GetBinCenter(startBin_far + y));
0690 }
0691 }
0692 }
0693 }
0694 }
0695 }
0696 if (plotOnline_) {
0697 for (auto &plt : potPlots_) {
0698 plt.second.planesWithDigis->Fill(plt.second.planesWithDigisSet.size());
0699 plt.second.planesWithDigisSet.clear();
0700 plt.second.planesWithTime->Fill(plt.second.planesWithTimeSet.size());
0701 plt.second.planesWithTimeSet.clear();
0702 }
0703
0704 for (const auto &plt : channelsPerPlane)
0705 planePlots_[plt.first].hitMultiplicity->Fill(plt.second);
0706
0707 for (const auto &plt : channelsPerPlaneWithTime)
0708 planePlots_[plt.first].hitMultiplicityWithTime->Fill(plt.second);
0709 }
0710 timeOfPreviousEvent_ = event.time().value();
0711 }
0712
0713
0714
0715 void DiamondSampicDQMSource::globalEndLuminosityBlock(const edm::LuminosityBlock &iLumi, const edm::EventSetup &) {
0716 auto lumiCache = luminosityBlockCache(iLumi.index());
0717 if (!perLSsaving_ && plotOnline_) {
0718 for (auto &plot : potPlots_)
0719 *(plot.second.hitDistribution2d_lumisection->getTH2F()) = *(lumiCache->hitDistribution2dMap[plot.first]);
0720 globalPlot_.digiSentPercentage->Reset();
0721 for (auto &plot : sectorPlots_)
0722 plot.second.digiSentPercentage->Reset();
0723 TH2F *hitHistoGlobalTmp = globalPlot_.digiSentPercentage->getTH2F();
0724 for (auto &plot : potPlots_) {
0725 TH2F *hitHistoTmp = plot.second.digiSentPercentage->getTH2F();
0726 TH2F *histoSent = plot.second.digiSent->getTH2F();
0727 TH2F *histoAll = plot.second.digiAll->getTH2F();
0728
0729 hitHistoTmp->Divide(histoSent, histoAll);
0730 hitHistoTmp->Scale(100);
0731 hitHistoGlobalTmp->Add(hitHistoTmp, 1);
0732
0733 plot.second.baseline->Reset();
0734 plot.second.noiseRMS->Reset();
0735 plot.second.meanAmplitude->Reset();
0736 plot.second.cellOfMax->Reset();
0737
0738 CTPPSDiamondDetId rpId(plot.first);
0739 TH2F *hitHistoSectorTmp = sectorPlots_[rpId.armId()].digiSentPercentage->getTH2F();
0740 hitHistoSectorTmp->Add(hitHistoTmp, 1);
0741
0742 for (auto &chPlot : channelPlots_) {
0743 CTPPSDiamondDetId chId(chPlot.first);
0744 if (chId.arm() == rpId.arm() && chId.rp() == rpId.rp()) {
0745 plot.second.baseline->Fill(chId.plane(), chId.channel(), chPlot.second.noiseSamples->getTH1F()->GetMean());
0746 plot.second.noiseRMS->Fill(chId.plane(), chId.channel(), chPlot.second.noiseSamples->getTH1F()->GetRMS());
0747 plot.second.meanAmplitude->Fill(chId.plane(), chId.channel(), chPlot.second.amplitude->getTH1F()->GetMean());
0748 plot.second.cellOfMax->Fill(chId.plane(), chId.channel(), chPlot.second.cellOfMax->getTH1F()->GetMean());
0749
0750
0751 }
0752 }
0753 }
0754
0755
0756
0757
0758
0759
0760
0761 }
0762 }
0763
0764 DEFINE_FWK_MODULE(DiamondSampicDQMSource);