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$
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 
58  std::vector<edm::InputTag> interestingDetIdCollections;
59 
60  //hcal severity ES
63  std::set<DetId> toBeKept;
64  template <typename CollectionType> void skim( const edm::Handle<CollectionType> & input, std::auto_ptr<CollectionType> & output);
65 
66  // ----------member data ---------------------------
67 };
68 
69 template <class CollectionType> void HcalHitSelection::skim( const edm::Handle<CollectionType> & input, std::auto_ptr<CollectionType> & output){
70  output->reserve(input->size());
71  typename CollectionType::const_iterator begin=input->begin();
72  typename CollectionType::const_iterator end=input->end();
73  typename CollectionType::const_iterator hit=begin;
74 
75  for (;hit!=end;++hit){
76  // edm::LogError("HcalHitSelection")<<"the hit pointer is"<<&(*hit);
77  const DetId & id = hit->detid();
78  const uint32_t & recHitFlag = hit->flags();
79  // edm::LogError("HcalHitSelection")<<"the hit id and flag are "<<id.rawId()<<" "<<recHitFlag;
80 
81  const uint32_t & dbStatusFlag = theHcalChStatus->getValues(id)->getValue();
82  int severityLevel = theHcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
83  //anything that is not "good" goes in
84  if (severityLevel!=0){
85  output->push_back(*hit);
86  }else{
87  //chek on the detid list
88  if (toBeKept.find(id)!=toBeKept.end())
89  output->push_back(*hit);
90  }
91  }
92 }
93 
94 //
95 // constants, enums and typedefs
96 //
97 
98 
99 //
100 // static data member definitions
101 //
102 
103 //
104 // constructors and destructor
105 //
107 {
108  hbheTag=iConfig.getParameter<edm::InputTag>("hbheTag");
109  hfTag=iConfig.getParameter<edm::InputTag>("hfTag");
110  hoTag=iConfig.getParameter<edm::InputTag>("hoTag");
111 
112  interestingDetIdCollections = iConfig.getParameter< std::vector<edm::InputTag> >("interestingDetIds");
113 
114  produces<HBHERecHitCollection>(hbheTag.label());
115  produces<HFRecHitCollection>(hfTag.label());
116  produces<HORecHitCollection>(hoTag.label());
117 
118 }
119 
120 
122 {
123 
124  // do anything here that needs to be done at desctruction time
125  // (e.g. close files, deallocate resources etc.)
126 
127 }
128 
129 
130 //
131 // member functions
132 //
133 
134 // ------------ method called to produce the data ------------
135 void
137 {
140 
144 
145  iEvent.getByLabel(hbheTag,hbhe);
146  iEvent.getByLabel(hfTag,hf);
147  iEvent.getByLabel(hoTag,ho);
148 
149  toBeKept.clear();
151  for( unsigned int t = 0; t < interestingDetIdCollections.size(); ++t )
152  {
154  if (!detId.isValid())
155  continue;
156  toBeKept.insert(detId->begin(),detId->end());
157  }
158 
159  std::auto_ptr<HBHERecHitCollection> hbhe_out(new HBHERecHitCollection());
160  skim(hbhe,hbhe_out);
161  iEvent.put(hbhe_out,hbheTag.label());
162 
163  std::auto_ptr<HFRecHitCollection> hf_out(new HFRecHitCollection());
164  skim(hf,hf_out);
165  iEvent.put(hf_out,hfTag.label());
166 
167  std::auto_ptr<HORecHitCollection> ho_out(new HORecHitCollection());
168  skim(ho,ho_out);
169  iEvent.put(ho_out,hoTag.label());
170 
171 
172 
173 }
174 
175 // ------------ method called once each job just before starting event loop ------------
176 void
178 {
179 }
180 
181 // ------------ method called once each job just after ending the event loop ------------
182 void
184 }
185 
186 //define this as a plug-in
void skim(const edm::Handle< CollectionType > &input, std::auto_ptr< CollectionType > &output)
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:84
#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:359
tuple input
Definition: collect_tpl.py:10
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 &)
edm::InputTag hbheTag
virtual void beginJob()
edm::SortedCollection< HORecHit > HORecHitCollection