Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add HGCAL Hit calibration to Phase2 HLT. #44684

Merged
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.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
43 changes: 28 additions & 15 deletions HLTriggerOffline/Common/python/HLTValidation_cff.py
Original file line number Diff line number Diff line change
Expand Up @@ -18,6 +18,17 @@
from Validation.HcalDigis.HLTHcalDigisParam_cfi import *
from Validation.HcalRecHits.HLTHcalRecHitParam_cfi import *

# HGCAL Rechit Calibration
from Validation.HGCalValidation.hgcalHitCalibrationDefault_cfi import hgcalHitCalibrationDefault as _hgcalHitCalibrationDefault
hgcalHitCalibrationHLT = _hgcalHitCalibrationDefault.clone()
hgcalHitCalibrationHLT.folder = "HGCalHitCalibrationHLT"
hgcalHitCalibrationHLT.recHitsEE = cms.InputTag("HGCalRecHit", "HGCEERecHits", "HLT")
hgcalHitCalibrationHLT.recHitsFH = cms.InputTag("HGCalRecHit", "HGCHEFRecHits", "HLT")
hgcalHitCalibrationHLT.recHitsBH = cms.InputTag("HGCalRecHit", "HGCHEBRecHits", "HLT")
hgcalHitCalibrationHLT.hgcalMultiClusters = cms.InputTag("None")
hgcalHitCalibrationHLT.electrons = cms.InputTag("None")
hgcalHitCalibrationHLT.photons = cms.InputTag("None")

# offline dqm:
# from DQMOffline.Trigger.DQMOffline_Trigger_cff.py import *
from DQMOffline.Trigger.HLTTauDQMOffline_cff import *
Expand Down Expand Up @@ -70,22 +81,24 @@
)

# Temporary Phase-2 config
# Exclude everything except Muon and JetMET for now
# Exclude everything except Muon and JetMET for now. Add HGCAL Hit Calibration
from Configuration.Eras.Modifier_phase2_common_cff import phase2_common
phase2_common.toReplaceWith(hltvalidationWithMC, hltvalidationWithMC.copyAndExclude([#HLTMuonVal,
HLTTauVal,
egammaValidationSequence,
heavyFlavorValidationSequence,
#HLTJetMETValSeq,
HLTSusyExoValSeq,
HiggsValidationSequence,
ExoticaValidationSequence,
b2gHLTriggerValidation,
SMPValidationSequence,
hltbtagValidationSequence,
hltHCALdigisAnalyzer,
hltHCALRecoAnalyzer,
hltHCALNoiseRates]))
_hltvalidationWithMC_Phase2 = hltvalidationWithMC.copyAndExclude([#HLTMuonVal,
HLTTauVal,
egammaValidationSequence,
heavyFlavorValidationSequence,
#HLTJetMETValSeq,
HLTSusyExoValSeq,
HiggsValidationSequence,
ExoticaValidationSequence,
b2gHLTriggerValidation,
SMPValidationSequence,
hltbtagValidationSequence,
hltHCALdigisAnalyzer,
hltHCALRecoAnalyzer,
hltHCALNoiseRates])
_hltvalidationWithMC_Phase2.insert(-1, hgcalHitCalibrationHLT)
phase2_common.toReplaceWith(hltvalidationWithMC, _hltvalidationWithMC_Phase2)

hltvalidationWithData = cms.Sequence(
)
Expand Down
96 changes: 64 additions & 32 deletions Validation/HGCalValidation/plugins/HGCalHitCalibration.cc
Original file line number Diff line number Diff line change
Expand Up @@ -34,6 +34,7 @@
#include <array>
#include <string>
#include <numeric>
#include <cmath>

class HGCalHitCalibration : public DQMEDAnalyzer {
public:
Expand All @@ -60,6 +61,7 @@ class HGCalHitCalibration : public DQMEDAnalyzer {
edm::ESGetToken<CaloGeometry, CaloGeometryRecord> caloGeomToken_;
bool rawRecHits_;
int debug_;
std::string folder_;
hgcal::RecHitTools recHitTools_;
static constexpr int depletion1_ = 200;
static constexpr int depletion2_ = 300;
Expand All @@ -79,7 +81,8 @@ class HGCalHitCalibration : public DQMEDAnalyzer {
HGCalHitCalibration::HGCalHitCalibration(const edm::ParameterSet& iConfig)
: caloGeomToken_(esConsumes<CaloGeometry, CaloGeometryRecord>()),
rawRecHits_(iConfig.getParameter<bool>("rawRecHits")),
debug_(iConfig.getParameter<int>("debug")) {
debug_(iConfig.getParameter<int>("debug")),
folder_(iConfig.getParameter<std::string>("folder")) {
auto detector = iConfig.getParameter<std::string>("detector");
auto recHitsEE = iConfig.getParameter<edm::InputTag>("recHitsEE");
auto recHitsFH = iConfig.getParameter<edm::InputTag>("recHitsFH");
Expand Down Expand Up @@ -121,7 +124,7 @@ void HGCalHitCalibration::bookHistograms(DQMStore::IBooker& ibooker,
edm::Run const& iRun,
edm::EventSetup const& /* iSetup */) {
ibooker.cd();
ibooker.setCurrentFolder("HGCalHitCalibration");
ibooker.setCurrentFolder(folder_);
h_EoP_CPene_calib_fraction_[depletion0_] = ibooker.book1D("h_EoP_CPene_100_calib_fraction", "", 1000, -0.5, 2.5);
h_EoP_CPene_calib_fraction_[depletion1_] = ibooker.book1D("h_EoP_CPene_200_calib_fraction", "", 1000, -0.5, 2.5);
h_EoP_CPene_calib_fraction_[depletion2_] = ibooker.book1D("h_EoP_CPene_300_calib_fraction", "", 1000, -0.5, 2.5);
Expand Down Expand Up @@ -179,14 +182,16 @@ void HGCalHitCalibration::fillWithRecHits(std::map<DetId, const HGCRecHit*>& hit
LayerOccupancy_->Fill(layer);
if (seedEnergy < hitmap[hitid]->energy()) {
seedEnergy = hitmap[hitid]->energy();
seedDet = recHitTools_.getSiThickness(hitid);
seedDet = std::rint(recHitTools_.getSiThickness(hitid));
if (hitid.det() == DetId::HGCalHSc) {
seedDet = scint_;
}
}
}

void HGCalHitCalibration::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
float constexpr max_dR2 = 0.0025;

const edm::ESHandle<CaloGeometry>& geom = iSetup.getHandle(caloGeomToken_);
recHitTools_.setGeometry(*geom);

Expand All @@ -203,9 +208,21 @@ void HGCalHitCalibration::analyze(const edm::Event& iEvent, const edm::EventSetu
std::map<DetId, const HGCRecHit*> hitmap;
switch (algo_) {
case 1: {
const auto& rechitsEE = iEvent.get(recHitsEE_);
const auto& rechitsFH = iEvent.get(recHitsFH_);
const auto& rechitsBH = iEvent.get(recHitsBH_);
const auto& rechitsEE_handle = iEvent.getHandle(recHitsEE_);
if (!rechitsEE_handle.isValid())
return;

const auto& rechitsFH_handle = iEvent.getHandle(recHitsFH_);
if (!rechitsFH_handle.isValid())
return;

const auto& rechitsBH_handle = iEvent.getHandle(recHitsBH_);
if (!rechitsBH_handle.isValid())
return;

auto const& rechitsEE = *rechitsEE_handle;
auto const& rechitsFH = *rechitsFH_handle;
auto const& rechitsBH = *rechitsBH_handle;
for (unsigned int i = 0; i < rechitsEE.size(); ++i) {
hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
}
Expand All @@ -218,15 +235,27 @@ void HGCalHitCalibration::analyze(const edm::Event& iEvent, const edm::EventSetu
break;
}
case 2: {
const HGCRecHitCollection& rechitsEE = iEvent.get(recHitsEE_);
const auto& rechitsEE_handle = iEvent.getHandle(recHitsEE_);
if (!rechitsEE_handle.isValid())
return;

auto const& rechitsEE = *rechitsEE_handle;
for (unsigned int i = 0; i < rechitsEE.size(); i++) {
hitmap[rechitsEE[i].detid()] = &rechitsEE[i];
}
break;
}
case 3: {
const auto& rechitsFH = iEvent.get(recHitsFH_);
const auto& rechitsBH = iEvent.get(recHitsBH_);
const auto& rechitsFH_handle = iEvent.getHandle(recHitsFH_);
if (!rechitsFH_handle.isValid())
return;

const auto& rechitsBH_handle = iEvent.getHandle(recHitsBH_);
if (!rechitsBH_handle.isValid())
return;

auto const& rechitsFH = *rechitsFH_handle;
auto const& rechitsBH = *rechitsBH_handle;
for (unsigned int i = 0; i < rechitsFH.size(); i++) {
hitmap[rechitsFH[i].detid()] = &rechitsFH[i];
}
Expand Down Expand Up @@ -286,29 +315,30 @@ void HGCalHitCalibration::analyze(const edm::Event& iEvent, const edm::EventSetu
h_EoP_CPene_calib_fraction_[seedDet]->Fill(sumCalibRecHitCalib_fraction / it_caloPart.energy());

// Loop on reconstructed SC.
const auto& clusters = *hgcalMultiClustersHandle;
float total_energy = 0.;
float max_dR2 = 0.0025;
auto closest =
std::min_element(clusters.begin(), clusters.end(), [&](const reco::PFCluster& a, const reco::PFCluster& b) {
auto dR2_a = reco::deltaR2(it_caloPart, a);
auto dR2_b = reco::deltaR2(it_caloPart, b);
auto ERatio_a = a.correctedEnergy() / it_caloPart.energy();
auto ERatio_b = b.correctedEnergy() / it_caloPart.energy();
// If both clusters are within 0.0025, mark as first (min) the
// element with the highest ratio against the SimCluster
if (dR2_a < max_dR2 && dR2_b < max_dR2)
return ERatio_a > ERatio_b;
return dR2_a < dR2_b;
});
if (closest != clusters.end() && reco::deltaR2(*closest, it_caloPart) < 0.01) {
total_energy = closest->correctedEnergy();
seedDet = recHitTools_.getSiThickness(closest->seed());
if (closest->seed().det() == DetId::HGCalHSc) {
seedDet = scint_;
}
if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy / it_caloPart.energy());
if (hgcalMultiClustersHandle.isValid()) {
const auto& clusters = *hgcalMultiClustersHandle;
float total_energy = 0.;
auto closest =
std::min_element(clusters.begin(), clusters.end(), [&](const reco::PFCluster& a, const reco::PFCluster& b) {
auto dR2_a = reco::deltaR2(it_caloPart, a);
auto dR2_b = reco::deltaR2(it_caloPart, b);
auto ERatio_a = a.correctedEnergy() / it_caloPart.energy();
auto ERatio_b = b.correctedEnergy() / it_caloPart.energy();
// If both clusters are within 0.0025, mark as first (min) the
// element with the highest ratio against the SimCluster
if (dR2_a < max_dR2 && dR2_b < max_dR2)
return ERatio_a > ERatio_b;
return dR2_a < dR2_b;
});
if (closest != clusters.end() && reco::deltaR2(*closest, it_caloPart) < 0.01) {
total_energy = closest->correctedEnergy();
seedDet = recHitTools_.getSiThickness(closest->seed());
if (closest->seed().det() == DetId::HGCalHSc) {
seedDet = scint_;
}
if (hgcal_EoP_CPene_calib_fraction_.count(seedDet)) {
hgcal_EoP_CPene_calib_fraction_[seedDet]->Fill(total_energy / it_caloPart.energy());
}
}
}

Expand All @@ -323,6 +353,7 @@ void HGCalHitCalibration::analyze(const edm::Event& iEvent, const edm::EventSetu
return ERatio_a > ERatio_b;
return dR2_a < dR2_b;
};

// ELECTRONS in HGCAL
if (PFElectronHandle.isValid()) {
auto const& ele = (*PFElectronHandle);
Expand Down Expand Up @@ -366,6 +397,7 @@ void HGCalHitCalibration::fillDescriptions(edm::ConfigurationDescriptions& descr
edm::ParameterSetDescription desc;
desc.add<int>("debug", 0);
desc.add<bool>("rawRecHits", true);
desc.add<std::string>("folder", "HGCalHitCalibration");
desc.add<std::string>("detector", "all");
desc.add<int>("depletionFine", 120);
desc.add<edm::InputTag>("caloParticles", edm::InputTag("mix", "MergedCaloTruth"));
Expand Down