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 Hcal reco producer with FACILE algorithm #1

Open
wants to merge 20 commits into
base: master
Choose a base branch
from
2 changes: 2 additions & 0 deletions RecoLocalCalo/HcalRecProducers/BuildFile.xml
Original file line number Diff line number Diff line change
Expand Up @@ -6,3 +6,5 @@
<use name="DataFormats/Common"/>
<use name="Geometry/Records"/>
<use name="boost"/>
<use name="HeterogeneousCore/SonicTriton"/>
<use name="HeterogeneousCore/SonicCore"/>
98 changes: 98 additions & 0 deletions RecoLocalCalo/HcalRecProducers/src/FacileHcalReconstructor.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,98 @@
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
#include "HeterogeneousCore/SonicTriton/interface/TritonClient.h"
#include "HeterogeneousCore/SonicCore/interface/SonicEDProducer.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "DataFormats/HcalRecHit/interface/HBHERecHit.h"
#include "Geometry/CaloTopology/interface/HcalTopology.h"
#include "Geometry/Records/interface/HcalRecNumberingRecord.h"

class FacileHcalReconstructor : public SonicEDProducer<TritonClient> {
public:
explicit FacileHcalReconstructor(edm::ParameterSet const& cfg)
: SonicEDProducer<TritonClient>(cfg),
fChannelInfoName_(cfg.getParameter<edm::InputTag>("ChannelInfoName")),
fTokChannelInfo_(consumes<HBHEChannelInfoCollection>(fChannelInfoName_)),
htopoToken_(esConsumes<HcalTopology, HcalRecNumberingRecord>()) {
produces<HBHERecHitCollection>();
setDebugName("FacileHcalReconstructor");
}

void acquire(edm::Event const& iEvent, edm::EventSetup const& iSetup, Input& iInput) override {
edm::Handle<HBHEChannelInfoCollection> hChannelInfo;
iEvent.getByToken(fTokChannelInfo_, hChannelInfo);

const HcalTopology* htopo = &iSetup.getData(htopoToken_);

auto& input1 = iInput.begin()->second;
auto data1 = std::make_shared<TritonInput<float>>();
data1->reserve(hChannelInfo->size());
client_.setBatchSize(hChannelInfo->size());

hcalIds_.clear();

for (const auto& pChannel : *hChannelInfo) {
std::vector<float> input;
const HcalDetId pDetId = pChannel.id();
hcalIds_.push_back(pDetId);

//inputs for Facile: iphi, gain, raw[8], depth (categorical), ieta (categorical)
input.push_back(pDetId.iphi());
input.push_back(pChannel.tsGain(0.));
for (unsigned int itTS = 0; itTS < pChannel.nSamples(); ++itTS) {
input.push_back(pChannel.tsRawCharge(itTS));
}

for (int itDepth = 1; itDepth <= htopo->maxDepth(); itDepth++) {
input.push_back(pDetId.depth() == itDepth);
}

for (int itIeta = 0; itIeta <= htopo->lastHERing(); itIeta++) {
input.push_back(pDetId.ietaAbs() == itIeta);
}

data1->push_back(input);
}

input1.toServer(data1);
}

void produce(edm::Event& iEvent, edm::EventSetup const& iSetup, Output const& iOutput) override {
std::unique_ptr<HBHERecHitCollection> out;
out = std::make_unique<HBHERecHitCollection>();
out->reserve(hcalIds_.size());

const auto& output1 = iOutput.begin()->second;
const auto& outputs = output1.fromServer<float>();

for (std::size_t iB = 0; iB < hcalIds_.size(); iB++) {
float rhE = outputs[iB][0];
if (rhE < 0. or std::isnan(rhE) or std::isinf(rhE))
rhE = 0;

HBHERecHit rh = HBHERecHit(hcalIds_[iB], rhE, 0.f, 0.f);
out->push_back(rh);
}
iEvent.put(std::move(out));
}

static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
TritonClient::fillPSetDescription(desc);
desc.add<edm::InputTag>("ChannelInfoName");
descriptions.add("FacileHcalReconstructor", desc);
}

private:
edm::InputTag fChannelInfoName_;
edm::EDGetTokenT<HBHEChannelInfoCollection> fTokChannelInfo_;
std::vector<HcalDetId> hcalIds_;
edm::ESGetToken<HcalTopology, HcalRecNumberingRecord> htopoToken_;
};

DEFINE_FWK_MODULE(FacileHcalReconstructor);
32 changes: 32 additions & 0 deletions RecoLocalCalo/HcalRecProducers/test/sonic_hlt_test.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,32 @@
import FWCore.ParameterSet.Config as cms
import os,sys

sys.path = sys.path + [os.path.expandvars("$CMSSW_BASE/src/HLTrigger/Configuration/test/"), os.path.expandvars("$CMSSW_RELEASE_BASE/src/HLTrigger/Configuration/test/")]

from OnLine_HLT_GRun import process

process.hltHbherecopre = process.hltHbhereco.clone(
makeRecHits = cms.bool(False),
saveInfos = cms.bool(True),
)

process.hltHbhereco = cms.EDProducer("FacileHcalReconstructor",
Client = cms.PSet(
batchSize = cms.untracked.uint32(16000),
address = cms.untracked.string("0.0.0.0"),
port = cms.untracked.uint32(8001),
timeout = cms.untracked.uint32(300),
modelName = cms.string("facile_all_v2"),
mode = cms.string("Async"),
modelVersion = cms.int32(-1),
verbose = cms.untracked.bool(False),
allowedTries = cms.untracked.uint32(5),
outputs = cms.untracked.vstring("output/BiasAdd"),
),
ChannelInfoName = cms.InputTag("hltHbherecopre")
)

process.HLTDoLocalHcalSequence = cms.Sequence( process.hltHcalDigis + process.hltHbherecopre + process.hltHbhereco + process.hltHfprereco + process.hltHfreco + process.hltHoreco )
process.HLTStoppedHSCPLocalHcalReco = cms.Sequence( process.hltHcalDigis + process.hltHbherecopre + process.hltHbhereco)

process.source.fileNames = cms.untracked.vstring("file:RelVal_Raw_GRun_MC.root")