File indexing completed on 2024-09-11 04:32:26
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include "FWCore/Framework/interface/ESHandle.h"
0020 #include "FWCore/Framework/interface/Event.h"
0021 #include "FWCore/Framework/interface/EventSetup.h"
0022 #include "FWCore/Framework/interface/MakerMacros.h"
0023 #include "FWCore/Framework/interface/Run.h"
0024 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0025 #include "FWCore/Utilities/interface/InputTag.h"
0026
0027 #include "DQMServices/Core/interface/DQMOneEDAnalyzer.h"
0028 #include "DQMServices/Core/interface/DQMStore.h"
0029
0030 #include "DataFormats/CTPPSDigi/interface/TotemFEDInfo.h"
0031 #include "DataFormats/CTPPSDigi/interface/TotemVFATStatus.h"
0032 #include "DataFormats/CTPPSReco/interface/TotemRPLocalTrack.h"
0033 #include "DataFormats/Common/interface/DetSetVector.h"
0034 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
0035 #include "DataFormats/Provenance/interface/EventRange.h"
0036
0037 #include "DataFormats/CTPPSDetId/interface/CTPPSDiamondDetId.h"
0038 #include "DataFormats/CTPPSDigi/interface/TotemTimingDigi.h"
0039 #include "DataFormats/CTPPSReco/interface/TotemTimingRecHit.h"
0040
0041 #include "DataFormats/CTPPSDetId/interface/TotemRPDetId.h"
0042 #include "DataFormats/CTPPSReco/interface/TotemTimingLocalTrack.h"
0043 #include "Geometry/Records/interface/VeryForwardRealGeometryRecord.h"
0044 #include "Geometry/VeryForwardGeometryBuilder/interface/CTPPSGeometry.h"
0045 #include "CondFormats/PPSObjects/interface/PPSTimingCalibration.h"
0046 #include "CondFormats/DataRecord/interface/PPSTimingCalibrationRcd.h"
0047
0048 #include <string>
0049
0050
0051
0052 class DiamondSampicCalibrationDQMSource : public DQMOneEDAnalyzer<> {
0053 public:
0054 DiamondSampicCalibrationDQMSource(const edm::ParameterSet &);
0055 ~DiamondSampicCalibrationDQMSource() override;
0056
0057 protected:
0058 void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override;
0059 void bookHistograms(DQMStore::IBooker &, const edm::Run &, const edm::EventSetup &) override;
0060 void analyze(const edm::Event &, const edm::EventSetup &) override;
0061
0062 private:
0063
0064 static const double DISPLAY_RESOLUTION_FOR_HITS_MM;
0065
0066
0067 static const double INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0068
0069 edm::EDGetTokenT<edm::DetSetVector<TotemTimingDigi>> totemTimingDigiToken_;
0070 edm::EDGetTokenT<edm::DetSetVector<TotemTimingRecHit>> tokenRecHit_;
0071 edm::ESGetToken<PPSTimingCalibration, PPSTimingCalibrationRcd> timingCalibrationToken_;
0072 edm::ESGetToken<CTPPSGeometry, VeryForwardRealGeometryRecord> geomEsToken_;
0073 unsigned int verbosity_;
0074 edm::TimeValue_t timeOfPreviousEvent_;
0075
0076 float verticalShiftBot_, verticalShiftTop_;
0077 std::unordered_map<unsigned int, double> horizontalShiftOfDiamond_;
0078
0079
0080 struct GlobalPlots {
0081 GlobalPlots() {}
0082 GlobalPlots(DQMStore::IBooker &ibooker);
0083 };
0084
0085 GlobalPlots globalPlot_;
0086
0087
0088 struct PotPlots {
0089
0090 MonitorElement *hitDistribution2d = nullptr;
0091 MonitorElement *recHitTime = nullptr;
0092
0093 PotPlots() {}
0094 PotPlots(DQMStore::IBooker &ibooker, unsigned int id);
0095 };
0096
0097 std::unordered_map<unsigned int, PotPlots> potPlots_;
0098
0099
0100 struct PlanePlots {
0101 PlanePlots() {}
0102 PlanePlots(DQMStore::IBooker &ibooker, unsigned int id);
0103 };
0104
0105 std::unordered_map<unsigned int, PlanePlots> planePlots_;
0106
0107
0108 struct ChannelPlots {
0109
0110 MonitorElement *recHitTime = nullptr;
0111
0112 ChannelPlots() {}
0113 ChannelPlots(DQMStore::IBooker &ibooker, unsigned int id);
0114 };
0115
0116 edm::ESHandle<PPSTimingCalibration> hTimingCalib_;
0117 std::unordered_map<unsigned int, ChannelPlots> channelPlots_;
0118 };
0119
0120
0121
0122
0123 const double DiamondSampicCalibrationDQMSource::DISPLAY_RESOLUTION_FOR_HITS_MM = 0.05;
0124 const double DiamondSampicCalibrationDQMSource::INV_DISPLAY_RESOLUTION_FOR_HITS_MM =
0125 1. / DISPLAY_RESOLUTION_FOR_HITS_MM;
0126
0127
0128
0129 DiamondSampicCalibrationDQMSource::GlobalPlots::GlobalPlots(DQMStore::IBooker &ibooker) {
0130 ibooker.setCurrentFolder("CTPPS/TimingFastSilicon");
0131 }
0132
0133
0134
0135 DiamondSampicCalibrationDQMSource::PotPlots::PotPlots(DQMStore::IBooker &ibooker, unsigned int id) {
0136 std::string path, title;
0137 CTPPSDiamondDetId(id).rpName(path, CTPPSDiamondDetId::nPath);
0138 ibooker.setCurrentFolder(path);
0139
0140 CTPPSDiamondDetId(id).rpName(title, CTPPSDiamondDetId::nFull);
0141
0142 hitDistribution2d = ibooker.book2D("hits in planes",
0143 title + " hits in planes;plane number;x (mm)",
0144 10,
0145 -0.5,
0146 4.5,
0147 19. * INV_DISPLAY_RESOLUTION_FOR_HITS_MM,
0148 -0.5,
0149 18.5);
0150
0151 recHitTime = ibooker.book1D("recHit time", title + " time in the recHits; t (ns)", 500, -25, 25);
0152 }
0153
0154
0155
0156 DiamondSampicCalibrationDQMSource::PlanePlots::PlanePlots(DQMStore::IBooker &ibooker, unsigned int id) {
0157 std::string path, title;
0158 CTPPSDiamondDetId(id).planeName(path, CTPPSDiamondDetId::nPath);
0159 ibooker.setCurrentFolder(path);
0160 }
0161
0162
0163
0164 DiamondSampicCalibrationDQMSource::ChannelPlots::ChannelPlots(DQMStore::IBooker &ibooker, unsigned int id) {
0165 std::string path, title;
0166 CTPPSDiamondDetId(id).channelName(path, CTPPSDiamondDetId::nPath);
0167 ibooker.setCurrentFolder(path);
0168
0169 CTPPSDiamondDetId(id).channelName(title, CTPPSDiamondDetId::nFull);
0170 recHitTime = ibooker.book1D("recHit Time", title + " recHit Time; t (ns)", 500, -25, 25);
0171 }
0172
0173
0174
0175 DiamondSampicCalibrationDQMSource::DiamondSampicCalibrationDQMSource(const edm::ParameterSet &ps)
0176 : totemTimingDigiToken_(
0177 consumes<edm::DetSetVector<TotemTimingDigi>>(ps.getUntrackedParameter<edm::InputTag>("totemTimingDigiTag"))),
0178 tokenRecHit_(
0179 consumes<edm::DetSetVector<TotemTimingRecHit>>(ps.getUntrackedParameter<edm::InputTag>("tagRecHits"))),
0180 timingCalibrationToken_(esConsumes<edm::Transition::BeginRun>()),
0181 geomEsToken_(esConsumes<edm::Transition::BeginRun>()),
0182 verbosity_(ps.getUntrackedParameter<unsigned int>("verbosity", 0)),
0183 timeOfPreviousEvent_(0) {}
0184
0185
0186
0187 DiamondSampicCalibrationDQMSource::~DiamondSampicCalibrationDQMSource() {}
0188
0189
0190
0191 void DiamondSampicCalibrationDQMSource::dqmBeginRun(const edm::Run &iRun, const edm::EventSetup &iSetup) {
0192
0193 const auto &geom = iSetup.getData(geomEsToken_);
0194 for (auto it = geom.beginSensor(); it != geom.endSensor(); it++) {
0195 if (!CTPPSDiamondDetId::check(it->first))
0196 continue;
0197 const CTPPSDiamondDetId detid(it->first);
0198
0199 const DetGeomDesc *det = geom.sensorNoThrow(detid);
0200 if (det)
0201 horizontalShiftOfDiamond_[detid.rpId()] = det->translation().x() - det->getDiamondDimensions().xHalfWidth;
0202 else
0203 edm::LogProblem("DiamondSampicCalibrationDQMSource") << "ERROR: no descriptor for detId";
0204 }
0205 }
0206
0207
0208
0209 void DiamondSampicCalibrationDQMSource::bookHistograms(DQMStore::IBooker &ibooker,
0210 const edm::Run &,
0211 const edm::EventSetup &iSetup) {
0212 ibooker.cd();
0213 ibooker.setCurrentFolder("CTPPS");
0214
0215 globalPlot_ = GlobalPlots(ibooker);
0216 const auto &geom = iSetup.getData(geomEsToken_);
0217 for (auto it = geom.beginSensor(); it != geom.endSensor(); ++it) {
0218 if (!CTPPSDiamondDetId::check(it->first))
0219 continue;
0220 const CTPPSDiamondDetId detid(it->first);
0221
0222 const CTPPSDiamondDetId rpId(detid.arm(), detid.station(), detid.rp());
0223 potPlots_[rpId] = PotPlots(ibooker, rpId);
0224
0225 const CTPPSDiamondDetId plId(detid.arm(), detid.station(), detid.rp(), detid.plane());
0226 planePlots_[plId] = PlanePlots(ibooker, plId);
0227
0228 const CTPPSDiamondDetId chId(detid.arm(), detid.station(), detid.rp(), detid.plane(), detid.channel());
0229 channelPlots_[chId] = ChannelPlots(ibooker, chId);
0230 }
0231 hTimingCalib_ = iSetup.getHandle(timingCalibrationToken_);
0232 }
0233
0234
0235
0236 void DiamondSampicCalibrationDQMSource::analyze(const edm::Event &event, const edm::EventSetup &eventSetup) {
0237 PPSTimingCalibration calib = *hTimingCalib_;
0238
0239 edm::Handle<edm::DetSetVector<TotemTimingRecHit>> timingRecHits;
0240 event.getByToken(tokenRecHit_, timingRecHits);
0241
0242 edm::Handle<edm::DetSetVector<TotemTimingDigi>> timingDigi;
0243 event.getByToken(totemTimingDigiToken_, timingDigi);
0244
0245 std::unordered_map<uint32_t, uint32_t> detIdToHw;
0246
0247 for (const auto &digis : *timingDigi) {
0248 const CTPPSDiamondDetId detId(digis.detId());
0249 for (const auto &digi : digis)
0250 detIdToHw[detId] = digi.hardwareId();
0251 }
0252
0253
0254 std::set<uint8_t> boardSet;
0255 std::unordered_map<unsigned int, unsigned int> channelsPerPlane;
0256 std::unordered_map<unsigned int, unsigned int> channelsPerPlaneWithTime;
0257
0258
0259
0260 for (const auto &rechits : *timingRecHits) {
0261 const CTPPSDiamondDetId detId(rechits.detId());
0262 CTPPSDiamondDetId detId_pot(rechits.detId());
0263 detId_pot.setPlane(0);
0264 detId_pot.setChannel(0);
0265 CTPPSDiamondDetId detId_plane(rechits.detId());
0266 detId_plane.setChannel(0);
0267
0268 for (const auto &rechit : rechits) {
0269 if (potPlots_.find(detId_pot) != potPlots_.end()) {
0270 float UFSDShift = 0.0;
0271 if (rechit.yWidth() < 3)
0272 UFSDShift = 0.5;
0273
0274 TH2F *hitHistoTmp = potPlots_[detId_pot].hitDistribution2d->getTH2F();
0275 TAxis *hitHistoTmpYAxis = hitHistoTmp->GetYaxis();
0276 int startBin =
0277 hitHistoTmpYAxis->FindBin(rechit.x() - horizontalShiftOfDiamond_[detId_pot] - 0.5 * rechit.xWidth());
0278 int numOfBins = rechit.xWidth() * INV_DISPLAY_RESOLUTION_FOR_HITS_MM;
0279 for (int i = 0; i < numOfBins; ++i)
0280 potPlots_[detId_pot].hitDistribution2d->Fill(detId.plane() + UFSDShift,
0281 hitHistoTmpYAxis->GetBinCenter(startBin + i));
0282
0283
0284 if (rechit.time() != TotemTimingRecHit::NO_T_AVAILABLE) {
0285 int db = (detIdToHw[detId] & 0xE0) >> 5;
0286 int sampic = (detIdToHw[detId] & 0x10) >> 4;
0287 int channel = (detIdToHw[detId] & 0x0F);
0288 double offset = calib.timeOffset(db, sampic, channel);
0289 potPlots_[detId_pot].recHitTime->Fill(rechit.time() + offset);
0290 if (channelPlots_.find(detId) != channelPlots_.end())
0291 channelPlots_[detId].recHitTime->Fill(rechit.time() + offset);
0292 }
0293 }
0294 }
0295 }
0296
0297
0298 timeOfPreviousEvent_ = event.time().value();
0299 }
0300
0301 DEFINE_FWK_MODULE(DiamondSampicCalibrationDQMSource);