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 public:
57  ~HFPhase1Reconstructor() override;
58 
59  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
60 
61 private:
62  void beginRun(const edm::Run&, const edm::EventSetup&) override;
63  void produce(edm::Event&, const edm::EventSetup&) override;
64 
65  // Configuration parameters
71 
72  // Other members
74  std::unique_ptr<AbsHFPhase1Algo> reco_;
75  std::unique_ptr<AbsHcalAlgoData> recoConfig_;
76 
77  // Noise cleanup algos
78  std::unique_ptr<HcalHF_S9S1algorithm> hfS9S1_;
79  std::unique_ptr<HcalHF_S9S1algorithm> hfS8S1_;
80  std::unique_ptr<HcalHF_PETalgorithm> hfPET_;
81  std::unique_ptr<HFStripFilter> hfStripFilter_;
82 };
83 
84 //
85 // constructors and destructor
86 //
88  : algoConfigClass_(conf.getParameter<std::string>("algoConfigClass")),
89  setNoiseFlags_(conf.getParameter<bool>("setNoiseFlags")),
90  runHFStripFilter_(conf.getParameter<bool>("runHFStripFilter")),
91  useChannelQualityFromDB_(conf.getParameter<bool>("useChannelQualityFromDB")),
92  checkChannelQualityForDepth3and4_(conf.getParameter<bool>("checkChannelQualityForDepth3and4")),
93  reco_(parseHFPhase1AlgoDescription(conf.getParameter<edm::ParameterSet>("algorithm"))) {
94  // Check that the reco algorithm has been successfully configured
95  if (!reco_.get())
96  throw cms::Exception("HFPhase1BadConfig") << "Invalid HFPhase1Reconstructor algorithm configuration" << std::endl;
97 
98  // Configure the noise cleanup algorithms
99  if (setNoiseFlags_) {
100  const edm::ParameterSet& psS9S1 = conf.getParameter<edm::ParameterSet>("S9S1stat");
101  hfS9S1_ = std::make_unique<HcalHF_S9S1algorithm>(psS9S1.getParameter<std::vector<double> >("short_optimumSlope"),
102  psS9S1.getParameter<std::vector<double> >("shortEnergyParams"),
103  psS9S1.getParameter<std::vector<double> >("shortETParams"),
104  psS9S1.getParameter<std::vector<double> >("long_optimumSlope"),
105  psS9S1.getParameter<std::vector<double> >("longEnergyParams"),
106  psS9S1.getParameter<std::vector<double> >("longETParams"),
107  psS9S1.getParameter<int>("HcalAcceptSeverityLevel"),
108  psS9S1.getParameter<bool>("isS8S1"));
109 
110  const edm::ParameterSet& psS8S1 = conf.getParameter<edm::ParameterSet>("S8S1stat");
111  hfS8S1_ = std::make_unique<HcalHF_S9S1algorithm>(psS8S1.getParameter<std::vector<double> >("short_optimumSlope"),
112  psS8S1.getParameter<std::vector<double> >("shortEnergyParams"),
113  psS8S1.getParameter<std::vector<double> >("shortETParams"),
114  psS8S1.getParameter<std::vector<double> >("long_optimumSlope"),
115  psS8S1.getParameter<std::vector<double> >("longEnergyParams"),
116  psS8S1.getParameter<std::vector<double> >("longETParams"),
117  psS8S1.getParameter<int>("HcalAcceptSeverityLevel"),
118  psS8S1.getParameter<bool>("isS8S1"));
119 
120  const edm::ParameterSet& psPET = conf.getParameter<edm::ParameterSet>("PETstat");
121  hfPET_ = std::make_unique<HcalHF_PETalgorithm>(psPET.getParameter<std::vector<double> >("short_R"),
122  psPET.getParameter<std::vector<double> >("shortEnergyParams"),
123  psPET.getParameter<std::vector<double> >("shortETParams"),
124  psPET.getParameter<std::vector<double> >("long_R"),
125  psPET.getParameter<std::vector<double> >("longEnergyParams"),
126  psPET.getParameter<std::vector<double> >("longETParams"),
127  psPET.getParameter<int>("HcalAcceptSeverityLevel"),
128  psPET.getParameter<std::vector<double> >("short_R_29"),
129  psPET.getParameter<std::vector<double> >("long_R_29"));
130 
131  // Configure HFStripFilter
132  if (runHFStripFilter_) {
133  const edm::ParameterSet& psStripFilter = conf.getParameter<edm::ParameterSet>("HFStripFilter");
135  }
136  }
137 
138  // Describe consumed data
139  tok_PreRecHit_ = consumes<HFPreRecHitCollection>(conf.getParameter<edm::InputTag>("inputLabel"));
140 
141  // Register the product
142  produces<HFRecHitCollection>();
143 }
144 
146  // do anything here that needs to be done at destruction time
147  // (e.g. close files, deallocate resources etc.)
148 }
149 
151  if (reco_->isConfigurable()) {
153  if (!recoConfig_.get())
154  throw cms::Exception("HFPhase1BadConfig")
155  << "Invalid HFPhase1Reconstructor \"algoConfigClass\" parameter value \"" << algoConfigClass_ << '"'
156  << std::endl;
157  if (!reco_->configure(recoConfig_.get()))
158  throw cms::Exception("HFPhase1BadConfig")
159  << "Failed to configure HFPhase1Reconstructor algorithm from EventSetup" << std::endl;
160  }
161 }
162 
163 // ------------ method called to produce the data ------------
165  using namespace edm;
166 
167  // Fetch the calibrations
168  ESHandle<HcalDbService> conditions;
169  eventSetup.get<HcalDbRecord>().get(conditions);
170 
172  eventSetup.get<HcalChannelQualityRcd>().get("withTopo", p);
173  const HcalChannelQuality* myqual = p.product();
174 
176  eventSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
177  const HcalSeverityLevelComputer* mySeverity = mycomputer.product();
178 
179  // Get the input data
181  e.getByToken(tok_PreRecHit_, preRecHits);
182 
183  // Create a new output collection
184  std::unique_ptr<HFRecHitCollection> rec(std::make_unique<HFRecHitCollection>());
185  rec->reserve(preRecHits->size());
186 
187  // Iterate over the input and fill the output collection
188  for (HFPreRecHitCollection::const_iterator it = preRecHits->begin(); it != preRecHits->end(); ++it) {
189  // The check whether this PMT is single-anode one should go here.
190  // Fix this piece of code if/when mixed-anode readout configurations
191  // become available.
192  const bool thisIsSingleAnodePMT = false;
193 
194  // Check if the anodes were tagged bad in the database
195  bool taggedBadByDb[2] = {false, false};
197  if (checkChannelQualityForDepth3and4_ && !thisIsSingleAnodePMT) {
198  HcalDetId anodeIds[2];
199  anodeIds[0] = it->id();
200  anodeIds[1] = anodeIds[0].secondAnodeId();
201  for (unsigned i = 0; i < 2; ++i) {
202  const HcalChannelStatus* mydigistatus = myqual->getValues(anodeIds[i].rawId());
203  taggedBadByDb[i] = mySeverity->dropChannel(mydigistatus->getValue());
204  }
205  } else {
206  const HcalChannelStatus* mydigistatus = myqual->getValues(it->id().rawId());
207  const bool b = mySeverity->dropChannel(mydigistatus->getValue());
208  taggedBadByDb[0] = b;
209  taggedBadByDb[1] = b;
210  }
211  }
212 
213  // Reconstruct the rechit
214  const HFRecHit& rh =
215  reco_->reconstruct(*it, conditions->getHcalCalibrations(it->id()), taggedBadByDb, thisIsSingleAnodePMT);
216 
217  // The rechit will have the id of 0 if the algorithm
218  // decides that it should be dropped
219  if (rh.id().rawId())
220  rec->push_back(rh);
221  }
222 
223  // At this point, the rechits contain energy, timing,
224  // as well as proper auxiliary words. We do still need
225  // to set certain flags using the noise cleanup algoritms.
226 
227  // The following flags require the full set of rechits.
228  // These flags need to be set consecutively.
229  if (setNoiseFlags_) {
230  // Step 1: Set PET flag (short fibers of |ieta|==29)
231  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
232  int depth = i->id().depth();
233  int ieta = i->id().ieta();
234  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
235  if (depth == 2 || abs(ieta) == 29)
236  hfPET_->HFSetFlagFromPET(*i, *rec, myqual, mySeverity);
237  }
238 
239  // Step 2: Set S8S1 flag (short fibers or |ieta|==29)
240  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
241  int depth = i->id().depth();
242  int ieta = i->id().ieta();
243  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
244  if (depth == 2 || abs(ieta) == 29)
245  hfS8S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
246  }
247 
248  // Step 3: Set S9S1 flag (long fibers)
249  for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
250  int depth = i->id().depth();
251  int ieta = i->id().ieta();
252  // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
253  if (depth == 1 && abs(ieta) != 29)
254  hfS9S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
255  }
256 
257  // Step 4: Run HFStripFilter if requested
258  if (runHFStripFilter_)
259  hfStripFilter_->runFilter(*rec, myqual);
260  }
261 
262  // Add the output collection to the event record
263  e.put(std::move(rec));
264 }
265 
266 #define add_param_set(name) \
267  edm::ParameterSetDescription name; \
268  name.setAllowAnything(); \
269  desc.add<edm::ParameterSetDescription>(#name, name)
270 
271 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
274 
275  desc.add<edm::InputTag>("inputLabel");
276  desc.add<std::string>("algoConfigClass");
277  desc.add<bool>("setNoiseFlags");
278  desc.add<bool>("runHFStripFilter", false);
279  desc.add<bool>("useChannelQualityFromDB");
280  desc.add<bool>("checkChannelQualityForDepth3and4");
283 
287 
288  descriptions.addDefault(desc);
289 }
290 
291 //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:131
void beginRun(const edm::Run &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
static edm::ParameterSetDescription fillDescription()
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
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:118
HFPhase1Reconstructor(const edm::ParameterSet &)
HLT enums.
size_type size() const
T get() const
Definition: EventSetup.h:73
edm::ParameterSetDescription fillDescriptionForParseHFPhase1AlgoDescription()
HcalDetId id() const
Definition: HFRecHit.h:26
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:241
def move(src, dest)
Definition: eostools.py:511
const_iterator begin() const
Definition: Run.h:45