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