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
32 changes: 32 additions & 0 deletions HLTrigger/Configuration/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("HcalReconstructor",
Client = cms.PSet(
batchSize = cms.untracked.uint32(16000),
address = cms.untracked.string("ailab01.fnal.gov"),
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
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")
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"/>
118 changes: 118 additions & 0 deletions RecoLocalCalo/HcalRecProducers/src/HcalReconstructor.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,118 @@
#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"

class HcalReconstructor : public SonicEDProducer<TritonClient>
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
{
public:
explicit HcalReconstructor(edm::ParameterSet const& cfg) :
SonicEDProducer<TritonClient>(cfg),
fChannelInfoName_(cfg.getParameter<edm::InputTag>("ChannelInfoName")),
fTokChannelInfo_(this->consumes<HBHEChannelInfoCollection>(fChannelInfoName_))
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved

{
this->produces<HBHERecHitCollection>();
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
this->setDebugName("HcalReconstructor");
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
}

void acquire(edm::Event const& iEvent, edm::EventSetup const& iSetup, Input& iInput) override {

const HBHEChannelInfoCollection *channelInfo = 0;
edm::Handle<HBHEChannelInfoCollection> hChannelInfo;
iEvent.getByToken(fTokChannelInfo_, hChannelInfo);
channelInfo = hChannelInfo.product();

auto& input1 = iInput.begin()->second;
auto data1 = std::make_shared<TritonInput<float>>();
data1->reserve(std::distance(channelInfo->begin(), channelInfo->end()));
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved

hcalIds_.clear();

for(HBHEChannelInfoCollection::const_iterator itC = channelInfo->begin(); itC != channelInfo->end(); itC++){
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved

std::vector<float> input;
const HBHEChannelInfo& pChannel(*itC);
const HcalDetId pDetId = pChannel.id();
hcalIds_.push_back(pDetId);

//FACILE uses iphi as a continuous variable
input.push_back((float)pDetId.iphi());
input.push_back((float)pChannel.tsGain(0.));
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
for (unsigned int itTS=0; itTS < pChannel.nSamples(); ++itTS) {
input.push_back((float)pChannel.tsRawCharge(itTS));
}

//FACILE considers 7 Hcal depths as binary variables
for (int itDepth=1; itDepth < 8; itDepth++){
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
if (pDetId.depth() == itDepth) input.push_back(1.f);
else input.push_back(0.f);
}

//ieta is also encoded as a binary variable
for (int itIeta = 0; itIeta < 30; itIeta++){
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
if (std::abs(pDetId.ieta()) == itIeta) input.push_back(1.f);
else input.push_back(0.f);
}

data1->push_back(input);
}

//set batch at maximum RHcollsize; pad the remainder
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
unsigned last = std::distance(channelInfo->begin(), channelInfo->end());
if(last < input1.batchSize()){
std::vector<float> pad(47,0.f);
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
for(unsigned int iP = last; iP != input1.batchSize(); iP++){
data1->push_back(pad);
}
}

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.) rhE = 0.;//shouldn't be necessary, relu activation function
jeffkrupa marked this conversation as resolved.
Show resolved Hide resolved
//exception?
if(std::isnan(rhE)) rhE = 0;
if(std::isinf(rhE)) rhE = 0;

//FACILE does no time reco
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) {
Copy link

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

reactivate this and provide the necessary parameters
(just ChannelInfoName based on the constructor?)

edm::ParameterSetDescription desc;
TritonClient::fillPSetDescription(desc);
//add producer-specific parameters
desc.add<edm::InputTag>("ChannelInfoName","hbheprereco");
descriptions.add("HcalPhase1Reconstructor_FACILE",desc);
}*/


private:
edm::InputTag fChannelInfoName_;
edm::EDGetTokenT<HBHEChannelInfoCollection> fTokChannelInfo_;
std::vector<HcalDetId> hcalIds_;
};

DEFINE_FWK_MODULE(HcalReconstructor);