CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HcalHitSelection.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalHitSelection
4 // Class: HcalHitSelection
5 //
13 //
14 // Original Author: Jean-Roch Vlimant,40 3-A28,+41227671209,
15 // Created: Thu Nov 4 22:17:56 CET 2010
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
25 
28 
30 
38 
41 
44 
45 //
46 // class declaration
47 //
48 
50 public:
51  explicit HcalHitSelection(const edm::ParameterSet&);
52  ~HcalHitSelection() override;
53 
54 private:
55  void produce(edm::Event&, const edm::EventSetup&) override;
56 
61  std::vector<edm::EDGetTokenT<DetIdCollection> > toks_did_;
63  std::vector<edm::InputTag> interestingDetIdCollections;
65 
66  //hcal severity ES
69  std::set<DetId> toBeKept;
70  template <typename CollectionType>
71  void skim(const edm::Handle<CollectionType>& input, CollectionType& output, int severityThreshold = 0) const;
72 
73  // ES tokens
77 };
78 
79 template <class CollectionType>
81  CollectionType& output,
82  int severityThreshold) const {
83  output.reserve(input->size());
84  typename CollectionType::const_iterator begin = input->begin();
85  typename CollectionType::const_iterator end = input->end();
86  typename CollectionType::const_iterator hit = begin;
87 
88  for (; hit != end; ++hit) {
89  // edm::LogError("HcalHitSelection")<<"the hit pointer is"<<&(*hit);
90  HcalDetId id = hit->detid();
91  if (theHcalTopology_->getMergePositionFlag() && id.subdet() == HcalEndcap) {
92  id = theHcalTopology_->idFront(id);
93  }
94  const uint32_t& recHitFlag = hit->flags();
95  // edm::LogError("HcalHitSelection")<<"the hit id and flag are "<<id.rawId()<<" "<<recHitFlag;
96 
97  const uint32_t& dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
98  int severityLevel = theHcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
99  //anything that is not "good" goes in
100  if (severityLevel > severityThreshold) {
101  output.push_back(*hit);
102  } else {
103  //chek on the detid list
104  if (toBeKept.find(id) != toBeKept.end())
105  output.push_back(*hit);
106  }
107  }
108 }
109 
110 //
111 // constants, enums and typedefs
112 //
113 
114 //
115 // static data member definitions
116 //
117 
118 //
119 // constructors and destructor
120 //
122  : hbheTag(iConfig.getParameter<edm::InputTag>("hbheTag")),
123  hoTag(iConfig.getParameter<edm::InputTag>("hoTag")),
124  hfTag(iConfig.getParameter<edm::InputTag>("hfTag")),
125  theHcalTopology_(nullptr),
126  theHcalChStatus(nullptr),
127  theHcalSevLvlComputer(nullptr) {
128  // register for data access
129  tok_hbhe_ = consumes<HBHERecHitCollection>(hbheTag);
130  tok_hf_ = consumes<HFRecHitCollection>(hfTag);
131  tok_ho_ = consumes<HORecHitCollection>(hoTag);
132 
133  interestingDetIdCollections = iConfig.getParameter<std::vector<edm::InputTag> >("interestingDetIds");
134 
135  const unsigned nLabels = interestingDetIdCollections.size();
136  for (unsigned i = 0; i != nLabels; i++)
137  toks_did_.push_back(consumes<DetIdCollection>(interestingDetIdCollections[i]));
138 
139  hoSeverityLevel = iConfig.getParameter<int>("hoSeverityLevel");
140 
141  produces<HBHERecHitCollection>(hbheTag.label());
142  produces<HFRecHitCollection>(hfTag.label());
143  produces<HORecHitCollection>(hoTag.label());
144 
145  // ES tokens
146  htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord>();
147  qualToken_ = esConsumes<HcalChannelQuality, HcalChannelQualityRcd>(edm::ESInputTag("", "withTopo"));
148  sevToken_ = esConsumes<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd>();
149 }
150 
152  // do anything here that needs to be done at desctruction time
153  // (e.g. close files, deallocate resources etc.)
154 }
155 
156 //
157 // member functions
158 //
159 
160 // ------------ method called to produce the data ------------
165 
169 
170  iEvent.getByToken(tok_hbhe_, hbhe);
171  iEvent.getByToken(tok_hf_, hf);
172  iEvent.getByToken(tok_ho_, ho);
173 
174  toBeKept.clear();
176  for (unsigned int t = 0; t < toks_did_.size(); ++t) {
177  iEvent.getByToken(toks_did_[t], detId);
178  if (!detId.isValid()) {
179  edm::LogError("MissingInput") << "the collection of interesting detIds:" << interestingDetIdCollections[t]
180  << " is not found.";
181  continue;
182  }
183  toBeKept.insert(detId->begin(), detId->end());
184  }
185 
186  auto hbhe_out = std::make_unique<HBHERecHitCollection>();
187  skim(hbhe, *hbhe_out);
188  iEvent.put(std::move(hbhe_out), hbheTag.label());
189 
190  auto hf_out = std::make_unique<HFRecHitCollection>();
191  skim(hf, *hf_out);
192  iEvent.put(std::move(hf_out), hfTag.label());
193 
194  auto ho_out = std::make_unique<HORecHitCollection>();
195  skim(ho, *ho_out, hoSeverityLevel);
196  iEvent.put(std::move(ho_out), hoTag.label());
197 }
198 
199 //define this as a plug-in
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< edm::InputTag > interestingDetIdCollections
~HcalHitSelection() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::set< DetId > toBeKept
const HcalChannelQuality * theHcalChStatus
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< HORecHitCollection > tok_ho_
const Item * getValues(DetId fId, bool throwOnFail=true) const
std::vector< edm::EDGetTokenT< DetIdCollection > > toks_did_
Log< level::Error, false > LogError
edm::ESGetToken< HcalChannelQuality, HcalChannelQualityRcd > qualToken_
edm::InputTag hoTag
bool getMergePositionFlag() const
Definition: HcalTopology.h:167
edm::InputTag hfTag
static std::string const input
Definition: EdmProvDump.cc:47
bool getData(T &iHolder) const
Definition: EventSetup.h:128
const HcalSeverityLevelComputer * theHcalSevLvlComputer
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
void skim(const edm::Handle< CollectionType > &input, CollectionType &output, int severityThreshold=0) const
bool isValid() const
Definition: HandleBase.h:70
edm::ESGetToken< HcalTopology, HcalRecNumberingRecord > htopoToken_
void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< HFRecHitCollection > tok_hf_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int getSeverityLevel(const DetId &myid, const uint32_t &myflag, const uint32_t &mystatus) const
const HcalTopology * theHcalTopology_
std::string const & label() const
Definition: InputTag.h:36
string end
Definition: dataset.py:937
HcalHitSelection(const edm::ParameterSet &)
HcalDetId idFront(const HcalDetId &id) const
Definition: HcalTopology.h:170
edm::InputTag hbheTag
uint32_t getValue() const
edm::ESGetToken< HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd > sevToken_