CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
13 
15 
16 // Reco candidates (might not need)
23 
24 
25 // Geometry and topology
30 
32 
36 
40 
42 
44 
45 /**************************************************************
46 / purpose: enable filtering of calo objects in eta/phi or deltaR
47 / regions around generic objects
48 /
49 / operation : accepts all objects with
50 / (dEta <dEtaMax && dPhi < dPhiMax) || dR < dRMax
51 / so the OR of a rectangular region and cone region
52 ****************************************************************/
53 
54 //this is a struct which contains all the eta/phi regions
55 //from which to filter the calo objs
57 private:
58  float centreEta_;
59  float centrePhi_;
60  float maxDeltaR2_;
61  float maxDEta_;
62  float maxDPhi_;
63  public:
64  EtaPhiRegion(float iEta,float iPhi,float iDR,float iDEta,float iDPhi):
65  centreEta_(iEta),centrePhi_(iPhi),maxDeltaR2_(iDR*iDR),
66  maxDEta_(iDEta),maxDPhi_(iDPhi){}
68  bool operator()(float eta,float phi)const{
71 
72 };
73 
75 public:
77  virtual void getEtaPhiRegions(const edm::Event&,std::vector<EtaPhiRegion>&)const=0;
78 };
79 
80 
81 //this class stores the tokens to access the objects around which we wish to filter
82 //it makes a vector of EtaPhiRegions which are then used to filter the CaloObjs
83 template<typename T1> class EtaPhiRegionData : public EtaPhiRegionDataBase{
84 private:
85  float minEt_;
86  float maxEt_;
87  float maxDeltaR_;
88  float maxDEta_;
89  float maxDPhi_;
91 public:
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"))){}
99 
100  void getEtaPhiRegions(const edm::Event&,std::vector<EtaPhiRegion>&)const override;
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();}
103  template<typename T2>static typename BXVector<T2>::const_iterator beginIt(const BXVector<T2>& coll){return coll.begin(0);}
104  template<typename T2>static typename BXVector<T2>::const_iterator endIt(const BXVector<T2>& coll){return coll.end(0);}
105 
106 };
107 
108 
109 template<typename CaloObjType,
110  typename CaloObjCollType =edm::SortedCollection<CaloObjType> >
112 
113 
114  public:
115 
116 
117 
120 
121  void produce(edm::Event&, const edm::EventSetup&) override;
122  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
123 
124  private:
126  static std::auto_ptr<CaloObjCollType>
128  const edm::ESHandle<CaloGeometry>& caloGeomHandle,
129  const std::vector<EtaPhiRegion>& regions);
130  std::vector<std::string> outputProductNames_;
131  std::vector<edm::InputTag> inputCollTags_;
132  std::vector<edm::EDGetTokenT<CaloObjCollType>> inputTokens_;
133  std::vector<std::unique_ptr<EtaPhiRegionDataBase>> etaPhiRegionData_;
134 };
135 
136 
137 template<typename CaloObjType,typename CaloObjCollType>
139 {
140  const std::vector<edm::ParameterSet> etaPhiRegions = para.getParameter<std::vector<edm::ParameterSet>>("etaPhiRegions");
141  for(auto& pset : etaPhiRegions){
142  const std::string type=pset.getParameter<std::string>("type");
143  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
144  }
145 
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]);
151  }
152 }
153 
154 template<typename CaloObjType,typename CaloObjCollType>
156 {
158  std::vector<std::string> outputProductNames;
159  outputProductNames.push_back("EcalRegionalRecHitsEB");
160  desc.add<std::vector<std::string>>("outputProductNames", outputProductNames);
161  std::vector<edm::InputTag> inputColls;
162  inputColls.push_back(edm::InputTag("hltHcalDigis"));
163  desc.add<std::vector<edm::InputTag>>("inputCollTags", inputColls);
164  std::vector<edm::ParameterSet> etaPhiRegions;
165 
166  edm::ParameterSet ecalCandPSet;
167  ecalCandPSet.addParameter<std::string>("type","RecoEcalCandidate");
168  ecalCandPSet.addParameter<double>("minEt",-1);
169  ecalCandPSet.addParameter<double>("maxEt",-1);
170  ecalCandPSet.addParameter<double>("maxDeltaR",0.5);
171  ecalCandPSet.addParameter<double>("maxDEta",0.);
172  ecalCandPSet.addParameter<double>("maxDPhi",0.);
173  ecalCandPSet.addParameter<edm::InputTag>("inputColl",edm::InputTag("hltEgammaCandidates"));
174  etaPhiRegions.push_back(ecalCandPSet);
175 
176 
177  edm::ParameterSetDescription etaPhiRegionDesc;
178  etaPhiRegionDesc.add<std::string>("type");
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");
184  etaPhiRegionDesc.add<edm::InputTag>("inputColl");
185  desc.addVPSet("etaPhiRegions",etaPhiRegionDesc,etaPhiRegions);
186 
188 }
189 
190 
191 template<typename CaloObjType,typename CaloObjCollType>
193 
194  // get the collection geometry:
195  edm::ESHandle<CaloGeometry> caloGeomHandle;
196  setup.get<CaloGeometryRecord>().get(caloGeomHandle);
197 
198  std::vector<EtaPhiRegion> regions;
199  std::for_each(etaPhiRegionData_.begin(),etaPhiRegionData_.end(),
200  [&event,&regions](const std::unique_ptr<EtaPhiRegionDataBase>& input)
201  {input->getEtaPhiRegions(event,regions);}
202  );
203 
204  for(size_t inputCollNr=0;inputCollNr<inputTokens_.size();inputCollNr++){
206  event.getByToken(inputTokens_[inputCollNr],inputColl);
207 
208  if (!(inputColl.isValid())) {
209  edm::LogError("ProductNotFound")<< "could not get a handle on the "<<typeid(CaloObjCollType).name() <<" named "<< inputCollTags_[inputCollNr].encode() << std::endl;
210  continue;
211  }
212  auto outputColl = makeFilteredColl(inputColl,caloGeomHandle,regions);
213  event.put(outputColl,outputProductNames_[inputCollNr]);
214  }
215 }
216 
217 template<typename CaloObjType,typename CaloObjCollType>
218 std::auto_ptr<CaloObjCollType>
221  const edm::ESHandle<CaloGeometry>& caloGeomHandle,
222  const std::vector<EtaPhiRegion>& regions)
223 {
224 
225  std::auto_ptr<CaloObjCollType> outputColl(new CaloObjCollType);
226  if(!inputColl->empty()){
227  const CaloSubdetectorGeometry* subDetGeom=caloGeomHandle->getSubdetectorGeometry(inputColl->front().id());
228  if(!regions.empty()){
229  for(const CaloObjType& obj : *inputColl){
230  const CaloCellGeometry* objGeom = subDetGeom->getGeometry(obj.id());
231  if(objGeom==nullptr){
232  //wondering what to do here
233  //something is very very wrong
234  //given HLT should never crash or throw, decided to log an error and
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);
237  }
238  float eta = objGeom->getPosition().eta();
239  float phi = objGeom->getPosition().phi();
240 
241  for(const auto& region : regions){
242  if(region(eta,phi)) {
243  outputColl->push_back(obj);
244  break;
245  }
246  }
247  }
248  }//end check of empty regions
249  }//end check of empty rec-hits
250  return outputColl;
251 
252 }
253 
254 
255 
256 
257 
258 template<typename CaloObjType,typename CaloObjCollType>
262  edm::ConsumesCollector && consumesColl)
263 {
264  if(type=="L1EGamma"){
265  return new EtaPhiRegionData<l1t::EGammaBxCollection>(para,consumesColl);
266  }else if(type=="L1Jet"){
267  return new EtaPhiRegionData<l1t::JetBxCollection>(para,consumesColl);
268  }else if(type=="L1Muon"){
269  return new EtaPhiRegionData<l1t::MuonBxCollection>(para,consumesColl);
270  }else if(type=="L1Tau"){
271  return new EtaPhiRegionData<l1t::TauBxCollection>(para,consumesColl);
272  }else if(type=="RecoEcalCandidate"){
273  return new EtaPhiRegionData<reco::RecoEcalCandidateCollection>(para,consumesColl);
274  }else if(type=="RecoChargedCandidate"){
275  return new EtaPhiRegionData<reco::RecoChargedCandidateCollection>(para,consumesColl);
276  }else if(type=="Electron"){
277  return new EtaPhiRegionData<reco::Electron>(para,consumesColl);
278  }else{
279  //this is a major issue and could lead to rather subtle efficiency losses, so if its incorrectly configured, we're aborting the job!
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;
281  }
282 
283 }
284 
285 template<typename CandCollType>
286 void EtaPhiRegionData<CandCollType>::getEtaPhiRegions(const edm::Event& event,std::vector<EtaPhiRegion>&regions)const
287 {
289  event.getByToken(token_,cands);
290 
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_));
295  }
296 
297  }
298 }
299 
300 
303 
304 #endif
305 
306 
const_iterator end(int bx) const
type
Definition: HCALResponse.h:21
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 &&)
ProductID id() const
Definition: HandleBase.cc:15
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
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_
HLTCaloObjInRegionsProducer< HBHEDataFrame > HLTHcalDigisInRegionsProducer
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)
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
bool isValid() const
Definition: HandleBase.h:75
HLTCaloObjInRegionsProducer(const edm::ParameterSet &ps)
static std::auto_ptr< CaloObjCollType > makeFilteredColl(const edm::Handle< CaloObjCollType > &inputColl, const edm::ESHandle< CaloGeometry > &caloGeomHandle, const std::vector< EtaPhiRegion > &regions)
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:22
SubDetector subDetGeom[18]
JetCorrectorParametersCollection coll
Definition: classes.h:10
edm::EDGetTokenT< T1 > token_
const T & get() const
Definition: EventSetup.h:56
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Geom::Phi< T > phi() const
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
T eta() const
Definition: PV3DBase.h:76
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