1 #ifndef RecoEgamma_EgammayHLTProducers_HLTRecHitInAllL1RegionsProducer_h_ 2 #define RecoEgamma_EgammayHLTProducers_HLTRecHitInAllL1RegionsProducer_h_ 67 minEt_(para.getParameter<double>(
"minEt")),
68 maxEt_(para.getParameter<double>(
"maxEt")),
69 regionEtaMargin_(para.getParameter<double>(
"regionEtaMargin")),
70 regionPhiMargin_(para.getParameter<double>(
"regionPhiMargin")),
71 token_(consumesColl.consumes<T1>(para.getParameter<
edm::InputTag>(
"inputColl"))){}
74 template<
typename T2>
static typename T2::const_iterator
beginIt(
const T2&
coll){
return coll.begin();}
75 template<
typename T2>
static typename T2::const_iterator
endIt(
const T2&
coll){
return coll.end();}
84 template<
typename RecHitType>
111 template<
typename RecHitType>
114 const std::vector<edm::ParameterSet> l1InputRegions = para.
getParameter<std::vector<edm::ParameterSet>>(
"l1InputRegions");
115 for(
auto&
pset : l1InputRegions){
117 l1RegionData_.emplace_back(createL1RegionData(type,
pset,consumesCollector()));
119 recHitLabels_ =para.
getParameter<std::vector<edm::InputTag>>(
"recHitLabels");
120 productLabels_=para.
getParameter<std::vector<std::string>>(
"productLabels");
122 for (
unsigned int collNr=0; collNr<recHitLabels_.size(); collNr++) {
123 recHitTokens_.push_back(consumes<RecHitCollectionType>(recHitLabels_[collNr]));
124 produces<RecHitCollectionType> (productLabels_[collNr]);
127 template<
typename RecHitType>
131 std::vector<std::string> productTags;
132 productTags.push_back(
"EcalRegionalRecHitsEB");
133 productTags.push_back(
"EcalRegionalRecHitsEE");
134 desc.
add<std::vector<std::string>>(
"productLabels", productTags);
135 std::vector<edm::InputTag> recHitLabels;
136 recHitLabels.push_back(
edm::InputTag(
"hltEcalRegionalEgammaRecHit:EcalRecHitsEB"));
137 recHitLabels.push_back(
edm::InputTag(
"hltEcalRegionalEgammaRecHit:EcalRecHitsEE"));
138 recHitLabels.push_back(
edm::InputTag(
"hltESRegionalEgammaRecHit:EcalRecHitsES"));
139 desc.
add<std::vector<edm::InputTag>>(
"recHitLabels", recHitLabels);
140 std::vector<edm::ParameterSet> l1InputRegions;
149 l1InputRegions.push_back(emIsoPSet);
154 emNonIsoPSet.
addParameter<
double>(
"regionEtaMargin",0.14);
155 emNonIsoPSet.
addParameter<
double>(
"regionPhiMargin",0.4);
157 l1InputRegions.push_back(emNonIsoPSet);
169 l1InputRegions.push_back(egPSet);
178 l1InputRegions.push_back(jetPSet);
183 l1InputRegionDesc.
add<
double>(
"minEt");
184 l1InputRegionDesc.
add<
double>(
"maxEt");
185 l1InputRegionDesc.
add<
double>(
"regionEtaMargin");
186 l1InputRegionDesc.
add<
double>(
"regionPhiMargin");
188 desc.
addVPSet(
"l1InputRegions",l1InputRegionDesc,l1InputRegions);
194 template<
typename RecHitType>
205 std::vector<EcalEtaPhiRegion> regions;
206 std::for_each(l1RegionData_.begin(),l1RegionData_.end(),
207 [&
event,®ions,l1CaloGeom](
const std::unique_ptr<L1RegionDataBase>&
input)
208 {
input->getEtaPhiRegions(event,regions,*l1CaloGeom);}
211 for(
size_t recHitCollNr=0;recHitCollNr<recHitTokens_.size();recHitCollNr++){
213 event.getByToken(recHitTokens_[recHitCollNr],recHits);
220 auto filteredRecHits = std::make_unique<RecHitCollectionType>();
222 if(!recHits->empty()){
224 if(!regions.empty()){
228 for(
const auto& region : regions){
229 if (region.inRegion(this_cell.
etaPos(),this_cell.
phiPos())) {
230 filteredRecHits->push_back(
recHit);
238 event.put(
std::move(filteredRecHits),productLabels_[recHitCollNr]);
247 template<
typename RecHitType>
250 if(type==
"L1EmParticle"){
252 }
else if(type==
"L1JetParticle"){
254 }
else if(type==
"L1MuonParticle"){
256 }
else if(type==
"EGamma"){
258 }
else if(type==
"Jet"){
260 }
else if(type==
"Muon"){
262 }
else if(type==
"Tau"){
266 throw cms::Exception(
"InvalidConfig") <<
" type "<<type<<
" is not recognised, this means the rec-hit you think you are keeping may not be and you should fix this error as it can lead to hard to find efficiency loses"<<std::endl;
272 template<
typename L1CollType>
276 event.getByToken(token_,l1Cands);
278 for(
auto l1CandIt = beginIt(*l1Cands);l1CandIt!=endIt(*l1Cands);++l1CandIt){
279 if(l1CandIt->et() >= minEt_ && l1CandIt->et() < maxEt_){
281 double etaLow = l1CandIt->eta() - regionEtaMargin_;
282 double etaHigh = l1CandIt->eta() + regionEtaMargin_;
283 double phiLow = l1CandIt->phi() - regionPhiMargin_;
284 double phiHigh = l1CandIt->phi() + regionPhiMargin_;
295 event.getByToken(token_,l1Cands);
297 for(
const auto& l1Cand : *l1Cands){
298 if(l1Cand.et() >= minEt_ && l1Cand.et() < maxEt_){
301 int etaIndex = l1Cand.gctJetCand()->etaIndex();
302 int phiIndex = l1Cand.gctJetCand()->phiIndex();
310 etaLow -= regionEtaMargin_;
311 etaHigh += regionEtaMargin_;
312 phiLow -= regionPhiMargin_;
313 phiHigh += regionPhiMargin_;
325 event.getByToken(token_,l1Cands);
327 for(
const auto& l1Cand : *l1Cands){
328 if(l1Cand.et() >= minEt_ && l1Cand.et() < maxEt_){
331 int etaIndex = l1Cand.gctEmCand()->etaIndex();
332 int phiIndex = l1Cand.gctEmCand()->phiIndex();
340 etaLow -= regionEtaMargin_;
341 etaHigh += regionEtaMargin_;
342 phiLow -= regionPhiMargin_;
343 phiHigh += regionPhiMargin_;
const_iterator end(int bx) const
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
std::string defaultModuleLabel()
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
double etaBinHighEdge(unsigned int etaIndex, bool central=true) const
HLTRecHitInAllL1RegionsProducer< EcalRecHit > HLTEcalRecHitInAllL1RegionsProducer
double etaBinLowEdge(unsigned int etaIndex, bool central=true) const
std::vector< std::unique_ptr< L1RegionDataBase > > l1RegionData_
void produce(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
L1RegionData(const edm::ParameterSet ¶, edm::ConsumesCollector &consumesColl)
virtual ~L1RegionDataBase()
def setup(process, global_tag, zero_tesla=False)
static BXVector< T2 >::const_iterator beginIt(const BXVector< T2 > &coll)
HLTRecHitInAllL1RegionsProducer(const edm::ParameterSet &ps)
std::vector< edm::InputTag > recHitLabels_
virtual const CaloCellGeometry * getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
static std::string const input
void getEtaPhiRegions(const edm::Event &, std::vector< EcalEtaPhiRegion > &, const L1CaloGeometry &) const override
std::vector< std::string > productLabels_
void addParameter(std::string const &name, T const &value)
HLTRecHitInAllL1RegionsProducer< EcalUncalibratedRecHit > HLTEcalUncalibratedRecHitInAllL1RegionsProducer
SubDetector subDetGeom[21]
ParameterDescriptionBase * add(U const &iLabel, T const &value)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
edm::EDGetTokenT< T1 > token_
static T2::const_iterator endIt(const T2 &coll)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
~HLTRecHitInAllL1RegionsProducer() override
static BXVector< T2 >::const_iterator endIt(const BXVector< T2 > &coll)
static T2::const_iterator beginIt(const T2 &coll)
virtual void getEtaPhiRegions(const edm::Event &, std::vector< EcalEtaPhiRegion > &, const L1CaloGeometry &) const =0
double emJetPhiBinLowEdge(unsigned int phiIndex) const
std::vector< edm::EDGetTokenT< RecHitCollectionType > > recHitTokens_
const_iterator begin(int bx) const
L1RegionDataBase * createL1RegionData(const std::string &, const edm::ParameterSet &, edm::ConsumesCollector &&)
double emJetPhiBinHighEdge(unsigned int phiIndex) const