File indexing completed on 2024-04-06 12:32:43
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include <string>
0015
0016 #include "FWCore/Framework/interface/Frameworkfwd.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020
0021 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0022 #include "DQMServices/Core/interface/DQMStore.h"
0023
0024 #include "DataFormats/Common/interface/ValidHandle.h"
0025 #include "DataFormats/ForwardDetId/interface/BTLDetId.h"
0026 #include "DataFormats/FTLDigi/interface/FTLDigiCollections.h"
0027
0028 #include "Geometry/Records/interface/MTDDigiGeometryRecord.h"
0029 #include "Geometry/Records/interface/MTDTopologyRcd.h"
0030 #include "Geometry/MTDGeometryBuilder/interface/MTDGeometry.h"
0031 #include "Geometry/MTDNumberingBuilder/interface/MTDTopology.h"
0032
0033 #include "Geometry/MTDGeometryBuilder/interface/ProxyMTDTopology.h"
0034 #include "Geometry/MTDGeometryBuilder/interface/RectangularMTDTopology.h"
0035
0036 #include "Geometry/MTDCommonData/interface/MTDTopologyMode.h"
0037
0038 class BtlDigiHitsValidation : public DQMEDAnalyzer {
0039 public:
0040 explicit BtlDigiHitsValidation(const edm::ParameterSet&);
0041 ~BtlDigiHitsValidation() override;
0042
0043 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0044
0045 private:
0046 void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0047
0048 void analyze(const edm::Event&, const edm::EventSetup&) override;
0049
0050
0051
0052 const std::string folder_;
0053 const bool optionalPlots_;
0054
0055 edm::EDGetTokenT<BTLDigiCollection> btlDigiHitsToken_;
0056
0057 edm::ESGetToken<MTDGeometry, MTDDigiGeometryRecord> mtdgeoToken_;
0058 edm::ESGetToken<MTDTopology, MTDTopologyRcd> mtdtopoToken_;
0059
0060
0061
0062 MonitorElement* meNhits_[2];
0063
0064 MonitorElement* meHitCharge_[2];
0065 MonitorElement* meHitTime_[2];
0066
0067 MonitorElement* meOccupancy_[2];
0068
0069
0070 MonitorElement* meLocalOccupancy_[2];
0071 MonitorElement* meHitXlocal_[2];
0072 MonitorElement* meHitYlocal_[2];
0073 MonitorElement* meHitZlocal_[2];
0074
0075 MonitorElement* meHitX_[2];
0076 MonitorElement* meHitY_[2];
0077 MonitorElement* meHitZ_[2];
0078 MonitorElement* meHitPhi_[2];
0079 MonitorElement* meHitEta_[2];
0080
0081 MonitorElement* meHitTvsQ_[2];
0082 MonitorElement* meHitQvsPhi_[2];
0083 MonitorElement* meHitQvsEta_[2];
0084 MonitorElement* meHitQvsZ_[2];
0085 MonitorElement* meHitTvsPhi_[2];
0086 MonitorElement* meHitTvsEta_[2];
0087 MonitorElement* meHitTvsZ_[2];
0088 };
0089
0090
0091 BtlDigiHitsValidation::BtlDigiHitsValidation(const edm::ParameterSet& iConfig)
0092 : folder_(iConfig.getParameter<std::string>("folder")),
0093 optionalPlots_(iConfig.getParameter<bool>("optionalPlots")) {
0094 btlDigiHitsToken_ = consumes<BTLDigiCollection>(iConfig.getParameter<edm::InputTag>("inputTag"));
0095 mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
0096 mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
0097 }
0098
0099 BtlDigiHitsValidation::~BtlDigiHitsValidation() {}
0100
0101
0102 void BtlDigiHitsValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0103 using namespace edm;
0104
0105 auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
0106 const MTDGeometry* geom = geometryHandle.product();
0107
0108 auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
0109 const MTDTopology* topology = topologyHandle.product();
0110
0111 auto btlDigiHitsHandle = makeValid(iEvent.getHandle(btlDigiHitsToken_));
0112
0113
0114
0115 unsigned int n_digi_btl[2] = {0, 0};
0116
0117 for (const auto& dataFrame : *btlDigiHitsHandle) {
0118 BTLDetId detId = dataFrame.id();
0119 DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
0120 const MTDGeomDet* thedet = geom->idToDet(geoId);
0121 if (thedet == nullptr)
0122 throw cms::Exception("BtlDigiHitsValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
0123 << detId.rawId() << ") is invalid!" << std::dec << std::endl;
0124 const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
0125 const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
0126
0127 Local3DPoint local_point(0., 0., 0.);
0128 local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
0129 const auto& global_point = thedet->toGlobal(local_point);
0130
0131 const auto& sample_L = dataFrame.sample(0);
0132 const auto& sample_R = dataFrame.sample(1);
0133
0134 uint32_t adc[2] = {sample_L.data(), sample_R.data()};
0135 uint32_t tdc[2] = {sample_L.toa(), sample_R.toa()};
0136
0137 for (int iside = 0; iside < 2; ++iside) {
0138 if (adc[iside] == 0)
0139 continue;
0140
0141 meHitCharge_[iside]->Fill(adc[iside]);
0142 meHitTime_[iside]->Fill(tdc[iside]);
0143
0144 meOccupancy_[iside]->Fill(global_point.z(), global_point.phi());
0145
0146 if (optionalPlots_) {
0147 meLocalOccupancy_[iside]->Fill(local_point.x(), local_point.y());
0148 meHitXlocal_[iside]->Fill(local_point.x());
0149 meHitYlocal_[iside]->Fill(local_point.y());
0150 meHitZlocal_[iside]->Fill(local_point.z());
0151 }
0152
0153 meHitX_[iside]->Fill(global_point.x());
0154 meHitY_[iside]->Fill(global_point.y());
0155 meHitZ_[iside]->Fill(global_point.z());
0156 meHitPhi_[iside]->Fill(global_point.phi());
0157 meHitEta_[iside]->Fill(global_point.eta());
0158
0159 meHitTvsQ_[iside]->Fill(adc[iside], tdc[iside]);
0160 meHitQvsPhi_[iside]->Fill(global_point.phi(), adc[iside]);
0161 meHitQvsEta_[iside]->Fill(global_point.eta(), adc[iside]);
0162 meHitQvsZ_[iside]->Fill(global_point.z(), adc[iside]);
0163 meHitTvsPhi_[iside]->Fill(global_point.phi(), tdc[iside]);
0164 meHitTvsEta_[iside]->Fill(global_point.eta(), tdc[iside]);
0165 meHitTvsZ_[iside]->Fill(global_point.z(), tdc[iside]);
0166
0167 n_digi_btl[iside]++;
0168
0169 }
0170
0171 }
0172
0173 if (n_digi_btl[0] > 0)
0174 meNhits_[0]->Fill(log10(n_digi_btl[0]));
0175 if (n_digi_btl[1] > 0)
0176 meNhits_[1]->Fill(log10(n_digi_btl[1]));
0177 }
0178
0179
0180 void BtlDigiHitsValidation::bookHistograms(DQMStore::IBooker& ibook,
0181 edm::Run const& run,
0182 edm::EventSetup const& iSetup) {
0183 ibook.setCurrentFolder(folder_);
0184
0185
0186 meNhits_[0] = ibook.book1D("BtlNhitsL", "Number of BTL DIGI hits (L);log_{10}(N_{DIGI})", 100, 0., 5.25);
0187 meNhits_[1] = ibook.book1D("BtlNhitsR", "Number of BTL DIGI hits (R);log_{10}(N_{DIGI})", 100, 0., 5.25);
0188
0189 meHitCharge_[0] = ibook.book1D("BtlHitChargeL", "BTL DIGI hits charge (L);Q_{DIGI} [ADC counts]", 100, 0., 1024.);
0190 meHitCharge_[1] = ibook.book1D("BtlHitChargeR", "BTL DIGI hits charge (R);Q_{DIGI} [ADC counts]", 100, 0., 1024.);
0191 meHitTime_[0] = ibook.book1D("BtlHitTimeL", "BTL DIGI hits ToA (L);ToA_{DIGI} [TDC counts]", 100, 0., 1024.);
0192 meHitTime_[1] = ibook.book1D("BtlHitTimeR", "BTL DIGI hits ToA (R);ToA_{DIGI} [TDC counts]", 100, 0., 1024.);
0193 meOccupancy_[0] = ibook.book2D("BtlOccupancyL",
0194 "BTL DIGI hits occupancy (L);Z_{DIGI} [cm]; #phi_{DIGI} [rad]",
0195 65,
0196 -260.,
0197 260.,
0198 126,
0199 -3.15,
0200 3.15);
0201 meOccupancy_[1] = ibook.book2D("BtlOccupancyR",
0202 "BTL DIGI hits occupancy (R);Z_{DIGI} [cm]; #phi_{DIGI} [rad]",
0203 65,
0204 -260.,
0205 260.,
0206 126,
0207 -3.15,
0208 3.15);
0209 if (optionalPlots_) {
0210 meLocalOccupancy_[0] = ibook.book2D("BtlLocalOccupancyL",
0211 "BTL DIGI hits local occupancy (L);X_{DIGI} [cm]; Y_{DIGI} [cm]",
0212 100,
0213 -10.,
0214 10,
0215 60,
0216 -3.,
0217 3.);
0218 meLocalOccupancy_[1] = ibook.book2D(
0219 "BtlLocalOccupancyR", "BTL DIGI hits occupancy (R);X_{DIGI} [cm]; Y_{DIGI} [cm]", 100, -10., 10., 60, -3., 3.);
0220 meHitXlocal_[0] = ibook.book1D("BtlHitXlocalL", "BTL DIGI local X (L);X_{DIGI}^{LOC} [cm]", 100, -10., 10.);
0221 meHitXlocal_[1] = ibook.book1D("BtlHitXlocalR", "BTL DIGI local X (R);X_{DIGI}^{LOC} [cm]", 100, -10., 10.);
0222 meHitYlocal_[0] = ibook.book1D("BtlHitYlocalL", "BTL DIGI local Y (L);Y_{DIGI}^{LOC} [cm]", 60, -3., 3.);
0223 meHitYlocal_[1] = ibook.book1D("BtlHitYlocalR", "BTL DIGI local Y (R);Y_{DIGI}^{LOC} [cm]", 60, -3., 3.);
0224 meHitZlocal_[0] = ibook.book1D("BtlHitZlocalL", "BTL DIGI local z (L);z_{DIGI}^{LOC} [cm]", 10, -1, 1);
0225 meHitZlocal_[1] = ibook.book1D("BtlHitZlocalR", "BTL DIGI local z (R);z_{DIGI}^{LOC} [cm]", 10, -1, 1);
0226 }
0227
0228 meHitX_[0] = ibook.book1D("BtlHitXL", "BTL DIGI hits X (L);X_{DIGI} [cm]", 60, -120., 120.);
0229 meHitX_[1] = ibook.book1D("BtlHitXR", "BTL DIGI hits X (R);X_{DIGI} [cm]", 60, -120., 120.);
0230 meHitY_[0] = ibook.book1D("BtlHitYL", "BTL DIGI hits Y (L);Y_{DIGI} [cm]", 60, -120., 120.);
0231 meHitY_[1] = ibook.book1D("BtlHitYR", "BTL DIGI hits Y (R);Y_{DIGI} [cm]", 60, -120., 120.);
0232 meHitZ_[0] = ibook.book1D("BtlHitZL", "BTL DIGI hits Z (L);Z_{DIGI} [cm]", 100, -260., 260.);
0233 meHitZ_[1] = ibook.book1D("BtlHitZR", "BTL DIGI hits Z (R);Z_{DIGI} [cm]", 100, -260., 260.);
0234 meHitPhi_[0] = ibook.book1D("BtlHitPhiL", "BTL DIGI hits #phi (L);#phi_{DIGI} [rad]", 126, -3.15, 3.15);
0235 meHitPhi_[1] = ibook.book1D("BtlHitPhiR", "BTL DIGI hits #phi (R);#phi_{DIGI} [rad]", 126, -3.15, 3.15);
0236 meHitEta_[0] = ibook.book1D("BtlHitEtaL", "BTL DIGI hits #eta (L);#eta_{DIGI}", 100, -1.55, 1.55);
0237 meHitEta_[1] = ibook.book1D("BtlHitEtaR", "BTL DIGI hits #eta (R);#eta_{DIGI}", 100, -1.55, 1.55);
0238
0239 meHitTvsQ_[0] = ibook.bookProfile("BtlHitTvsQL",
0240 "BTL DIGI ToA vs charge (L);Q_{DIGI} [ADC counts];ToA_{DIGI} [TDC counts]",
0241 50,
0242 0.,
0243 1024.,
0244 0.,
0245 1024.);
0246 meHitTvsQ_[1] = ibook.bookProfile("BtlHitTvsQR",
0247 "BTL DIGI ToA vs charge (R);Q_{DIGI} [ADC counts];ToA_{DIGI} [TDC counts]",
0248 50,
0249 0.,
0250 1024.,
0251 0.,
0252 1024.);
0253 meHitQvsPhi_[0] = ibook.bookProfile("BtlHitQvsPhiL",
0254 "BTL DIGI charge vs #phi (L);#phi_{DIGI} [rad];Q_{DIGI} [ADC counts]",
0255 50,
0256 -3.15,
0257 3.15,
0258 0.,
0259 1024.);
0260 meHitQvsPhi_[1] = ibook.bookProfile("BtlHitQvsPhiR",
0261 "BTL DIGI charge vs #phi (R);#phi_{DIGI} [rad];Q_{DIGI} [ADC counts]",
0262 50,
0263 -3.15,
0264 3.15,
0265 0.,
0266 1024.);
0267 meHitQvsEta_[0] = ibook.bookProfile(
0268 "BtlHitQvsEtaL", "BTL DIGI charge vs #eta (L);#eta_{DIGI};Q_{DIGI} [ADC counts]", 50, -1.55, 1.55, 0., 1024.);
0269 meHitQvsEta_[1] = ibook.bookProfile(
0270 "BtlHitQvsEtaR", "BTL DIGI charge vs #eta (R);#eta_{DIGI};Q_{DIGI} [ADC counts]", 50, -1.55, 1.55, 0., 1024.);
0271 meHitQvsZ_[0] = ibook.bookProfile(
0272 "BtlHitQvsZL", "BTL DIGI charge vs Z (L);Z_{DIGI} [cm];Q_{DIGI} [ADC counts]", 50, -260., 260., 0., 1024.);
0273 meHitQvsZ_[1] = ibook.bookProfile(
0274 "BtlHitQvsZR", "BTL DIGI charge vs Z (R);Z_{DIGI} [cm];Q_{DIGI} [ADC counts]", 50, -260., 260., 0., 1024.);
0275 meHitTvsPhi_[0] = ibook.bookProfile(
0276 "BtlHitTvsPhiL", "BTL DIGI ToA vs #phi (L);#phi_{DIGI} [rad];ToA_{DIGI} [TDC counts]", 50, -3.15, 3.15, 0., 1024.);
0277 meHitTvsPhi_[1] = ibook.bookProfile(
0278 "BtlHitTvsPhiR", "BTL DIGI ToA vs #phi (R);#phi_{DIGI} [rad];ToA_{DIGI} [TDC counts]", 50, -3.15, 3.15, 0., 1024.);
0279 meHitTvsEta_[0] = ibook.bookProfile(
0280 "BtlHitTvsEtaL", "BTL DIGI ToA vs #eta (L);#eta_{DIGI};ToA_{DIGI} [TDC counts]", 50, -1.55, 1.55, 0., 1024.);
0281 meHitTvsEta_[1] = ibook.bookProfile(
0282 "BtlHitTvsEtaR", "BTL DIGI ToA vs #eta (R);#eta_{DIGI};ToA_{DIGI} [TDC counts]", 50, -1.55, 1.55, 0., 1024.);
0283 meHitTvsZ_[0] = ibook.bookProfile(
0284 "BtlHitTvsZL", "BTL DIGI ToA vs Z (L);Z_{DIGI} [cm];ToA_{DIGI} [TDC counts]", 50, -260., 260., 0., 1024.);
0285 meHitTvsZ_[1] = ibook.bookProfile(
0286 "BtlHitTvsZR", "BTL DIGI ToA vs Z (R);Z_{DIGI} [cm];ToA_{DIGI} [TDC counts]", 50, -260., 260., 0., 1024.);
0287 }
0288
0289
0290 void BtlDigiHitsValidation::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0291 edm::ParameterSetDescription desc;
0292
0293 desc.add<std::string>("folder", "MTD/BTL/DigiHits");
0294 desc.add<edm::InputTag>("inputTag", edm::InputTag("mix", "FTLBarrel"));
0295 desc.add<bool>("optionalPlots", false);
0296
0297 descriptions.add("btlDigiHitsDefaultValid", desc);
0298 }
0299
0300 DEFINE_FWK_MODULE(BtlDigiHitsValidation);