CMS 3D CMS Logo

HLTCaloObjInRegionsProducer.h
Go to the documentation of this file.
1 #ifndef RecoEgamma_EgammaHLTProducers_HLTCaloObjInRegionsProducer_h_
2 #define RecoEgamma_EgammaHLTProducers_HLTCaloObjInRegionsProducer_h_
3 
4 #include <memory>
5 
16 
23 
24 
29 
31 
35 
38 
39 /**************************************************************
40 / purpose: enable filtering of calo objects in eta/phi or deltaR
41 / regions around generic objects
42 /
43 / operation : accepts all objects with
44 / (dEta <dEtaMax && dPhi < dPhiMax) || dR < dRMax
45 / so the OR of a rectangular region and cone region
46 ****************************************************************/
47 
48 //this is a struct which contains all the eta/phi regions
49 //from which to filter the calo objs
51 private:
52  float centreEta_;
53  float centrePhi_;
54  float maxDeltaR2_;
55  float maxDEta_;
56  float maxDPhi_;
57  public:
58  EtaPhiRegion(float iEta,float iPhi,float iDR,float iDEta,float iDPhi):
59  centreEta_(iEta),centrePhi_(iPhi),maxDeltaR2_(iDR*iDR),
60  maxDEta_(iDEta),maxDPhi_(iDPhi){}
62  bool operator()(float eta,float phi)const{
63  return reco::deltaR2(eta,phi,centreEta_,centrePhi_)<maxDeltaR2_ ||
64  (std::abs(eta-centreEta_)<maxDEta_ && std::abs(reco::deltaPhi(phi,centrePhi_))<maxDPhi_);}
65 
66 };
67 
69 public:
71  virtual ~EtaPhiRegionDataBase() = default;
72  virtual void getEtaPhiRegions(const edm::Event&,std::vector<EtaPhiRegion>&)const=0;
73 };
74 
75 
76 //this class stores the tokens to access the objects around which we wish to filter
77 //it makes a vector of EtaPhiRegions which are then used to filter the CaloObjs
78 template<typename T1> class EtaPhiRegionData : public EtaPhiRegionDataBase{
79 private:
80  float minEt_;
81  float maxEt_;
82  float maxDeltaR_;
83  float maxDEta_;
84  float maxDPhi_;
86 public:
88  minEt_(para.getParameter<double>("minEt")),
89  maxEt_(para.getParameter<double>("maxEt")),
90  maxDeltaR_(para.getParameter<double>("maxDeltaR")),
91  maxDEta_(para.getParameter<double>("maxDEta")),
92  maxDPhi_(para.getParameter<double>("maxDPhi")),
93  token_(consumesColl.consumes<T1>(para.getParameter<edm::InputTag>("inputColl"))){}
94 
95  void getEtaPhiRegions(const edm::Event&,std::vector<EtaPhiRegion>&)const override;
96  template<typename T2>static typename T2::const_iterator beginIt(const T2& coll){return coll.begin();}
97  template<typename T2>static typename T2::const_iterator endIt(const T2& coll){return coll.end();}
98  template<typename T2>static typename BXVector<T2>::const_iterator beginIt(const BXVector<T2>& coll){return coll.begin(0);}
99  template<typename T2>static typename BXVector<T2>::const_iterator endIt(const BXVector<T2>& coll){return coll.end(0);}
100 
101 };
102 
103 
104 template<typename CaloObjType,
105  typename CaloObjCollType =edm::SortedCollection<CaloObjType> >
107 
108 
109  public:
110 
111 
112 
115 
116  void produce(edm::Event&, const edm::EventSetup&) override;
117  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
118 
119  private:
120  EtaPhiRegionDataBase* createEtaPhiRegionData(const std::string&,const edm::ParameterSet&,edm::ConsumesCollector &&); //calling function owns this
121  static std::unique_ptr<CaloObjCollType>
122  makeFilteredColl(const edm::Handle<CaloObjCollType>& inputColl,
123  const edm::ESHandle<CaloGeometry>& caloGeomHandle,
124  const std::vector<EtaPhiRegion>& regions);
125  std::vector<std::string> outputProductNames_;
126  std::vector<edm::InputTag> inputCollTags_;
127  std::vector<edm::EDGetTokenT<CaloObjCollType>> inputTokens_;
128  std::vector<std::unique_ptr<EtaPhiRegionDataBase>> etaPhiRegionData_;
129 };
130 
131 
132 template<typename CaloObjType,typename CaloObjCollType>
134 {
135  const std::vector<edm::ParameterSet> etaPhiRegions = para.getParameter<std::vector<edm::ParameterSet>>("etaPhiRegions");
136  for(auto& pset : etaPhiRegions){
137  const std::string type=pset.getParameter<std::string>("type");
138  etaPhiRegionData_.emplace_back(createEtaPhiRegionData(type,pset,consumesCollector())); //meh I was going to use a factory but it was going to be overly complex for my needs
139  }
140 
141  outputProductNames_=para.getParameter<std::vector<std::string>>("outputProductNames");
142  inputCollTags_=para.getParameter<std::vector<edm::InputTag>>("inputCollTags");
143  if(outputProductNames_.size()!=inputCollTags_.size()){
144  throw cms::Exception("InvalidConfiguration") <<" error outputProductNames and inputCollTags must be the same size, they are "<<outputProductNames_.size()<<" vs "<<inputCollTags_.size();
145  }
146  for (unsigned int collNr=0; collNr<inputCollTags_.size(); collNr++) {
147  inputTokens_.push_back(consumes<CaloObjCollType>(inputCollTags_[collNr]));
148  produces<CaloObjCollType> (outputProductNames_[collNr]);
149  }
150 }
151 
152 template<typename CaloObjType,typename CaloObjCollType>
154 {
156  std::vector<std::string> outputProductNames;
157  outputProductNames.push_back("EcalRegionalRecHitsEB");
158  desc.add<std::vector<std::string>>("outputProductNames", outputProductNames);
159  std::vector<edm::InputTag> inputColls;
160  inputColls.push_back(edm::InputTag("hltHcalDigis"));
161  desc.add<std::vector<edm::InputTag>>("inputCollTags", inputColls);
162  std::vector<edm::ParameterSet> etaPhiRegions;
163 
164  edm::ParameterSet ecalCandPSet;
165  ecalCandPSet.addParameter<std::string>("type","RecoEcalCandidate");
166  ecalCandPSet.addParameter<double>("minEt",-1);
167  ecalCandPSet.addParameter<double>("maxEt",-1);
168  ecalCandPSet.addParameter<double>("maxDeltaR",0.5);
169  ecalCandPSet.addParameter<double>("maxDEta",0.);
170  ecalCandPSet.addParameter<double>("maxDPhi",0.);
171  ecalCandPSet.addParameter<edm::InputTag>("inputColl",edm::InputTag("hltEgammaCandidates"));
172  etaPhiRegions.push_back(ecalCandPSet);
173 
174 
175  edm::ParameterSetDescription etaPhiRegionDesc;
176  etaPhiRegionDesc.add<std::string>("type");
177  etaPhiRegionDesc.add<double>("minEt");
178  etaPhiRegionDesc.add<double>("maxEt");
179  etaPhiRegionDesc.add<double>("maxDeltaR");
180  etaPhiRegionDesc.add<double>("maxDEta");
181  etaPhiRegionDesc.add<double>("maxDPhi");
182  etaPhiRegionDesc.add<edm::InputTag>("inputColl");
183  desc.addVPSet("etaPhiRegions",etaPhiRegionDesc,etaPhiRegions);
184 
186 }
187 
188 
189 template<typename CaloObjType,typename CaloObjCollType>
191 
192  // get the collection geometry:
193  edm::ESHandle<CaloGeometry> caloGeomHandle;
194  setup.get<CaloGeometryRecord>().get(caloGeomHandle);
195 
196  std::vector<EtaPhiRegion> regions;
197  std::for_each(etaPhiRegionData_.begin(),etaPhiRegionData_.end(),
198  [&event,&regions](const std::unique_ptr<EtaPhiRegionDataBase>& input)
199  {input->getEtaPhiRegions(event,regions);}
200  );
201 
202  for(size_t inputCollNr=0;inputCollNr<inputTokens_.size();inputCollNr++){
204  event.getByToken(inputTokens_[inputCollNr],inputColl);
205 
206  if (!(inputColl.isValid())) {
207  edm::LogError("ProductNotFound")<< "could not get a handle on the "<<typeid(CaloObjCollType).name() <<" named "<< inputCollTags_[inputCollNr].encode() << std::endl;
208  continue;
209  }
210  auto outputColl = makeFilteredColl(inputColl,caloGeomHandle,regions);
211  event.put(std::move(outputColl),outputProductNames_[inputCollNr]);
212  }
213 }
214 
215 template<typename CaloObjType,typename CaloObjCollType>
216 std::unique_ptr<CaloObjCollType>
219  const edm::ESHandle<CaloGeometry>& caloGeomHandle,
220  const std::vector<EtaPhiRegion>& regions)
221 {
222 
223  auto outputColl = std::make_unique<CaloObjCollType>();
224  if(!inputColl->empty()){
225  const CaloSubdetectorGeometry* subDetGeom=caloGeomHandle->getSubdetectorGeometry(inputColl->begin()->id());
226  if(!regions.empty()){
227  for(const CaloObjType& obj : *inputColl){
228  const CaloCellGeometry* objGeom = subDetGeom->getGeometry(obj.id());
229  if(objGeom==nullptr){
230  //wondering what to do here
231  //something is very very wrong
232  //given HLT should never crash or throw, decided to log an error and
233  edm::LogError("HLTCaloObjInRegionsProducer") << "for an object of type "<<typeid(CaloObjType).name()<<" the geometry returned null for id "<<DetId(obj.id()).rawId()<<" in HLTCaloObjsInRegion, this shouldnt be possible and something has gone wrong, auto accepting hit";
234  outputColl->push_back(obj);
235  }
236  float eta = objGeom->getPosition().eta();
237  float phi = objGeom->getPosition().phi();
238 
239  for(const auto& region : regions){
240  if(region(eta,phi)) {
241  outputColl->push_back(obj);
242  break;
243  }
244  }
245  }
246  }//end check of empty regions
247  }//end check of empty rec-hits
248  return outputColl;
249 
250 }
251 
252 
253 
254 
255 
256 template<typename CaloObjType,typename CaloObjCollType>
260  edm::ConsumesCollector && consumesColl)
261 {
262  if(type=="L1EGamma"){
263  return new EtaPhiRegionData<l1t::EGammaBxCollection>(para,consumesColl);
264  }else if(type=="L1Jet"){
265  return new EtaPhiRegionData<l1t::JetBxCollection>(para,consumesColl);
266  }else if(type=="L1Muon"){
267  return new EtaPhiRegionData<l1t::MuonBxCollection>(para,consumesColl);
268  }else if(type=="L1Tau"){
269  return new EtaPhiRegionData<l1t::TauBxCollection>(para,consumesColl);
270  }else if(type=="RecoEcalCandidate"){
271  return new EtaPhiRegionData<reco::RecoEcalCandidateCollection>(para,consumesColl);
272  }else if(type=="RecoChargedCandidate"){
273  return new EtaPhiRegionData<reco::RecoChargedCandidateCollection>(para,consumesColl);
274  }else if(type=="Electron"){
275  return new EtaPhiRegionData<reco::Electron>(para,consumesColl);
276  }else{
277  //this is a major issue and could lead to rather subtle efficiency losses, so if its incorrectly configured, we're aborting the job!
278  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;
279  }
280 
281 }
282 
283 template<typename CandCollType>
284 void EtaPhiRegionData<CandCollType>::getEtaPhiRegions(const edm::Event& event,std::vector<EtaPhiRegion>&regions)const
285 {
287  event.getByToken(token_,cands);
288 
289  for(auto candIt = beginIt(*cands);candIt!=endIt(*cands);++candIt){
290  if(candIt->et() >= minEt_ && (maxEt_<0 || candIt->et() < maxEt_)){
291  regions.push_back(EtaPhiRegion(candIt->eta(),candIt->phi(),
292  maxDeltaR_,maxDEta_,maxDPhi_));
293  }
294 
295  }
296 }
297 
298 
299 
300 #endif
301 
302 
const_iterator end(int bx) const
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:45
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_
EtaPhiRegionDataBase * createEtaPhiRegionData(const std::string &, const edm::ParameterSet &, edm::ConsumesCollector &&)
ProductID id() const
Definition: HandleBase.cc:15
virtual const GlobalPoint & getPosition() const
Returns the position of reference for this cell.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
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
Definition: EdmProvDump.cc:44
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< std::string > outputProductNames_
static std::unique_ptr< CaloObjCollType > makeFilteredColl(const edm::Handle< CaloObjCollType > &inputColl, const edm::ESHandle< CaloGeometry > &caloGeomHandle, const std::vector< EtaPhiRegion > &regions)
static T2::const_iterator endIt(const T2 &coll)
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:144
bool operator()(float eta, float phi) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
EtaPhiRegionData(const edm::ParameterSet &para, edm::ConsumesCollector &consumesColl)
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
bool isValid() const
Definition: HandleBase.h:74
HLTCaloObjInRegionsProducer(const edm::ParameterSet &ps)
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
Definition: DetId.h:18
JetCorrectorParametersCollection coll
Definition: classes.h:10
edm::EDGetTokenT< T1 > token_
const T & get() const
Definition: EventSetup.h:55
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
T eta() const
Definition: PV3DBase.h:76
HLT enums.
static BXVector< T2 >::const_iterator endIt(const BXVector< T2 > &coll)
void getEtaPhiRegions(const edm::Event &, std::vector< EtaPhiRegion > &) const override
const_iterator begin(int bx) const
EtaPhiRegion(float iEta, float iPhi, float iDR, float iDEta, float iDPhi)
std::vector< std::unique_ptr< EtaPhiRegionDataBase > > etaPhiRegionData_
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1
void produce(edm::Event &, const edm::EventSetup &) override