CMS 3D CMS Logo

HLTHcalPhiSymFilter.cc
Go to the documentation of this file.
1 #include "HLTHcalPhiSymFilter.h"
4 
6 {
7  HBHEHits_ = iConfig.getParameter< edm::InputTag > ("HBHEHitCollection");
8  HOHits_ = iConfig.getParameter< edm::InputTag > ("HOHitCollection");
9  HFHits_= iConfig.getParameter<edm::InputTag>("HFHitCollection");
10  phiSymHBHEHits_ = iConfig.getParameter<std::string > ("phiSymHBHEHitCollection");
11  phiSymHOHits_ = iConfig.getParameter<std::string > ("phiSymHOHitCollection");
12  phiSymHFHits_ = iConfig.getParameter<std::string > ("phiSymHFHitCollection");
13 
14  eCut_HB_ = iConfig.getParameter< double > ("eCut_HB");
15  eCut_HE_ = iConfig.getParameter< double > ("eCut_HE");
16  eCut_HO_ = iConfig.getParameter<double>("eCut_HO");
17  eCut_HF_ = iConfig.getParameter<double>("eCut_HF");
18 
19  HBHEHitsToken_ = consumes<HBHERecHitCollection>(HBHEHits_);
20  HOHitsToken_ = consumes<HORecHitCollection>(HOHits_);
21  HFHitsToken_ = consumes<HFRecHitCollection>(HFHits_);
22 
23  //register your products
24  produces< HBHERecHitCollection >(phiSymHBHEHits_);
25  produces< HORecHitCollection >(phiSymHOHits_);
26  produces< HFRecHitCollection >(phiSymHFHits_);
27 }
28 
29 
31 
32 void
36  desc.add<edm::InputTag>("HBHEHitCollection",edm::InputTag("hbhereco"));
37  desc.add<edm::InputTag>("HOHitCollection",edm::InputTag("horeco"));
38  desc.add<edm::InputTag>("HFHitCollection",edm::InputTag("hfreco"));
39  desc.add<double>("eCut_HE",-10.);
40  desc.add<double>("eCut_HF",-10.);
41  desc.add<double>("eCut_HB",-10.);
42  desc.add<double>("eCut_HO",-10.);
43  desc.add<std::string>("phiSymHOHitCollection","phiSymHcalRecHitsHO");
44  desc.add<std::string>("phiSymHBHEHitCollection","phiSymHcalRecHitsHBHE");
45  desc.add<std::string>("phiSymHFHitCollection","phiSymHcalRecHitsHF");
46  descriptions.add("alCaHcalPhiSymStream",desc);
47 }
48 
49 // ------------ method called to produce the data ------------
50 bool
52 {
56 
57  iEvent.getByToken(HBHEHitsToken_,HBHERecHitsH);
58  iEvent.getByToken(HOHitsToken_,HORecHitsH);
59  iEvent.getByToken(HFHitsToken_,HFRecHitsH);
60 
61  //Create empty output collections
62  std::unique_ptr< HBHERecHitCollection > phiSymHBHERecHitCollection( new HBHERecHitCollection );
63  std::unique_ptr< HORecHitCollection > phiSymHORecHitCollection( new HORecHitCollection );
64  std::unique_ptr< HFRecHitCollection > phiSymHFRecHitCollection( new HFRecHitCollection );
65 
66  //Select interesting HBHERecHits
67  for (auto const & it : *HBHERecHitsH) {
68  if (it.energy()>eCut_HB_&&it.id().subdet()==1) {
69  phiSymHBHERecHitCollection->push_back(it);
70  }
71  if (it.energy()>eCut_HE_&&it.id().subdet()==2) {
72  phiSymHBHERecHitCollection->push_back(it);
73  }
74 
75  }
76 
77  //Select interesting HORecHits
78  for (auto const & it : *HORecHitsH) {
79  if (it.energy()>eCut_HO_&&it.id().subdet()==3) {
80  phiSymHORecHitCollection->push_back(it);
81  }
82  }
83 
84  //Select interesting HFRecHits
85  for (auto const & it : *HFRecHitsH) {
86  if (it.energy()>eCut_HF_&&it.id().subdet()==4) {
87  phiSymHFRecHitCollection->push_back(it);
88  }
89  }
90 
91  if ((phiSymHBHERecHitCollection->empty() ) && (phiSymHORecHitCollection->empty()) && (phiSymHFRecHitCollection->empty())) return false;
92 
93  //Put selected information in the event
94  if (!phiSymHBHERecHitCollection->empty()) iEvent.put(std::move(phiSymHBHERecHitCollection), phiSymHBHEHits_);
95  if (!phiSymHORecHitCollection->empty()) iEvent.put(std::move(phiSymHORecHitCollection), phiSymHOHits_);
96  if (!phiSymHFRecHitCollection->empty()) iEvent.put(std::move(phiSymHFRecHitCollection), phiSymHFHits_);
97 
98  return true;
99 }
100 
101 // declare this class as a framework plugin
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
edm::EDGetTokenT< HFRecHitCollection > HFHitsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
~HLTHcalPhiSymFilter() override
edm::EDGetTokenT< HORecHitCollection > HOHitsToken_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool hltFilter(edm::Event &, const edm::EventSetup &, trigger::TriggerFilterObjectWithRefs &filterproduct) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void makeHLTFilterDescription(edm::ParameterSetDescription &desc)
Definition: HLTFilter.cc:29
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLTHcalPhiSymFilter(const edm::ParameterSet &)
edm::EDGetTokenT< HBHERecHitCollection > HBHEHitsToken_
def move(src, dest)
Definition: eostools.py:511