CMS 3D CMS Logo

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 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
31 
39 
42 
45 
46 //
47 // class declaration
48 //
49 
51  public:
52  explicit HcalHitSelection(const edm::ParameterSet&);
53  ~HcalHitSelection() override;
54 
55  private:
56  void produce(edm::Event&, const edm::EventSetup&) override;
57 
62  std::vector<edm::EDGetTokenT<DetIdCollection> > toks_did_;
64  std::vector<edm::InputTag> interestingDetIdCollections;
66 
67  //hcal severity ES
70  std::set<DetId> toBeKept;
71  template <typename CollectionType> void skim( const edm::Handle<CollectionType> & input, CollectionType & output,int severityThreshold=0) const;
72 
73  // ----------member data ---------------------------
74 };
75 
76 template <class CollectionType> void HcalHitSelection::skim( const edm::Handle<CollectionType> & input, CollectionType & output,int severityThreshold) const {
77  output.reserve(input->size());
78  typename CollectionType::const_iterator begin=input->begin();
79  typename CollectionType::const_iterator end=input->end();
80  typename CollectionType::const_iterator hit=begin;
81 
82  for (;hit!=end;++hit){
83  // edm::LogError("HcalHitSelection")<<"the hit pointer is"<<&(*hit);
84  HcalDetId id = hit->detid();
85  if (theHcalTopology_->getMergePositionFlag() && id.subdet() == HcalEndcap) {
86  id = theHcalTopology_->idFront(id);
87  }
88  const uint32_t & recHitFlag = hit->flags();
89  // edm::LogError("HcalHitSelection")<<"the hit id and flag are "<<id.rawId()<<" "<<recHitFlag;
90 
91  const uint32_t & dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
92  int severityLevel = theHcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
93  //anything that is not "good" goes in
94  if (severityLevel>severityThreshold){
95  output.push_back(*hit);
96  }else{
97  //chek on the detid list
98  if (toBeKept.find(id)!=toBeKept.end())
99  output.push_back(*hit);
100  }
101  }
102 }
103 
104 //
105 // constants, enums and typedefs
106 //
107 
108 
109 //
110 // static data member definitions
111 //
112 
113 //
114 // constructors and destructor
115 //
117 {
118  hbheTag=iConfig.getParameter<edm::InputTag>("hbheTag");
119  hfTag=iConfig.getParameter<edm::InputTag>("hfTag");
120  hoTag=iConfig.getParameter<edm::InputTag>("hoTag");
121 
122  // register for data access
123  tok_hbhe_ = consumes<HBHERecHitCollection>(hbheTag);
124  tok_hf_ = consumes<HFRecHitCollection>(hfTag);
125  tok_ho_ = consumes<HORecHitCollection>(hoTag);
126 
127  interestingDetIdCollections = iConfig.getParameter< std::vector<edm::InputTag> >("interestingDetIds");
128 
129  const unsigned nLabels = interestingDetIdCollections.size();
130  for ( unsigned i=0; i != nLabels; i++ )
131  toks_did_.push_back(consumes<DetIdCollection>(interestingDetIdCollections[i]));
132 
133  hoSeverityLevel=iConfig.getParameter<int>("hoSeverityLevel");
134 
135  produces<HBHERecHitCollection>(hbheTag.label());
136  produces<HFRecHitCollection>(hfTag.label());
137  produces<HORecHitCollection>(hoTag.label());
138 
139 }
140 
141 
143 {
144 
145  // do anything here that needs to be done at desctruction time
146  // (e.g. close files, deallocate resources etc.)
147 
148 }
149 
150 
151 //
152 // member functions
153 //
154 
155 // ------------ method called to produce the data ------------
156 void
158 {
159  iSetup.get<HcalChannelQualityRcd>().get("withTopo", theHcalChStatus);
162  iSetup.get<HcalRecNumberingRecord>().get(topo);
163  theHcalTopology_ = topo.product();
164 
168 
169  iEvent.getByToken(tok_hbhe_,hbhe);
170  iEvent.getByToken(tok_hf_,hf);
171  iEvent.getByToken(tok_ho_,ho);
172 
173  toBeKept.clear();
175  for( unsigned int t = 0; t < toks_did_.size(); ++t )
176  {
177  iEvent.getByToken(toks_did_[t],detId);
178  if (!detId.isValid()){
179  edm::LogError("MissingInput")<<"the collection of interesting detIds:"<<interestingDetIdCollections[t]<<" is not found.";
180  continue;
181  }
182  toBeKept.insert(detId->begin(),detId->end());
183  }
184 
185  auto hbhe_out = std::make_unique<HBHERecHitCollection>();
186  skim(hbhe,*hbhe_out);
187  iEvent.put(std::move(hbhe_out),hbheTag.label());
188 
189  auto hf_out = std::make_unique<HFRecHitCollection>();
190  skim(hf,*hf_out);
191  iEvent.put(std::move(hf_out),hfTag.label());
192 
193  auto ho_out = std::make_unique<HORecHitCollection>();
194  skim(ho,*ho_out,hoSeverityLevel);
195  iEvent.put(std::move(ho_out),hoTag.label());
196 
197 }
198 
199 //define this as a plug-in
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
std::vector< edm::InputTag > interestingDetIdCollections
~HcalHitSelection() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::set< DetId > toBeKept
const_iterator end() const
Definition: EDCollection.h:153
edm::ESHandle< HcalChannelQuality > theHcalChStatus
edm::EDGetTokenT< HORecHitCollection > tok_ho_
const Item * getValues(DetId fId, bool throwOnFail=true) const
std::vector< edm::EDGetTokenT< DetIdCollection > > toks_did_
edm::InputTag hoTag
bool getMergePositionFlag() const
Definition: HcalTopology.h:171
edm::InputTag hfTag
static std::string const input
Definition: EdmProvDump.cc:48
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
void skim(const edm::Handle< CollectionType > &input, CollectionType &output, int severityThreshold=0) const
#define end
Definition: vmac.h:39
bool isValid() const
Definition: HandleBase.h:74
void produce(edm::Event &, const edm::EventSetup &) override
const_iterator begin() const
Definition: EDCollection.h:146
edm::EDGetTokenT< HFRecHitCollection > tok_hf_
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
edm::ESHandle< HcalSeverityLevelComputer > theHcalSevLvlComputer
#define begin
Definition: vmac.h:32
T get() const
Definition: EventSetup.h:71
HcalHitSelection(const edm::ParameterSet &)
HcalDetId idFront(const HcalDetId &id) const
Definition: HcalTopology.h:177
edm::InputTag hbheTag
uint32_t getValue() const
T const * product() const
Definition: ESHandle.h:86
def move(src, dest)
Definition: eostools.py:511