Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:52

0001 #include "HcalSimpleReconstructor.h"
0002 #include "DataFormats/Common/interface/EDCollection.h"
0003 #include "DataFormats/Common/interface/Handle.h"
0004 #include "FWCore/Framework/interface/ConsumesCollector.h"
0005 #include "FWCore/Framework/interface/EventSetup.h"
0006 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0008 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
0009 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
0010 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
0011 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
0012 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0013 
0014 #include <iostream>
0015 
0016 HcalSimpleReconstructor::HcalSimpleReconstructor(edm::ParameterSet const& conf)
0017     : reco_(conf.getParameter<bool>("correctForTimeslew"),
0018             conf.getParameter<bool>("correctForPhaseContainment"),
0019             conf.getParameter<double>("correctionPhaseNS"),
0020             consumesCollector()),
0021       det_(DetId::Hcal),
0022       inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
0023       dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
0024       firstSample_(conf.getParameter<int>("firstSample")),
0025       samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
0026       tsFromDB_(conf.getParameter<bool>("tsFromDB")) {
0027   // register for data access
0028   tok_hf_ = consumes<HFDigiCollection>(inputLabel_);
0029   tok_ho_ = consumes<HODigiCollection>(inputLabel_);
0030   tok_calib_ = consumes<HcalCalibDigiCollection>(inputLabel_);
0031 
0032   std::string subd = conf.getParameter<std::string>("Subdetector");
0033   if (!strcasecmp(subd.c_str(), "HO")) {
0034     subdet_ = HcalOuter;
0035     produces<HORecHitCollection>();
0036   } else if (!strcasecmp(subd.c_str(), "HF")) {
0037     subdet_ = HcalForward;
0038     produces<HFRecHitCollection>();
0039   } else {
0040     std::cout << "HcalSimpleReconstructor is not associated with a specific subdetector!" << std::endl;
0041   }
0042 
0043   // ES tokens
0044   if (tsFromDB_) {
0045     htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
0046     paramsToken_ = esConsumes<HcalRecoParams, HcalRecoParamsRcd, edm::Transition::BeginRun>();
0047   }
0048   conditionsToken_ = esConsumes<HcalDbService, HcalDbRecord>();
0049 }
0050 
0051 HcalSimpleReconstructor::~HcalSimpleReconstructor() {}
0052 
0053 void HcalSimpleReconstructor::beginRun(edm::Run const& r, edm::EventSetup const& es) {
0054   if (tsFromDB_) {
0055     const HcalTopology& htopo = es.getData(htopoToken_);
0056     const HcalRecoParams& p = es.getData(paramsToken_);
0057     paramTS_ = std::make_unique<HcalRecoParams>(p);
0058     paramTS_->setTopo(&htopo);
0059   }
0060   reco_.beginRun(es);
0061 }
0062 
0063 void HcalSimpleReconstructor::endRun(edm::Run const& r, edm::EventSetup const& es) { reco_.endRun(); }
0064 
0065 template <class DIGICOLL, class RECHITCOLL>
0066 void HcalSimpleReconstructor::process(edm::Event& e,
0067                                       const edm::EventSetup& eventSetup,
0068                                       const edm::EDGetTokenT<DIGICOLL>& tok) {
0069   // get conditions
0070   const HcalDbService* conditions = &eventSetup.getData(conditionsToken_);
0071 
0072   edm::Handle<DIGICOLL> digi;
0073   e.getByToken(tok, digi);
0074 
0075   // create empty output
0076   auto rec = std::make_unique<RECHITCOLL>();
0077   rec->reserve(digi->size());
0078   // run the algorithm
0079   int first = firstSample_;
0080   int toadd = samplesToAdd_;
0081   typename DIGICOLL::const_iterator i;
0082   for (i = digi->begin(); i != digi->end(); i++) {
0083     HcalDetId cell = i->id();
0084     DetId detcell = (DetId)cell;
0085     // rof 27.03.09: drop ZS marked and passed digis:
0086     if (dropZSmarkedPassed_)
0087       if (i->zsMarkAndPass())
0088         continue;
0089 
0090     const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0091     const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0092     const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0093     HcalCoderDb coder(*channelCoder, *shape);
0094 
0095     //>>> firstSample & samplesToAdd
0096     if (tsFromDB_) {
0097       const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0098       first = param_ts->firstSample();
0099       toadd = param_ts->samplesToAdd();
0100     }
0101     rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0102   }
0103   // return result
0104   e.put(std::move(rec));
0105 }
0106 
0107 void HcalSimpleReconstructor::produce(edm::Event& e, const edm::EventSetup& eventSetup) {
0108   if (det_ == DetId::Hcal) {
0109     if (subdet_ == HcalForward) {
0110       process<HFDigiCollection, HFRecHitCollection>(e, eventSetup, tok_hf_);
0111     } else if (subdet_ == HcalOuter) {
0112       process<HODigiCollection, HORecHitCollection>(e, eventSetup, tok_ho_);
0113     } else if (subdet_ == HcalOther && subdetOther_ == HcalCalibration) {
0114       process<HcalCalibDigiCollection, HcalCalibRecHitCollection>(e, eventSetup, tok_calib_);
0115     }
0116   }
0117 }
0118 
0119 void HcalSimpleReconstructor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0120   // horeco
0121   edm::ParameterSetDescription descHO;
0122   descHO.add<double>("correctionPhaseNS", 13.0);
0123   descHO.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
0124   descHO.add<bool>("tsFromDB", true);
0125   descHO.add<int>("samplesToAdd", 4);
0126   descHO.add<std::string>("Subdetector", "HO");
0127   descHO.add<bool>("correctForTimeslew", true);
0128   descHO.add<bool>("dropZSmarkedPassed", true);
0129   descHO.add<bool>("correctForPhaseContainment", true);
0130   descHO.add<int>("firstSample", 4);
0131   descriptions.add("hosimplereco", descHO);
0132 
0133   // hfreco
0134   edm::ParameterSetDescription descHF;
0135   descHF.add<double>("correctionPhaseNS", 0.0);
0136   descHF.add<edm::InputTag>("digiLabel", edm::InputTag("hcalDigis"));
0137   descHF.add<bool>("tsFromDB", true);
0138   descHF.add<int>("samplesToAdd", 2);
0139   descHF.add<std::string>("Subdetector", "HF");
0140   descHF.add<bool>("correctForTimeslew", false);
0141   descHF.add<bool>("dropZSmarkedPassed", true);
0142   descHF.add<bool>("correctForPhaseContainment", false);
0143   descHF.add<int>("firstSample", 4);
0144   descriptions.add("hfsimplereco", descHF);
0145 }