Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
5 changes: 5 additions & 0 deletions DQM/HLTEvF/plugins/ScoutingCollectionMonitor.cc
Original file line number Diff line number Diff line change
Expand Up @@ -434,6 +434,7 @@ class ScoutingCollectionMonitor : public DQMEDAnalyzer {
dqm::reco::MonitorElement* eeMinusRecHitsXYMap[2];
dqm::reco::MonitorElement* hbheRecHitsNumber_hist;
dqm::reco::MonitorElement* hbheRecHits_energy_hist;
dqm::reco::MonitorElement* hbheRecHits_time_hist;
dqm::reco::MonitorElement* hbheRecHitsEtaPhiMap;
dqm::reco::MonitorElement* hbRecHitsEtaPhiMap;
dqm::reco::MonitorElement* heRecHitsEtaPhiMap;
Expand Down Expand Up @@ -970,6 +971,7 @@ void ScoutingCollectionMonitor::analyze(const edm::Event& iEvent, const edm::Eve
hbheRecHitsNumber_hist->Fill(hbheRecHitsH->size());
for (const auto& hbheRecHit : *hbheRecHitsH) {
hbheRecHits_energy_hist->Fill(hbheRecHit.energy());
hbheRecHits_time_hist->Fill(hbheRecHit.time());
HcalDetId hcalid(hbheRecHit.detId());
hbheRecHitsEtaPhiMap->Fill(hcalid.ieta(), hcalid.iphi());
const auto& subdet = hcalid.subdetId();
Expand Down Expand Up @@ -1478,6 +1480,9 @@ void ScoutingCollectionMonitor::bookHistograms(DQMStore::IBooker& ibook,
hbheRecHits_energy_hist = ibook.book1D(
"hbheRechits_energy", "Energy spectrum of hbhe RecHits; Energy of HBHE recHits (GeV); Entries", 100, 0.0, 200.0);

hbheRecHits_time_hist =
ibook.book1D("hbheRechits time", "Time of HBHE RecHits; Time of HBHE recHits (ns); Entries", 100, 0.0, 30.0);

hbheRecHitsEtaPhiMap = ibook.book2D(
"hbheRecHitsEtaPhitMap", "Occupancy map of HBHE rechits;ieta;iphi;Entries", 61, -30.5, 30.5, 74, -0.5, 73.5);
hbheRecHitsEtaPhiMap->setOption("colz");
Expand Down
56 changes: 41 additions & 15 deletions HLTriggerOffline/Scouting/plugins/ScoutingRecHitAnalyzer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -101,13 +101,20 @@ class ScoutingRecHitAnalyzer : public DQMEDAnalyzer {
// histograms
std::vector<MonitorElement*> number_histograms_;
std::vector<MonitorElement*> energy_histograms_;
std::vector<MonitorElement*> time_histograms_;
std::vector<MonitorElement*> ieta_iphi_histograms_;
std::vector<MonitorElement*> eta_phi_histograms_;
std::vector<MonitorElement*> ieta_iphi_energy_profiles_;
std::vector<MonitorElement*> eta_phi_energy_profiles_;
std::vector<MonitorElement*> ieta_iphi_time_profiles_;
std::vector<MonitorElement*> eta_phi_time_profiles_;

// trigger names
std::vector<std::string> trigger_names_;

// named constants
static constexpr double kDegToRad = M_PI / 180.0;
static constexpr double kFiveDegToRad = 5.0 * M_PI / 180.0;
};

template <typename RecHitType>
Expand Down Expand Up @@ -168,22 +175,21 @@ void ScoutingRecHitAnalyzer<RecHitType>::bookHistograms(DQMStore::IBooker& ibook
if constexpr (std::is_same<RecHitType, Run3ScoutingEBRecHit>()) {
number_histograms_.push_back(ibooker.book1D("number", "Number;Events", 100, 0., 1000.));
energy_histograms_.push_back(ibooker.book1D("energy", "Energy (GeV);Events", 100, 0., 20.));
time_histograms_.push_back(ibooker.book1D("time", "Time (ps);Events", 100, 0., 1000.));
ieta_iphi_histograms_.push_back(
ibooker.book2D("ieta_iphi", "i#eta;i#phi;Entries", 171, -85.5, 85.5, 360, 0.5, 360.5));
eta_phi_histograms_.push_back(ibooker.book2D(
"eta_phi", "#eta;#phi;Entries", 170, -85 * 0.017453292519943295, 85 * 0.017453292519943295, 360, -M_PI, M_PI));
eta_phi_histograms_.push_back(
ibooker.book2D("eta_phi", "#eta;#phi;Entries", 170, -85 * kDegToRad, 85 * kDegToRad, 360, -M_PI, M_PI));
ieta_iphi_energy_profiles_.push_back(ibooker.bookProfile2D(
"ieta_iphi_energy", "i#eta;i#phi;mean Energy", 171, -85.5, 85.5, 360, 0.5, 360.5, 0, 20));
eta_phi_energy_profiles_.push_back(ibooker.bookProfile2D("eta_phi_energy",
"#eta;#phi;mean Energy",
170,
-85 * 0.017453292519943295,
85 * 0.017453292519943295,
360,
-M_PI,
M_PI,
0,
20));
eta_phi_energy_profiles_.push_back(ibooker.bookProfile2D(
"eta_phi_energy", "#eta;#phi;mean Energy", 170, -85 * kDegToRad, 85 * kDegToRad, 360, -M_PI, M_PI, 0, 20));

ieta_iphi_time_profiles_.push_back(
ibooker.bookProfile2D("ieta_iphi_time", "i#eta;i#phi;mean Time", 171, -85.5, 85.5, 360, 0.5, 360.5, 0, 20));
eta_phi_time_profiles_.push_back(ibooker.bookProfile2D(
"eta_phi_time", "#eta;#phi;mean Time", 170, -85 * kDegToRad, 85 * kDegToRad, 360, -M_PI, M_PI, 0, 20));

} else if constexpr (std::is_same<RecHitType, Run3ScoutingHBHERecHit>()) {
std::vector<double> eta_edges(59);
eta_edges[29] = 0.0;
Expand All @@ -210,11 +216,15 @@ void ScoutingRecHitAnalyzer<RecHitType>::bookHistograms(DQMStore::IBooker& ibook

number_histograms_.push_back(ibooker.book1D("number", "Number;Events", 100, 0., 2000.));
energy_histograms_.push_back(ibooker.book1D("energy", "Energy (GeV);Events", 100, 0., 20.));
time_histograms_.push_back(ibooker.book1D("time", "Time (ns);Events", 100, 0., 30.));
ieta_iphi_histograms_.push_back(
ibooker.book2D("ieta_iphi", "i#eta;i#phi;Entries", 59, -29.5, 29.5, 72, 0.5, 72.5));
ieta_iphi_energy_profiles_.push_back(
ibooker.bookProfile2D("ieta_iphi_energy", "i#eta;i#phi;mean Energy", 59, -29.5, 29.5, 72, 0.5, 72.5, 0, 20));

ieta_iphi_time_profiles_.push_back(
ibooker.bookProfile2D("ieta_iphi_time", "i#eta;i#phi;mean Time", 59, -29.5, 29.5, 72, 0.5, 72.5, 0, 20));

// demote edges from double to float
std::vector<float> eta_edges_f(eta_edges.begin(), eta_edges.end());
std::vector<float> phi_edges_f(phi_edges.begin(), phi_edges.end());
Expand All @@ -231,6 +241,13 @@ void ScoutingRecHitAnalyzer<RecHitType>::bookHistograms(DQMStore::IBooker& ibook
eta_edges.data(),
phi_edges.size() - 1,
phi_edges.data()));

eta_phi_time_profiles_.push_back(ibooker.bookProfile2D("eta_phi_time",
"#eta;#phi;mean Time",
eta_edges.size() - 1,
eta_edges.data(),
phi_edges.size() - 1,
phi_edges.data()));
}
}
}
Expand Down Expand Up @@ -271,38 +288,47 @@ void ScoutingRecHitAnalyzer<RecHitType>::analyze(const edm::Event& iEvent, const
float eta, phi;
decodeDetId(rechit, ieta, eta, iphi, phi);
const auto& energy = rechit.energy();
const auto& time = rechit.time();

if (cut_(rechit)) {
energy_histograms_[0]->Fill(energy);
time_histograms_[0]->Fill(time);
ieta_iphi_histograms_[0]->Fill(ieta, iphi);
eta_phi_histograms_[0]->Fill(eta, phi);
ieta_iphi_energy_profiles_[0]->Fill(ieta, iphi, energy);
eta_phi_energy_profiles_[0]->Fill(eta, phi, energy);
ieta_iphi_time_profiles_[0]->Fill(ieta, iphi, time);
eta_phi_time_profiles_[0]->Fill(eta, phi, time);
if (std::is_same<RecHitType, Run3ScoutingHBHERecHit>()) {
if (std::abs(ieta) >= 21) {
float phi2 = (iphi + 0.5) * 5 * M_PI / 180.;
if (phi2 > M_PI)
phi2 = phi2 - 2 * M_PI;
eta_phi_histograms_[0]->Fill(eta, phi2);
eta_phi_energy_profiles_[0]->Fill(eta, phi2, energy);
eta_phi_time_profiles_[0]->Fill(eta, phi2, time);
}
}
unsigned int itrigger = 0;
for (const bool& trigger_accept : trigger_decisions) {
itrigger++;
if (trigger_accept) {
energy_histograms_[itrigger]->Fill(energy);
time_histograms_[itrigger]->Fill(time);
ieta_iphi_histograms_[itrigger]->Fill(ieta, iphi);
eta_phi_histograms_[itrigger]->Fill(eta, phi);
ieta_iphi_energy_profiles_[itrigger]->Fill(ieta, iphi, energy);
eta_phi_energy_profiles_[itrigger]->Fill(eta, phi, energy);
ieta_iphi_time_profiles_[itrigger]->Fill(ieta, iphi, time);
eta_phi_time_profiles_[itrigger]->Fill(eta, phi, time);
if (std::is_same<RecHitType, Run3ScoutingHBHERecHit>()) {
if (std::abs(ieta) >= 21) {
float phi2 = (iphi + 0.5) * 5 * M_PI / 180.;
if (phi2 > M_PI)
phi2 = phi2 - 2 * M_PI;
eta_phi_histograms_[itrigger]->Fill(eta, phi2);
eta_phi_energy_profiles_[itrigger]->Fill(eta, phi2, energy);
eta_phi_time_profiles_[itrigger]->Fill(eta, phi2, time);
}
}
}
Expand All @@ -329,7 +355,7 @@ void ScoutingRecHitAnalyzer<RecHitType>::decodeDetId(
ieta = detId.ieta();
const auto& ietaAbs = detId.ietaAbs();
const auto& zside = detId.zside();
eta = zside * (ietaAbs - 0.5) * 0.017453292519943295;
eta = zside * (ietaAbs - 0.5) * kDegToRad;
iphi = detId.iphi();
phi = (iphi - 0.5) * M_PI / 180.;
if (phi > M_PI)
Expand All @@ -340,9 +366,9 @@ void ScoutingRecHitAnalyzer<RecHitType>::decodeDetId(
const auto& ietaAbs = detId.ietaAbs();
const auto& zside = detId.zside();
if (ietaAbs <= 20) {
eta = zside * (ietaAbs - 0.5) * 0.087266462599716475;
eta = zside * (ietaAbs - 0.5) * kFiveDegToRad;
} else {
eta = zside * ((20 * 0.087266462599716475) + (ietaAbs - 20 - 0.5) * 0.087266462599716475 * 2);
eta = zside * ((20 * kFiveDegToRad) + (ietaAbs - 20 - 0.5) * kFiveDegToRad * 2);
}
iphi = detId.iphi();
phi = (iphi - 0.5) * 5 * M_PI / 180.;
Expand Down