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 // $Id: HcalHitSelection.cc,v 1.4 2011/03/24 19:41:28 vlimant Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
27 
30 
32 
36 
39 
42 
43 //
44 // class declaration
45 //
46 
48  public:
49  explicit HcalHitSelection(const edm::ParameterSet&);
51 
52  private:
53  virtual void beginJob() ;
54  virtual void produce(edm::Event&, const edm::EventSetup&);
55  virtual void endJob() ;
56 
59  std::vector<edm::InputTag> interestingDetIdCollections;
60 
61  //hcal severity ES
64  std::set<DetId> toBeKept;
65  template <typename CollectionType> void skim( const edm::Handle<CollectionType> & input, std::auto_ptr<CollectionType> & output,int severityThreshold=0);
66 
67  // ----------member data ---------------------------
68 };
69 
70 template <class CollectionType> void HcalHitSelection::skim( const edm::Handle<CollectionType> & input, std::auto_ptr<CollectionType> & output,int severityThreshold){
71  output->reserve(input->size());
72  typename CollectionType::const_iterator begin=input->begin();
73  typename CollectionType::const_iterator end=input->end();
74  typename CollectionType::const_iterator hit=begin;
75 
76  for (;hit!=end;++hit){
77  // edm::LogError("HcalHitSelection")<<"the hit pointer is"<<&(*hit);
78  const DetId & id = hit->detid();
79  const uint32_t & recHitFlag = hit->flags();
80  // edm::LogError("HcalHitSelection")<<"the hit id and flag are "<<id.rawId()<<" "<<recHitFlag;
81 
82  const uint32_t & dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
83  int severityLevel = theHcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
84  //anything that is not "good" goes in
85  if (severityLevel>severityThreshold){
86  output->push_back(*hit);
87  }else{
88  //chek on the detid list
89  if (toBeKept.find(id)!=toBeKept.end())
90  output->push_back(*hit);
91  }
92  }
93 }
94 
95 //
96 // constants, enums and typedefs
97 //
98 
99 
100 //
101 // static data member definitions
102 //
103 
104 //
105 // constructors and destructor
106 //
108 {
109  hbheTag=iConfig.getParameter<edm::InputTag>("hbheTag");
110  hfTag=iConfig.getParameter<edm::InputTag>("hfTag");
111  hoTag=iConfig.getParameter<edm::InputTag>("hoTag");
112 
113  interestingDetIdCollections = iConfig.getParameter< std::vector<edm::InputTag> >("interestingDetIds");
114 
115  hoSeverityLevel=iConfig.getParameter<int>("hoSeverityLevel");
116 
117  produces<HBHERecHitCollection>(hbheTag.label());
118  produces<HFRecHitCollection>(hfTag.label());
119  produces<HORecHitCollection>(hoTag.label());
120 
121 }
122 
123 
125 {
126 
127  // do anything here that needs to be done at desctruction time
128  // (e.g. close files, deallocate resources etc.)
129 
130 }
131 
132 
133 //
134 // member functions
135 //
136 
137 // ------------ method called to produce the data ------------
138 void
140 {
143 
147 
148  iEvent.getByLabel(hbheTag,hbhe);
149  iEvent.getByLabel(hfTag,hf);
150  iEvent.getByLabel(hoTag,ho);
151 
152  toBeKept.clear();
154  for( unsigned int t = 0; t < interestingDetIdCollections.size(); ++t )
155  {
157  if (!detId.isValid()){
158  edm::LogError("MissingInput")<<"the collection of interesting detIds:"<<interestingDetIdCollections[t]<<" is not found.";
159  continue;
160  }
161  toBeKept.insert(detId->begin(),detId->end());
162  }
163 
164  std::auto_ptr<HBHERecHitCollection> hbhe_out(new HBHERecHitCollection());
165  skim(hbhe,hbhe_out);
166  iEvent.put(hbhe_out,hbheTag.label());
167 
168  std::auto_ptr<HFRecHitCollection> hf_out(new HFRecHitCollection());
169  skim(hf,hf_out);
170  iEvent.put(hf_out,hfTag.label());
171 
172  std::auto_ptr<HORecHitCollection> ho_out(new HORecHitCollection());
173  skim(ho,ho_out,hoSeverityLevel);
174  iEvent.put(ho_out,hoTag.label());
175 
176 }
177 
178 // ------------ method called once each job just before starting event loop ------------
179 void
181 {
182 }
183 
184 // ------------ method called once each job just after ending the event loop ------------
185 void
187 }
188 
189 //define this as a plug-in
T getParameter(std::string const &) const
std::vector< edm::InputTag > interestingDetIdCollections
std::set< DetId > toBeKept
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::ESHandle< HcalChannelQuality > theHcalChStatus
edm::InputTag hoTag
edm::InputTag hfTag
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
#define end
Definition: vmac.h:38
bool isValid() const
Definition: HandleBase.h:76
virtual void produce(edm::Event &, const edm::EventSetup &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
edm::SortedCollection< HBHERecHit > HBHERecHitCollection
Definition: DetId.h:20
const T & get() const
Definition: EventSetup.h:55
virtual void endJob()
std::string const & label() const
Definition: InputTag.h:25
edm::ESHandle< HcalSeverityLevelComputer > theHcalSevLvlComputer
#define begin
Definition: vmac.h:31
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
virtual void beginJob()
edm::SortedCollection< HORecHit > HORecHitCollection