Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:58:35

0001 #include "FWCore/Utilities/interface/EDGetToken.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/Framework/interface/stream/EDProducer.h"
0008 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
0009 #include "DataFormats/Common/interface/ValueMap.h"
0010 
0011 template <typename ObjType, typename ObjCollType>
0012 class HLTDQMObjSelector : public edm::stream::EDProducer<> {
0013 public:
0014   explicit HLTDQMObjSelector(const edm::ParameterSet& config);
0015   void produce(edm::Event&, edm::EventSetup const&) override;
0016   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0017 
0018 private:
0019   edm::EDGetTokenT<ObjCollType> token_;
0020   StringCutObjectSelector<ObjType, true> selection_;
0021 };
0022 
0023 template <typename ObjType, typename ObjCollType>
0024 HLTDQMObjSelector<ObjType, ObjCollType>::HLTDQMObjSelector(const edm::ParameterSet& config)
0025     : token_(consumes<ObjCollType>(config.getParameter<edm::InputTag>("objs"))),
0026       selection_(config.getParameter<std::string>("selection")) {
0027   produces<edm::ValueMap<bool> >();
0028 }
0029 
0030 template <typename ObjType, typename ObjCollType>
0031 void HLTDQMObjSelector<ObjType, ObjCollType>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0032   edm::ParameterSetDescription desc;
0033   desc.add<edm::InputTag>("objs", edm::InputTag(""));
0034   desc.add<std::string>("selection", "et > 5");
0035   descriptions.add("hltDQMObjSelector", desc);
0036 }
0037 
0038 template <typename ObjType, typename ObjCollType>
0039 void HLTDQMObjSelector<ObjType, ObjCollType>::produce(edm::Event& event, const edm::EventSetup& setup) {
0040   edm::Handle<ObjCollType> handle;
0041   event.getByToken(token_, handle);
0042 
0043   if (!handle.isValid())
0044     return;
0045 
0046   std::vector<bool> selResults;
0047   for (auto& obj : *handle) {
0048     selResults.push_back(selection_(obj));
0049   }
0050   auto valMap = std::make_unique<edm::ValueMap<bool> >();
0051   edm::ValueMap<bool>::Filler filler(*valMap);
0052   filler.insert(handle, selResults.begin(), selResults.end());
0053   filler.fill();
0054   event.put(std::move(valMap));
0055 }
0056 
0057 #include "FWCore/Framework/interface/MakerMacros.h"
0058 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0059 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
0060 using HLTDQMGsfEleSelector = HLTDQMObjSelector<reco::GsfElectron, reco::GsfElectronCollection>;
0061 DEFINE_FWK_MODULE(HLTDQMGsfEleSelector);