CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
31 
35 
38 
41 
42 //
43 // class declaration
44 //
45 
47  public:
48  explicit HcalHitSelection(const edm::ParameterSet&);
50 
51  private:
52  virtual void beginJob() override ;
53  virtual void produce(edm::Event&, const edm::EventSetup&) override;
54  virtual void endJob() override ;
55 
60  std::vector<edm::EDGetTokenT<DetIdCollection> > toks_did_;
62  std::vector<edm::InputTag> interestingDetIdCollections;
63 
64  //hcal severity ES
67  std::set<DetId> toBeKept;
68  template <typename CollectionType> void skim( const edm::Handle<CollectionType> & input, std::auto_ptr<CollectionType> & output,int severityThreshold=0);
69 
70  // ----------member data ---------------------------
71 };
72 
73 template <class CollectionType> void HcalHitSelection::skim( const edm::Handle<CollectionType> & input, std::auto_ptr<CollectionType> & output,int severityThreshold){
74  output->reserve(input->size());
75  typename CollectionType::const_iterator begin=input->begin();
76  typename CollectionType::const_iterator end=input->end();
77  typename CollectionType::const_iterator hit=begin;
78 
79  for (;hit!=end;++hit){
80  // edm::LogError("HcalHitSelection")<<"the hit pointer is"<<&(*hit);
81  const DetId & id = hit->detid();
82  const uint32_t & recHitFlag = hit->flags();
83  // edm::LogError("HcalHitSelection")<<"the hit id and flag are "<<id.rawId()<<" "<<recHitFlag;
84 
85  const uint32_t & dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
86  int severityLevel = theHcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
87  //anything that is not "good" goes in
88  if (severityLevel>severityThreshold){
89  output->push_back(*hit);
90  }else{
91  //chek on the detid list
92  if (toBeKept.find(id)!=toBeKept.end())
93  output->push_back(*hit);
94  }
95  }
96 }
97 
98 //
99 // constants, enums and typedefs
100 //
101 
102 
103 //
104 // static data member definitions
105 //
106 
107 //
108 // constructors and destructor
109 //
111 {
112  hbheTag=iConfig.getParameter<edm::InputTag>("hbheTag");
113  hfTag=iConfig.getParameter<edm::InputTag>("hfTag");
114  hoTag=iConfig.getParameter<edm::InputTag>("hoTag");
115 
116  // register for data access
117  tok_hbhe_ = consumes<HBHERecHitCollection>(hbheTag);
118  tok_hf_ = consumes<HFRecHitCollection>(hfTag);
119  tok_ho_ = consumes<HORecHitCollection>(hoTag);
120 
121  interestingDetIdCollections = iConfig.getParameter< std::vector<edm::InputTag> >("interestingDetIds");
122 
123  const unsigned nLabels = interestingDetIdCollections.size();
124  for ( unsigned i=0; i != nLabels; i++ )
125  toks_did_.push_back(consumes<DetIdCollection>(interestingDetIdCollections[i]));
126 
127  hoSeverityLevel=iConfig.getParameter<int>("hoSeverityLevel");
128 
129  produces<HBHERecHitCollection>(hbheTag.label());
130  produces<HFRecHitCollection>(hfTag.label());
131  produces<HORecHitCollection>(hoTag.label());
132 
133 }
134 
135 
137 {
138 
139  // do anything here that needs to be done at desctruction time
140  // (e.g. close files, deallocate resources etc.)
141 
142 }
143 
144 
145 //
146 // member functions
147 //
148 
149 // ------------ method called to produce the data ------------
150 void
152 {
155 
159 
160  iEvent.getByToken(tok_hbhe_,hbhe);
161  iEvent.getByToken(tok_hf_,hf);
162  iEvent.getByToken(tok_ho_,ho);
163 
164  toBeKept.clear();
166  for( unsigned int t = 0; t < toks_did_.size(); ++t )
167  {
168  iEvent.getByToken(toks_did_[t],detId);
169  if (!detId.isValid()){
170  edm::LogError("MissingInput")<<"the collection of interesting detIds:"<<interestingDetIdCollections[t]<<" is not found.";
171  continue;
172  }
173  toBeKept.insert(detId->begin(),detId->end());
174  }
175 
176  std::auto_ptr<HBHERecHitCollection> hbhe_out(new HBHERecHitCollection());
177  skim(hbhe,hbhe_out);
178  iEvent.put(hbhe_out,hbheTag.label());
179 
180  std::auto_ptr<HFRecHitCollection> hf_out(new HFRecHitCollection());
181  skim(hf,hf_out);
182  iEvent.put(hf_out,hfTag.label());
183 
184  std::auto_ptr<HORecHitCollection> ho_out(new HORecHitCollection());
185  skim(ho,ho_out,hoSeverityLevel);
186  iEvent.put(ho_out,hoTag.label());
187 
188 }
189 
190 // ------------ method called once each job just before starting event loop ------------
191 void
193 {
194 }
195 
196 // ------------ method called once each job just after ending the event loop ------------
197 void
199 }
200 
201 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< edm::InputTag > interestingDetIdCollections
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::set< DetId > toBeKept
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::ESHandle< HcalChannelQuality > theHcalChStatus
edm::EDGetTokenT< HORecHitCollection > tok_ho_
std::vector< edm::EDGetTokenT< DetIdCollection > > toks_did_
edm::InputTag hoTag
edm::InputTag hfTag
static std::string const input
Definition: EdmProvDump.cc:44
int iEvent
Definition: GenABIO.cc:243
virtual void beginJob() override
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
virtual void endJob() override
#define end
Definition: vmac.h:37
bool isValid() const
Definition: HandleBase.h:76
edm::SortedCollection< HBHERecHit > HBHERecHitCollection
Definition: DetId.h:18
virtual void produce(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< HFRecHitCollection > tok_hf_
const T & get() const
Definition: EventSetup.h:55
std::string const & label() const
Definition: InputTag.h:42
edm::ESHandle< HcalSeverityLevelComputer > theHcalSevLvlComputer
#define begin
Definition: vmac.h:30
edm::SortedCollection< HFRecHit > HFRecHitCollection
HcalHitSelection(const edm::ParameterSet &)
void skim(const edm::Handle< CollectionType > &input, std::auto_ptr< CollectionType > &output, int severityThreshold=0)
edm::InputTag hbheTag
edm::SortedCollection< HORecHit > HORecHitCollection