1 #ifndef RecoEgamma_EgammaHLTProducers_HLTCaloObjInRegionsProducer_h_
2 #define RecoEgamma_EgammaHLTProducers_HLTCaloObjInRegionsProducer_h_
64 EtaPhiRegion(
float iEta,
float iPhi,
float iDR,
float iDEta,
float iDPhi):
93 minEt_(para.getParameter<double>(
"minEt")),
94 maxEt_(para.getParameter<double>(
"maxEt")),
95 maxDeltaR_(para.getParameter<double>(
"maxDeltaR")),
96 maxDEta_(para.getParameter<double>(
"maxDEta")),
97 maxDPhi_(para.getParameter<double>(
"maxDPhi")),
98 token_(consumesColl.consumes<T1>(para.getParameter<edm::
InputTag>(
"inputColl"))){}
101 template<
typename T2>
static typename T2::const_iterator
beginIt(
const T2&
coll){
return coll.begin();}
102 template<
typename T2>
static typename T2::const_iterator
endIt(
const T2&
coll){
return coll.end();}
109 template<
typename CaloObjType,
126 static std::auto_ptr<CaloObjCollType>
129 const std::vector<EtaPhiRegion>& regions);
137 template<
typename CaloObjType,
typename CaloObjCollType>
141 for(
auto&
pset : etaPhiRegions){
143 etaPhiRegionData_.emplace_back(createEtaPhiRegionData(type,
pset,consumesCollector()));
146 outputProductNames_=para.
getParameter<std::vector<std::string>>(
"outputProductNames");
147 inputCollTags_=para.
getParameter<std::vector<edm::InputTag>>(
"inputCollTags");
148 for (
unsigned int collNr=0; collNr<inputCollTags_.size(); collNr++) {
149 inputTokens_.push_back(consumes<CaloObjCollType>(inputCollTags_[collNr]));
150 produces<CaloObjCollType> (outputProductNames_[collNr]);
154 template<
typename CaloObjType,
typename CaloObjCollType>
159 outputProductNames.push_back(
"EcalRegionalRecHitsEB");
161 std::vector<edm::InputTag> inputColls;
163 desc.
add<std::vector<edm::InputTag>>(
"inputCollTags", inputColls);
174 etaPhiRegions.push_back(ecalCandPSet);
179 etaPhiRegionDesc.
add<
double>(
"minEt");
180 etaPhiRegionDesc.
add<
double>(
"maxEt");
181 etaPhiRegionDesc.
add<
double>(
"maxDeltaR");
182 etaPhiRegionDesc.
add<
double>(
"maxDEta");
183 etaPhiRegionDesc.
add<
double>(
"maxDPhi");
185 desc.
addVPSet(
"etaPhiRegions",etaPhiRegionDesc,etaPhiRegions);
191 template<
typename CaloObjType,
typename CaloObjCollType>
198 std::vector<EtaPhiRegion> regions;
199 std::for_each(etaPhiRegionData_.begin(),etaPhiRegionData_.end(),
200 [&
event,®ions](
const std::unique_ptr<EtaPhiRegionDataBase>&
input)
201 {
input->getEtaPhiRegions(event,regions);}
204 for(
size_t inputCollNr=0;inputCollNr<inputTokens_.size();inputCollNr++){
206 event.getByToken(inputTokens_[inputCollNr],inputColl);
209 edm::LogError(
"ProductNotFound")<<
"could not get a handle on the "<<
typeid(CaloObjCollType).
name() <<
" named "<< inputCollTags_[inputCollNr].encode() << std::endl;
212 auto outputColl = makeFilteredColl(inputColl,caloGeomHandle,regions);
213 event.put(
outputColl,outputProductNames_[inputCollNr]);
217 template<
typename CaloObjType,
typename CaloObjCollType>
218 std::auto_ptr<CaloObjCollType>
222 const std::vector<EtaPhiRegion>& regions)
225 std::auto_ptr<CaloObjCollType>
outputColl(
new CaloObjCollType);
226 if(!inputColl->empty()){
228 if(!regions.empty()){
229 for(
const CaloObjType&
obj : *inputColl){
231 if(objGeom==
nullptr){
235 edm::LogError(
"HLTCaloObjInRegionsProducer") <<
"for an object of type "<<
typeid(CaloObjType).
name()<<
" the geometry returned null for id "<<
obj.id()<<
" in HLTCaloObjsInRegion, this shouldnt be possible and something has gone wrong, auto accepting hit";
236 outputColl->push_back(
obj);
241 for(
const auto&
region : regions){
243 outputColl->push_back(
obj);
258 template<
typename CaloObjType,
typename CaloObjCollType>
264 if(type==
"L1EGamma"){
266 }
else if(type==
"L1Jet"){
268 }
else if(type==
"L1Muon"){
270 }
else if(type==
"L1Tau"){
272 }
else if(type==
"RecoEcalCandidate"){
274 }
else if(type==
"RecoChargedCandidate"){
276 }
else if(type==
"Electron"){
280 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;
285 template<
typename CandCollType>
289 event.getByToken(token_,cands);
291 for(
auto candIt = beginIt(*cands);candIt!=endIt(*cands);++candIt){
292 if(candIt->et() >= minEt_ && (maxEt_<0 || candIt->et() < maxEt_)){
293 regions.push_back(
EtaPhiRegion(candIt->eta(),candIt->phi(),
294 maxDeltaR_,maxDEta_,maxDPhi_));
const_iterator end(int bx) const
T getParameter(std::string const &) const
std::string defaultModuleLabel()
ParameterDescriptionBase * addVPSet(U const &iLabel, ParameterSetDescription const &validator, std::vector< ParameterSet > const &defaults)
std::vector< edm::EDGetTokenT< CaloObjCollType > > inputTokens_
std::vector< edm::InputTag > inputCollTags_
void getEtaPhiRegions(const edm::Event &, std::vector< EtaPhiRegion > &) const override
EtaPhiRegionDataBase * createEtaPhiRegionData(const std::string &, const edm::ParameterSet &, edm::ConsumesCollector &&)
#define DEFINE_FWK_MODULE(type)
Geom::Phi< T > phi() const
static BXVector< T2 >::const_iterator beginIt(const BXVector< T2 > &coll)
static T2::const_iterator beginIt(const T2 &coll)
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
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< std::string > outputProductNames_
HLTCaloObjInRegionsProducer< HBHEDataFrame > HLTHcalDigisInRegionsProducer
static T2::const_iterator endIt(const T2 &coll)
void addParameter(std::string const &name, T const &value)
bool operator()(float eta, float phi) const
Abs< T >::type abs(const T &t)
EtaPhiRegionData(const edm::ParameterSet ¶, edm::ConsumesCollector &consumesColl)
virtual void getEtaPhiRegions(const edm::Event &, std::vector< EtaPhiRegion > &) const =0
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
HLTCaloObjInRegionsProducer(const edm::ParameterSet &ps)
~HLTCaloObjInRegionsProducer()
static std::auto_ptr< CaloObjCollType > makeFilteredColl(const edm::Handle< CaloObjCollType > &inputColl, const edm::ESHandle< CaloGeometry > &caloGeomHandle, const std::vector< EtaPhiRegion > ®ions)
double deltaPhi(double phi1, double phi2)
SubDetector subDetGeom[18]
edm::EDGetTokenT< T1 > token_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Geom::Phi< T > phi() const
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
static BXVector< T2 >::const_iterator endIt(const BXVector< T2 > &coll)
const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
const_iterator begin(int bx) const
EtaPhiRegion(float iEta, float iPhi, float iDR, float iDEta, float iDPhi)
std::vector< std::unique_ptr< EtaPhiRegionDataBase > > etaPhiRegionData_
void produce(edm::Event &, const edm::EventSetup &) override