File indexing completed on 2025-06-03 00:12:44
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/Math/interface/GeantUnits.h"
0026 #include "DataFormats/ForwardDetId/interface/BTLDetId.h"
0027
0028 #include "SimDataFormats/CrossingFrame/interface/CrossingFrame.h"
0029 #include "SimDataFormats/CrossingFrame/interface/MixCollection.h"
0030 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0031
0032 #include "Geometry/Records/interface/MTDDigiGeometryRecord.h"
0033 #include "Geometry/Records/interface/MTDTopologyRcd.h"
0034 #include "Geometry/MTDGeometryBuilder/interface/MTDGeometry.h"
0035 #include "Geometry/MTDGeometryBuilder/interface/MTDTopology.h"
0036
0037 #include "Geometry/MTDGeometryBuilder/interface/ProxyMTDTopology.h"
0038 #include "Geometry/MTDGeometryBuilder/interface/RectangularMTDTopology.h"
0039
0040 #include "Geometry/MTDCommonData/interface/MTDTopologyMode.h"
0041
0042 #include "MTDHit.h"
0043
0044 class BtlSimHitsValidation : public DQMEDAnalyzer {
0045 public:
0046 explicit BtlSimHitsValidation(const edm::ParameterSet&);
0047 ~BtlSimHitsValidation() override;
0048
0049 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0050
0051 private:
0052 void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0053
0054 void analyze(const edm::Event&, const edm::EventSetup&) override;
0055
0056
0057
0058 const std::string folder_;
0059 const bool optionalPlots_;
0060 const float hitMinEnergy_;
0061
0062 static constexpr float BXTime_ = 25.;
0063
0064 edm::EDGetTokenT<CrossingFrame<PSimHit>> btlSimHitsToken_;
0065
0066 edm::ESGetToken<MTDGeometry, MTDDigiGeometryRecord> mtdgeoToken_;
0067 edm::ESGetToken<MTDTopology, MTDTopologyRcd> mtdtopoToken_;
0068
0069
0070
0071 MonitorElement* meNevents_;
0072
0073 MonitorElement* meNhits_;
0074 MonitorElement* meNtrkPerCell_;
0075
0076 MonitorElement* meHitEnergy_;
0077
0078 static constexpr int nRU_ = BTLDetId::kRUPerRod;
0079 static constexpr double logE_min = -2;
0080 static constexpr double logE_max = 2;
0081 static constexpr double n_bin_logE = 40;
0082
0083 MonitorElement* meHitLogEnergy_;
0084 MonitorElement* meHitLogEnergyRUSlice_[nRU_];
0085 MonitorElement* meHitMultCell_;
0086 MonitorElement* meHitMultCellRUSlice_[nRU_];
0087 MonitorElement* meHitMultSM_;
0088 MonitorElement* meHitMultSMRUSlice_[nRU_];
0089 MonitorElement* meHitMultCellperSMvsEth_;
0090 MonitorElement* meHitMultCellperSMvsEthRUSlice_[nRU_];
0091
0092 MonitorElement* meHitTime_;
0093
0094 MonitorElement* meHitXlocal_;
0095 MonitorElement* meHitYlocal_;
0096 MonitorElement* meHitZlocal_;
0097
0098 MonitorElement* meOccupancy_;
0099
0100 MonitorElement* meHitX_;
0101 MonitorElement* meHitY_;
0102 MonitorElement* meHitZ_;
0103 MonitorElement* meHitPhi_;
0104 MonitorElement* meHitEta_;
0105
0106 MonitorElement* meHitTvsE_;
0107 MonitorElement* meHitEvsPhi_;
0108 MonitorElement* meHitEvsEta_;
0109 MonitorElement* meHitEvsZ_;
0110 MonitorElement* meHitTvsPhi_;
0111 MonitorElement* meHitTvsEta_;
0112 MonitorElement* meHitTvsZ_;
0113
0114 MonitorElement* meHitTrkID1_;
0115 MonitorElement* meHitTrkID2_;
0116 MonitorElement* meHitTrkID3_;
0117 MonitorElement* meHitTrkID4_;
0118
0119 MonitorElement* meHitDeltaT2_;
0120 MonitorElement* meHitDeltaT3_;
0121 MonitorElement* meHitDeltaT4_;
0122
0123 MonitorElement* meHitDeltaE12vsE1_;
0124
0125 static constexpr float cellEneCut_ = 10.;
0126
0127 MonitorElement* meHitE1overEcellBulk1_;
0128 MonitorElement* meHitE1overEcellBulk2_;
0129 MonitorElement* meHitE1overEcellBulk3_;
0130 MonitorElement* meHitE1overEcellBulk4_;
0131 MonitorElement* meHitE1overEcellTail1_;
0132 MonitorElement* meHitE1overEcellTail2_;
0133 MonitorElement* meHitE1overEcellTail3_;
0134 MonitorElement* meHitE1overEcellTail4_;
0135 };
0136
0137
0138 BtlSimHitsValidation::BtlSimHitsValidation(const edm::ParameterSet& iConfig)
0139 : folder_(iConfig.getParameter<std::string>("folder")),
0140 optionalPlots_(iConfig.getParameter<bool>("optionalPlots")),
0141 hitMinEnergy_(iConfig.getParameter<double>("hitMinimumEnergy")) {
0142 btlSimHitsToken_ = consumes<CrossingFrame<PSimHit>>(iConfig.getParameter<edm::InputTag>("inputTag"));
0143 mtdgeoToken_ = esConsumes<MTDGeometry, MTDDigiGeometryRecord>();
0144 mtdtopoToken_ = esConsumes<MTDTopology, MTDTopologyRcd>();
0145 }
0146
0147 BtlSimHitsValidation::~BtlSimHitsValidation() {}
0148
0149
0150 void BtlSimHitsValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0151 using namespace edm;
0152 using namespace geant_units::operators;
0153
0154 auto geometryHandle = iSetup.getTransientHandle(mtdgeoToken_);
0155 const MTDGeometry* geom = geometryHandle.product();
0156
0157 auto topologyHandle = iSetup.getTransientHandle(mtdtopoToken_);
0158 const MTDTopology* topology = topologyHandle.product();
0159
0160 auto btlSimHitsHandle = makeValid(iEvent.getHandle(btlSimHitsToken_));
0161 MixCollection<PSimHit> btlSimHits(btlSimHitsHandle.product());
0162
0163
0164
0165
0166
0167
0168
0169
0170
0171 std::unordered_map<uint32_t, std::unordered_map<uint32_t, std::unordered_map<uint64_t, MTDHit>>>
0172 m_btlHitsPerCellAndModule;
0173
0174
0175 for (auto const& simHit : btlSimHits) {
0176
0177 if (simHit.tof() < 0. || simHit.tof() > BXTime_) {
0178 continue;
0179 }
0180
0181 DetId id = simHit.detUnitId();
0182
0183
0184 uint64_t globalTrkID = ((uint64_t)simHit.eventId().rawId() << 32) | simHit.trackId();
0185
0186
0187 BTLDetId detId(id.rawId());
0188 DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
0189
0190
0191 m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].energy += convertGeVToMeV(simHit.energyLoss());
0192
0193
0194 if (m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].time == 0. ||
0195 simHit.tof() < m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].time) {
0196 m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].time = simHit.tof();
0197
0198 auto hit_pos = simHit.localPosition();
0199 m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].x = hit_pos.x();
0200 m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].y = hit_pos.y();
0201 m_btlHitsPerCellAndModule[geoId.rawId()][id.rawId()][globalTrkID].z = hit_pos.z();
0202 }
0203
0204 }
0205
0206
0207
0208
0209
0210 int Nhits = 0;
0211 for (const auto& module : m_btlHitsPerCellAndModule) {
0212 Nhits += module.second.size();
0213 }
0214 if (Nhits > 0)
0215 meNhits_->Fill(log10(Nhits));
0216
0217
0218 for (const auto& module : m_btlHitsPerCellAndModule) {
0219
0220 for (auto const& cell : module.second) {
0221
0222 const auto& m_hits = cell.second;
0223
0224
0225 if (m_hits.empty()) {
0226 continue;
0227 }
0228
0229
0230 std::vector<uint64_t> v_hitID;
0231 float ene_tot_cell = 0.;
0232
0233 for (auto const& hit : m_hits) {
0234 ene_tot_cell += hit.second.energy;
0235 v_hitID.push_back(hit.first);
0236 }
0237
0238 meHitLogEnergy_->Fill(log10(ene_tot_cell));
0239
0240
0241 BTLDetId detId(cell.first);
0242 DetId geoId = detId.geographicalId(MTDTopologyMode::crysLayoutFromTopoMode(topology->getMTDTopologyMode()));
0243
0244 if (optionalPlots_)
0245 meHitLogEnergyRUSlice_[detId.runit() - 1]->Fill(log10(ene_tot_cell));
0246
0247
0248 if (ene_tot_cell < hitMinEnergy_) {
0249 continue;
0250 }
0251
0252
0253 bool swapped;
0254 for (unsigned int ihit = 0; ihit < v_hitID.size() - 1; ++ihit) {
0255 swapped = false;
0256 for (unsigned int jhit = 0; jhit < v_hitID.size() - ihit - 1; ++jhit) {
0257 if (m_hits.at(v_hitID[jhit]).time > m_hits.at(v_hitID[jhit + 1]).time) {
0258 std::swap(v_hitID[jhit], v_hitID[jhit + 1]);
0259 swapped = true;
0260 }
0261 }
0262 if (swapped == false) {
0263 break;
0264 }
0265 }
0266
0267
0268 float deltaT_max = 0.;
0269 for (unsigned int ihit = 0; ihit < v_hitID.size() - 1; ++ihit) {
0270 float deltaT = m_hits.at(v_hitID[ihit + 1]).time - m_hits.at(v_hitID[ihit]).time;
0271
0272 if (deltaT > deltaT_max) {
0273 deltaT_max = deltaT;
0274 }
0275 }
0276
0277
0278 const MTDGeomDet* thedet = geom->idToDet(geoId);
0279 if (thedet == nullptr) {
0280 throw cms::Exception("BtlSimHitsValidation") << "GeographicalID: " << std::hex << geoId.rawId() << " ("
0281 << detId.rawId() << ") is invalid!" << std::dec << std::endl;
0282 }
0283 const ProxyMTDTopology& topoproxy = static_cast<const ProxyMTDTopology&>(thedet->topology());
0284 const RectangularMTDTopology& topo = static_cast<const RectangularMTDTopology&>(topoproxy.specificTopology());
0285
0286 Local3DPoint local_point(convertMmToCm(m_hits.at(v_hitID[0]).x),
0287 convertMmToCm(m_hits.at(v_hitID[0]).y),
0288 convertMmToCm(m_hits.at(v_hitID[0]).z));
0289
0290 local_point = topo.pixelToModuleLocalPoint(local_point, detId.row(topo.nrows()), detId.column(topo.nrows()));
0291 const auto& global_point = thedet->toGlobal(local_point);
0292
0293
0294 meHitEnergy_->Fill(ene_tot_cell);
0295 meHitTime_->Fill(m_hits.at(v_hitID[0]).time);
0296
0297 meHitXlocal_->Fill(m_hits.at(v_hitID[0]).x);
0298 meHitYlocal_->Fill(m_hits.at(v_hitID[0]).y);
0299 meHitZlocal_->Fill(m_hits.at(v_hitID[0]).z);
0300
0301 meOccupancy_->Fill(global_point.z(), global_point.phi());
0302
0303 meHitX_->Fill(global_point.x());
0304 meHitY_->Fill(global_point.y());
0305 meHitZ_->Fill(global_point.z());
0306 meHitPhi_->Fill(global_point.phi());
0307 meHitEta_->Fill(global_point.eta());
0308
0309 meHitTvsE_->Fill(ene_tot_cell, m_hits.at(v_hitID[0]).time);
0310 meHitEvsPhi_->Fill(global_point.phi(), ene_tot_cell);
0311 meHitEvsEta_->Fill(global_point.eta(), ene_tot_cell);
0312 meHitEvsZ_->Fill(global_point.z(), ene_tot_cell);
0313 meHitTvsPhi_->Fill(global_point.phi(), m_hits.at(v_hitID[0]).time);
0314 meHitTvsEta_->Fill(global_point.eta(), m_hits.at(v_hitID[0]).time);
0315 meHitTvsZ_->Fill(global_point.z(), m_hits.at(v_hitID[0]).time);
0316
0317 meNtrkPerCell_->Fill(v_hitID.size());
0318
0319
0320 int trackID = (int)(v_hitID[0] & 0xFFFFFFFF) / 100000000;
0321 meHitTrkID1_->Fill(trackID);
0322
0323
0324 if (ene_tot_cell < cellEneCut_) {
0325 meHitE1overEcellBulk1_->Fill(m_hits.at(v_hitID[0]).energy / ene_tot_cell);
0326 }
0327
0328 else {
0329 meHitE1overEcellTail1_->Fill(m_hits.at(v_hitID[0]).energy / ene_tot_cell);
0330 }
0331
0332
0333 if (v_hitID.size() == 2) {
0334 trackID = (int)(v_hitID[1] & 0xFFFFFFFF) / 100000000;
0335 meHitTrkID2_->Fill(trackID);
0336
0337 meHitDeltaT2_->Fill(deltaT_max);
0338
0339
0340 if (ene_tot_cell < cellEneCut_) {
0341 meHitE1overEcellBulk2_->Fill(m_hits.at(v_hitID[1]).energy / ene_tot_cell);
0342 }
0343
0344 else {
0345 meHitE1overEcellTail2_->Fill(m_hits.at(v_hitID[1]).energy / ene_tot_cell);
0346 }
0347
0348 meHitDeltaE12vsE1_->Fill(m_hits.at(v_hitID[0]).energy,
0349 m_hits.at(v_hitID[0]).energy - m_hits.at(v_hitID[1]).energy);
0350
0351 }
0352
0353 else if (v_hitID.size() == 3) {
0354 trackID = (int)(v_hitID[2] & 0xFFFFFFFF) / 100000000;
0355 meHitTrkID3_->Fill(trackID);
0356
0357 meHitDeltaT3_->Fill(deltaT_max);
0358
0359
0360 if (ene_tot_cell < cellEneCut_) {
0361 meHitE1overEcellBulk3_->Fill(m_hits.at(v_hitID[2]).energy / ene_tot_cell);
0362 }
0363
0364 else {
0365 meHitE1overEcellTail3_->Fill(m_hits.at(v_hitID[2]).energy / ene_tot_cell);
0366 }
0367
0368 }
0369
0370 else if (v_hitID.size() >= 4) {
0371 for (unsigned int ihit = 3; ihit < v_hitID.size(); ++ihit) {
0372 trackID = (int)(v_hitID[ihit] & 0xFFFFFFFF) / 100000000;
0373 meHitTrkID4_->Fill(trackID);
0374
0375
0376 if (ene_tot_cell < cellEneCut_) {
0377 meHitE1overEcellBulk4_->Fill(m_hits.at(v_hitID[ihit]).energy / ene_tot_cell);
0378 }
0379
0380 else {
0381 meHitE1overEcellTail4_->Fill(m_hits.at(v_hitID[ihit]).energy / ene_tot_cell);
0382 }
0383 }
0384
0385 meHitDeltaT4_->Fill(deltaT_max);
0386 }
0387
0388 }
0389
0390 }
0391
0392
0393
0394 double bin_w_logE = (logE_max - logE_min) / n_bin_logE;
0395 for (int i = 0; i < n_bin_logE; i++) {
0396 double th_logE = logE_min + i * bin_w_logE;
0397
0398 for (const auto& module : m_btlHitsPerCellAndModule) {
0399 bool SM_bool = false;
0400 int crystal_count = 0;
0401 int SM_globalRunit = -1;
0402 for (const auto& cell : module.second) {
0403 float ene_tot_cell = 0.;
0404 for (const auto& hit : cell.second) {
0405 ene_tot_cell += hit.second.energy;
0406 }
0407 BTLDetId detId(cell.first);
0408 SM_globalRunit = detId.runit() - 1;
0409 if (log10(ene_tot_cell) > th_logE) {
0410 SM_bool = true;
0411 crystal_count++;
0412 meHitMultCell_->Fill(th_logE + bin_w_logE / 2);
0413 if (optionalPlots_)
0414 meHitMultCellRUSlice_[SM_globalRunit]->Fill(th_logE + bin_w_logE / 2);
0415 }
0416 }
0417 meHitMultCellperSMvsEth_->Fill(th_logE + bin_w_logE / 2, crystal_count);
0418 if (optionalPlots_)
0419 meHitMultCellperSMvsEthRUSlice_[SM_globalRunit]->Fill(th_logE + bin_w_logE / 2, crystal_count);
0420
0421 if (SM_bool && SM_globalRunit != -1) {
0422 meHitMultSM_->Fill(th_logE + bin_w_logE / 2);
0423 if (optionalPlots_)
0424 meHitMultSMRUSlice_[SM_globalRunit]->Fill(th_logE + bin_w_logE / 2);
0425 }
0426 }
0427 }
0428
0429
0430 meNevents_->Fill(0.5);
0431 }
0432
0433
0434 void BtlSimHitsValidation::bookHistograms(DQMStore::IBooker& ibook,
0435 edm::Run const& run,
0436 edm::EventSetup const& iSetup) {
0437 ibook.setCurrentFolder(folder_);
0438
0439
0440
0441 meNevents_ = ibook.book1D("BtlNevents", "Number of events", 1, 0., 1.);
0442 meNhits_ = ibook.book1D("BtlNhits", "Number of BTL cells with SIM hits;log_{10}(N_{BTL cells})", 100, 0., 5.25);
0443 meNtrkPerCell_ = ibook.book1D("BtlNtrkPerCell", "Number of tracks per BTL cell;N_{trk}", 10, 0., 10.);
0444
0445 meHitEnergy_ = ibook.book1D("BtlHitEnergy", "BTL SIM hits energy;E_{SIM} [MeV]", 100, 0., 20.);
0446 meHitLogEnergy_ = ibook.book1D("BtlHitLogEnergy", "BTL SIM hits energy;log_{10}(E_{SIM} [MeV])", 200, -6., 3.);
0447 meHitMultCell_ = ibook.book1D("BtlHitMultCell",
0448 "BTL cell multiplicity vs energy threshold ;log_{10}(E_{th} [MeV])",
0449 n_bin_logE,
0450 logE_min,
0451 logE_max);
0452 meHitMultSM_ = ibook.book1D(
0453 "BtlHitMultSM", "BTL SM multiplicity vs energy threshold; log_{10}(E_{th} [MeV])", n_bin_logE, logE_min, logE_max);
0454 meHitMultCellperSMvsEth_ =
0455 ibook.bookProfile("BtlHitMultCellperSMvsEth",
0456 "BTL cell per SM vs energy threshold;log_{10}(E_{th} [MeV]); cells per SM",
0457 n_bin_logE,
0458 logE_min,
0459 logE_max,
0460 0,
0461 16);
0462
0463 if (optionalPlots_) {
0464 for (unsigned int ihistoRU = 0; ihistoRU < nRU_; ++ihistoRU) {
0465 std::string name_LogEnergy = "BtlHitLogEnergyRUSlice" + std::to_string(ihistoRU);
0466 std::string title_LogEnergy = "BTL SIM hits energy (RU " + std::to_string(ihistoRU) + ");log_{10}(E_{SIM} [MeV])";
0467 meHitLogEnergyRUSlice_[ihistoRU] = ibook.book1D(name_LogEnergy, title_LogEnergy, 200, -6., 3.);
0468 std::string name_Cell = "BtlHitMultCellRUSlice" + std::to_string(ihistoRU);
0469 std::string title_Cell =
0470 "BTL cell multiplicity vs energy threshold (RU " + std::to_string(ihistoRU) + ");log_{10}(E_{th} [MeV])";
0471 meHitMultCellRUSlice_[ihistoRU] = ibook.book1D(name_Cell, title_Cell, n_bin_logE, logE_min, logE_max);
0472 std::string name_SM = "BtlHitMultSMRUSlice" + std::to_string(ihistoRU);
0473 std::string title_SM =
0474 "BTL SM multiplicity vs energy threshold (RU " + std::to_string(ihistoRU) + ");log_{10}(E_{th} [MeV])";
0475 meHitMultSMRUSlice_[ihistoRU] = ibook.book1D(name_SM, title_SM, n_bin_logE, logE_min, logE_max);
0476 std::string name_cellperSM = "BtlHitMultCellperSMvsEthRUSlice" + std::to_string(ihistoRU);
0477 std::string title_cellperSM = "BTL cell per SM vs energy threshold (RU " + std::to_string(ihistoRU) +
0478 ");log_{10}(E_{th} [MeV]);cells per SM";
0479 meHitMultCellperSMvsEthRUSlice_[ihistoRU] =
0480 ibook.bookProfile(name_cellperSM, title_cellperSM, n_bin_logE, logE_min, logE_max, 0, 16);
0481 }
0482 }
0483
0484 meHitTime_ = ibook.book1D("BtlHitTime", "BTL SIM hits ToA;ToA_{SIM} [ns]", 100, 0., 25.);
0485
0486 meHitXlocal_ = ibook.book1D("BtlHitXlocal", "BTL SIM local X;X_{SIM}^{LOC} [mm]", 100, -30., 30.);
0487 meHitYlocal_ = ibook.book1D("BtlHitYlocal", "BTL SIM local Y;Y_{SIM}^{LOC} [mm]", 100, -1.65, 1.65);
0488 meHitZlocal_ = ibook.book1D("BtlHitZlocal", "BTL SIM local Z;Z_{SIM}^{LOC} [mm]", 100, -2., 2.);
0489
0490 meOccupancy_ = ibook.book2D(
0491 "BtlOccupancy", "BTL SIM hits occupancy;z_{SIM} [cm];#phi_{SIM} [rad]", 130, -260., 260., 200, -3.15, 3.15);
0492
0493 meHitX_ = ibook.book1D("BtlHitX", "BTL SIM hits X;X_{SIM} [cm]", 100, -120., 120.);
0494 meHitY_ = ibook.book1D("BtlHitY", "BTL SIM hits Y;Y_{SIM} [cm]", 100, -120., 120.);
0495 meHitZ_ = ibook.book1D("BtlHitZ", "BTL SIM hits Z;Z_{SIM} [cm]", 100, -260., 260.);
0496 meHitPhi_ = ibook.book1D("BtlHitPhi", "BTL SIM hits #phi;#phi_{SIM} [rad]", 200, -3.15, 3.15);
0497 meHitEta_ = ibook.book1D("BtlHitEta", "BTL SIM hits #eta;#eta_{SIM}", 100, -1.55, 1.55);
0498
0499 meHitTvsE_ =
0500 ibook.bookProfile("BtlHitTvsE", "BTL SIM time vs energy;E_{SIM} [MeV];T_{SIM} [ns]", 50, 0., 20., 0., 100.);
0501 meHitEvsPhi_ = ibook.bookProfile(
0502 "BtlHitEvsPhi", "BTL SIM energy vs #phi;#phi_{SIM} [rad];E_{SIM} [MeV]", 50, -3.15, 3.15, 0., 100.);
0503 meHitEvsEta_ =
0504 ibook.bookProfile("BtlHitEvsEta", "BTL SIM energy vs #eta;#eta_{SIM};E_{SIM} [MeV]", 50, -1.55, 1.55, 0., 100.);
0505 meHitEvsZ_ =
0506 ibook.bookProfile("BtlHitEvsZ", "BTL SIM energy vs Z;Z_{SIM} [cm];E_{SIM} [MeV]", 50, -260., 260., 0., 100.);
0507 meHitTvsPhi_ = ibook.bookProfile(
0508 "BtlHitTvsPhi", "BTL SIM time vs #phi;#phi_{SIM} [rad];T_{SIM} [ns]", 50, -3.15, 3.15, 0., 100.);
0509 meHitTvsEta_ =
0510 ibook.bookProfile("BtlHitTvsEta", "BTL SIM time vs #eta;#eta_{SIM};T_{SIM} [ns]", 50, -1.55, 1.55, 0., 100.);
0511 meHitTvsZ_ =
0512 ibook.bookProfile("BtlHitTvsZ", "BTL SIM time vs Z;Z_{SIM} [cm];T_{SIM} [ns]", 50, -260., 260., 0., 100.);
0513
0514 meHitTrkID1_ = ibook.book1I("BtlHitTrkID1", "Category of the 1^{st} hit in time;Hit category", 8, 0, 8);
0515 meHitTrkID2_ = ibook.book1I("BtlHitTrkID2", "Category of the 2^{nd} hit in time;Hit category", 8, 0, 8);
0516 meHitTrkID3_ = ibook.book1I("BtlHitTrkID3", "Category of the 3^{rd} hit in time;Hit category", 8, 0, 8);
0517 meHitTrkID4_ = ibook.book1I("BtlHitTrkID4", "Category of the #geq4^{th} hit in time;Hit category", 8, 0, 8);
0518
0519 meHitDeltaT2_ = ibook.book1D(
0520 "BtlHitDeltaT2", "Time interval between hits in the same cell (2 hits);#DeltaT_{2} [ns]", 100., 0., 25.);
0521 meHitDeltaT3_ = ibook.book1D(
0522 "BtlHitDeltaT3", "Max time interval between hits in the same cell (3 hits);#DeltaT_{3} [ns]", 100., 0., 25.);
0523 meHitDeltaT4_ = ibook.book1D("BtlHitDeltaT4",
0524 "Max time interval between hits in the same cell (#geq4 hits);#DeltaT_{#geq4} [ns]",
0525 100.,
0526 0.,
0527 25.);
0528
0529 meHitDeltaE12vsE1_ = ibook.bookProfile(
0530 "BtlHitDeltaE12vsE", "E_{1}-E_{2} vs E_{1};E_{1} [MeV];E_{1}-E_{2} [MeV]", 100, 0., 100., 0., 100.);
0531
0532 meHitE1overEcellBulk1_ = ibook.book1D("BtlHitE1overEtotBulk1",
0533 "Energy fraction of the 1^{st} hit in time (E_{cell}<10 MeV);E_{1} / E_{cell}",
0534 100,
0535 0.,
0536 1.);
0537 meHitE1overEcellBulk2_ = ibook.book1D("BtlHitE1overEtotBulk2",
0538 "Energy fraction of the 2^{nd} hit in time (E_{cell}<10 MeV);E_{2} / E_{cell}",
0539 100,
0540 0.,
0541 1.);
0542 meHitE1overEcellBulk3_ = ibook.book1D("BtlHitE1overEtotBulk3",
0543 "Energy fraction of the 3^{rd} hit in time (E_{cell}<10 MeV);E_{3} / E_{cell}",
0544 100,
0545 0.,
0546 1.);
0547 meHitE1overEcellBulk4_ =
0548 ibook.book1D("BtlHitE1overEtotBulk4",
0549 "Energy fraction of the #geq4^{th} hit in time (E_{cell}<10 MeV);E_{#geq4} / E_{cell}",
0550 100,
0551 0.,
0552 1.);
0553 meHitE1overEcellTail1_ = ibook.book1D("BtlHitE1overEtotTail1",
0554 "Energy fraction of the 1^{st} hit in time (E_{cell}>10 MeV);E_{1} / E_{cell}",
0555 100,
0556 0.,
0557 1.);
0558 meHitE1overEcellTail2_ = ibook.book1D("BtlHitE1overEtotTail2",
0559 "Energy fraction of the 2^{nd} hit in time (E_{cell}>10 MeV);E_{2} / E_{cell}",
0560 100,
0561 0.,
0562 1.);
0563 meHitE1overEcellTail3_ = ibook.book1D("BtlHitE1overEtotTail3",
0564 "Energy fraction of the 3^{rd} hit in time (E_{cell}>10 MeV);E_{3} / E_{cell}",
0565 100,
0566 0.,
0567 1.);
0568 meHitE1overEcellTail4_ =
0569 ibook.book1D("BtlHitE1overEtotTail4",
0570 "Energy fraction of the #geq4^{th} hit in time (E_{cell}>10 MeV);E_{#geq4} / E_{cell}",
0571 100,
0572 0.,
0573 1.);
0574 }
0575
0576
0577 void BtlSimHitsValidation::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0578 edm::ParameterSetDescription desc;
0579
0580 desc.add<std::string>("folder", "MTD/BTL/SimHits");
0581 desc.add<bool>("optionalPlots", false);
0582 desc.add<edm::InputTag>("inputTag", edm::InputTag("mix", "g4SimHitsFastTimerHitsBarrel"));
0583 desc.add<double>("hitMinimumEnergy", 1.);
0584
0585 descriptions.add("btlSimHitsValid", desc);
0586 }
0587
0588 DEFINE_FWK_MODULE(BtlSimHitsValidation);