CMS 3D CMS Logo

HFPhase1Reconstructor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoLocalCalo/HcalRecProducers
4 // Class: HFPhase1Reconstructor
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Thu, 25 May 2016 00:17:51 GMT
16 //
17 //
18 
19 // system include files
20 
21 // user include files
28 
31 
33 
36 
39 
40 // Noise cleanup algos
44 
45 // Parser for Phase 1 HF reco algorithms
47 
48 // Fetcher for reco algorithm data
50 
51 //
52 // class declaration
53 //
55 {
56 public:
58  ~HFPhase1Reconstructor() override;
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
61 
62 private:
63  void beginRun(const edm::Run&, const edm::EventSetup&) override;
64  void produce(edm::Event&, const edm::EventSetup&) override;
65 
66  // Configuration parameters
72 
73  // Other members
75  std::unique_ptr<AbsHFPhase1Algo> reco_;
76  std::unique_ptr<AbsHcalAlgoData> recoConfig_;
77 
78  // Noise cleanup algos
79  std::unique_ptr<HcalHF_S9S1algorithm> hfS9S1_;
80  std::unique_ptr<HcalHF_S9S1algorithm> hfS8S1_;
81  std::unique_ptr<HcalHF_PETalgorithm> hfPET_;
82  std::unique_ptr<HFStripFilter> hfStripFilter_;
83 };
84 
85 //
86 // constructors and destructor
87 //
89  : algoConfigClass_(conf.getParameter<std::string>("algoConfigClass")),
90  setNoiseFlags_(conf.getParameter<bool>("setNoiseFlags")),
91  runHFStripFilter_(conf.getParameter<bool>("runHFStripFilter")),
92  useChannelQualityFromDB_(conf.getParameter<bool>("useChannelQualityFromDB")),
93  checkChannelQualityForDepth3and4_(conf.getParameter<bool>("checkChannelQualityForDepth3and4")),
94  reco_(parseHFPhase1AlgoDescription(conf.getParameter<edm::ParameterSet>("algorithm")))
95 {
96  // Check that the reco algorithm has been successfully configured
97  if (!reco_.get())
98  throw cms::Exception("HFPhase1BadConfig")
99  << "Invalid HFPhase1Reconstructor algorithm configuration"
100  << std::endl;
101 
102  // Configure the noise cleanup algorithms
103  if (setNoiseFlags_)
104  {
105  const edm::ParameterSet& psS9S1 = conf.getParameter<edm::ParameterSet>("S9S1stat");
106  hfS9S1_ = std::make_unique<HcalHF_S9S1algorithm>(
107  psS9S1.getParameter<std::vector<double> >("short_optimumSlope"),
108  psS9S1.getParameter<std::vector<double> >("shortEnergyParams"),
109  psS9S1.getParameter<std::vector<double> >("shortETParams"),
110  psS9S1.getParameter<std::vector<double> >("long_optimumSlope"),
111  psS9S1.getParameter<std::vector<double> >("longEnergyParams"),
112  psS9S1.getParameter<std::vector<double> >("longETParams"),
113  psS9S1.getParameter<int>("HcalAcceptSeverityLevel"),
114  psS9S1.getParameter<bool>("isS8S1") );
115 
116  const edm::ParameterSet& psS8S1 = conf.getParameter<edm::ParameterSet>("S8S1stat");
117  hfS8S1_ = std::make_unique<HcalHF_S9S1algorithm>(
118  psS8S1.getParameter<std::vector<double> >("short_optimumSlope"),
119  psS8S1.getParameter<std::vector<double> >("shortEnergyParams"),
120  psS8S1.getParameter<std::vector<double> >("shortETParams"),
121  psS8S1.getParameter<std::vector<double> >("long_optimumSlope"),
122  psS8S1.getParameter<std::vector<double> >("longEnergyParams"),
123  psS8S1.getParameter<std::vector<double> >("longETParams"),
124  psS8S1.getParameter<int>("HcalAcceptSeverityLevel"),
125  psS8S1.getParameter<bool>("isS8S1") );
126 
127  const edm::ParameterSet& psPET = conf.getParameter<edm::ParameterSet>("PETstat");
128  hfPET_ = std::make_unique<HcalHF_PETalgorithm>(
129  psPET.getParameter<std::vector<double> >("short_R"),
130  psPET.getParameter<std::vector<double> >("shortEnergyParams"),
131  psPET.getParameter<std::vector<double> >("shortETParams"),
132  psPET.getParameter<std::vector<double> >("long_R"),
133  psPET.getParameter<std::vector<double> >("longEnergyParams"),
134  psPET.getParameter<std::vector<double> >("longETParams"),
135  psPET.getParameter<int>("HcalAcceptSeverityLevel"),
136  psPET.getParameter<std::vector<double> >("short_R_29"),
137  psPET.getParameter<std::vector<double> >("long_R_29") );
138 
139  // Configure HFStripFilter
140  if (runHFStripFilter_)
141  {
142  const edm::ParameterSet& psStripFilter = conf.getParameter<edm::ParameterSet>("HFStripFilter");
144  }
145  }
146 
147  // Describe consumed data
148  tok_PreRecHit_ = consumes<HFPreRecHitCollection>(
149  conf.getParameter<edm::InputTag>("inputLabel"));
150 
151  // Register the product
152  produces<HFRecHitCollection>();
153 }
154 
155 
157 {
158 
159  // do anything here that needs to be done at destruction time
160  // (e.g. close files, deallocate resources etc.)
161 
162 }
163 
164 
165 void
167 {
168  if (reco_->isConfigurable())
169  {
171  if (!recoConfig_.get())
172  throw cms::Exception("HFPhase1BadConfig")
173  << "Invalid HFPhase1Reconstructor \"algoConfigClass\" parameter value \""
174  << algoConfigClass_ << '"' << std::endl;
175  if (!reco_->configure(recoConfig_.get()))
176  throw cms::Exception("HFPhase1BadConfig")
177  << "Failed to configure HFPhase1Reconstructor algorithm from EventSetup"
178  << std::endl;
179  }
180 }
181 
182 // ------------ method called to produce the data ------------
183 void
185 {
186  using namespace edm;
187 
188  // Fetch the calibrations
189  ESHandle<HcalDbService> conditions;
190  eventSetup.get<HcalDbRecord>().get(conditions);
191 
193  eventSetup.get<HcalChannelQualityRcd>().get("withTopo", p);
194  const HcalChannelQuality* myqual = p.product();
195 
197  eventSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
198  const HcalSeverityLevelComputer* mySeverity = mycomputer.product();
199 
200  // Get the input data
202  e.getByToken(tok_PreRecHit_, preRecHits);
203 
204  // Create a new output collection
205  std::unique_ptr<HFRecHitCollection> rec(std::make_unique<HFRecHitCollection>());
206  rec->reserve(preRecHits->size());
207 
208  // Iterate over the input and fill the output collection
209  for (HFPreRecHitCollection::const_iterator it = preRecHits->begin();
210  it != preRecHits->end(); ++it)
211  {
212  // The check whether this PMT is single-anode one should go here.
213  // Fix this piece of code if/when mixed-anode readout configurations
214  // become available.
215  const bool thisIsSingleAnodePMT = false;
216 
217  // Check if the anodes were tagged bad in the database
218  bool taggedBadByDb[2] = {false, false};
220  {
221  if (checkChannelQualityForDepth3and4_ && !thisIsSingleAnodePMT)
222  {
223  HcalDetId anodeIds[2];
224  anodeIds[0] = it->id();
225  anodeIds[1] = anodeIds[0].secondAnodeId();
226  for (unsigned i=0; i<2; ++i)
227  {
228  const HcalChannelStatus* mydigistatus = myqual->getValues(anodeIds[i].rawId());
229  taggedBadByDb[i] = mySeverity->dropChannel(mydigistatus->getValue());
230  }
231  }
232  else
233  {
234  const HcalChannelStatus* mydigistatus = myqual->getValues(it->id().rawId());
235  const bool b = mySeverity->dropChannel(mydigistatus->getValue());
236  taggedBadByDb[0] = b;
237  taggedBadByDb[1] = b;
238  }
239  }
240 
241  // Reconstruct the rechit
242  const HFRecHit& rh = reco_->reconstruct(
243  *it, conditions->getHcalCalibrations(it->id()),
244  taggedBadByDb, thisIsSingleAnodePMT);
245 
246  // The rechit will have the id of 0 if the algorithm
247  // decides that it should be dropped
248  if (rh.id().rawId())
249  rec->push_back(rh);
250  }
251 
252  // At this point, the rechits contain energy, timing,
253  // as well as proper auxiliary words. We do still need
254  // to set certain flags using the noise cleanup algoritms.
255 
256  // The following flags require the full set of rechits.
257  // These flags need to be set consecutively.
258  if (setNoiseFlags_)
259  {
260  // Step 1: Set PET flag (short fibers of |ieta|==29)
261  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i)
262  {
263  int depth=i->id().depth();
264  int ieta=i->id().ieta();
265  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
266  if (depth==2 || abs(ieta)==29 )
267  hfPET_->HFSetFlagFromPET(*i, *rec, myqual, mySeverity);
268  }
269 
270  // Step 2: Set S8S1 flag (short fibers or |ieta|==29)
271  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i)
272  {
273  int depth=i->id().depth();
274  int ieta=i->id().ieta();
275  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
276  if (depth==2 || abs(ieta)==29 )
277  hfS8S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
278  }
279 
280  // Step 3: Set S9S1 flag (long fibers)
281  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i)
282  {
283  int depth=i->id().depth();
284  int ieta=i->id().ieta();
285  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
286  if (depth==1 && abs(ieta)!=29 )
287  hfS9S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
288  }
289 
290  // Step 4: Run HFStripFilter if requested
291  if (runHFStripFilter_)
292  hfStripFilter_->runFilter(*rec, myqual);
293  }
294 
295  // Add the output collection to the event record
296  e.put(std::move(rec));
297 }
298 
299 #define add_param_set(name) \
300  edm::ParameterSetDescription name; \
301  name.setAllowAnything(); \
302  desc.add<edm::ParameterSetDescription>(#name, name)
303 
304 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
305 void
307 {
309 
310  desc.add<edm::InputTag>("inputLabel");
311  desc.add<std::string>("algoConfigClass");
312  desc.add<bool>("setNoiseFlags");
313  desc.add<bool>("runHFStripFilter", false);
314  desc.add<bool>("useChannelQualityFromDB");
315  desc.add<bool>("checkChannelQualityForDepth3and4");
318 
322 
323  descriptions.addDefault(desc);
324 }
325 
326 //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
void beginRun(const edm::Run &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
static edm::ParameterSetDescription fillDescription()
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
std::vector< T >::const_iterator const_iterator
const Item * getValues(DetId fId, bool throwOnFail=true) const
std::unique_ptr< AbsHcalAlgoData > fetchHcalAlgoData(const std::string &className, const edm::EventSetup &es)
std::unique_ptr< HcalHF_PETalgorithm > hfPET_
#define add_param_set(name)
void produce(edm::Event &, const edm::EventSetup &) override
std::unique_ptr< HFStripFilter > hfStripFilter_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
std::unique_ptr< AbsHcalAlgoData > recoConfig_
std::unique_ptr< AbsHFPhase1Algo > reco_
std::unique_ptr< AbsHFPhase1Algo > parseHFPhase1AlgoDescription(const edm::ParameterSet &ps)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool dropChannel(const uint32_t &mystatus) const
static std::unique_ptr< HFStripFilter > parseParameterSet(const edm::ParameterSet &ps)
std::unique_ptr< HcalHF_S9S1algorithm > hfS8S1_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< T >::iterator iterator
const_iterator end() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double b
Definition: hdecay.h:120
HFPhase1Reconstructor(const edm::ParameterSet &)
HLT enums.
size_type size() const
T get() const
Definition: EventSetup.h:71
edm::ParameterSetDescription fillDescriptionForParseHFPhase1AlgoDescription()
HcalDetId id() const
Definition: HFRecHit.h:31
uint32_t getValue() const
const HcalCalibrations & getHcalCalibrations(const HcalGenericDetId &fId) const
edm::EDGetTokenT< HFPreRecHitCollection > tok_PreRecHit_
std::unique_ptr< HcalHF_S9S1algorithm > hfS9S1_
T const * product() const
Definition: ESHandle.h:86
HcalDetId secondAnodeId() const
second PMT anode detId for HF dual channels
Definition: HcalDetId.h:240
def move(src, dest)
Definition: eostools.py:511
const_iterator begin() const
Definition: Run.h:45