Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "HcalHitReconstructor.h"
0002 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0003 #include "DataFormats/Common/interface/EDCollection.h"
0004 #include "DataFormats/Common/interface/Handle.h"
0005 #include "FWCore/Framework/interface/ConsumesCollector.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
0008 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
0009 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
0010 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
0011 #include "CalibCalorimetry/HcalAlgos/interface/HcalDbASCIIIO.h"
0012 #include "Geometry/CaloTopology/interface/HcalTopology.h"
0013 #include "Geometry/Records/interface/HcalRecNumberingRecord.h"
0014 #include "CondFormats/DataRecord/interface/HcalFrontEndMapRcd.h"
0015 #include "CondFormats/DataRecord/interface/HcalOOTPileupCorrectionRcd.h"
0016 #include "CondFormats/DataRecord/interface/HcalOOTPileupCompatibilityRcd.h"
0017 #include "CondFormats/DataRecord/interface/HBHENegativeEFilterRcd.h"
0018 #include "CondFormats/HcalObjects/interface/HcalFrontEndMap.h"
0019 #include "CondFormats/HcalObjects/interface/OOTPileupCorrectionColl.h"
0020 #include "CondFormats/HcalObjects/interface/OOTPileupCorrData.h"
0021 #include <iostream>
0022 #include <fstream>
0023 
0024 /*  Hcal Hit reconstructor allows for CaloRecHits with status words */
0025 
0026 HcalHitReconstructor::HcalHitReconstructor(edm::ParameterSet const& conf)
0027     : reco_(conf.getParameter<bool>("correctForTimeslew"),
0028             conf.getParameter<bool>("correctForPhaseContainment"),
0029             conf.getParameter<double>("correctionPhaseNS"),
0030             consumesCollector()),
0031       det_(DetId::Hcal),
0032       inputLabel_(conf.getParameter<edm::InputTag>("digiLabel")),
0033       correctTiming_(conf.getParameter<bool>("correctTiming")),
0034       setNoiseFlags_(conf.getParameter<bool>("setNoiseFlags")),
0035       setHSCPFlags_(conf.getParameter<bool>("setHSCPFlags")),
0036       setSaturationFlags_(conf.getParameter<bool>("setSaturationFlags")),
0037       setTimingTrustFlags_(conf.getParameter<bool>("setTimingTrustFlags")),
0038       setPulseShapeFlags_(conf.getParameter<bool>("setPulseShapeFlags")),
0039       setNegativeFlags_(false),
0040       dropZSmarkedPassed_(conf.getParameter<bool>("dropZSmarkedPassed")),
0041       firstAuxTS_(conf.getParameter<int>("firstAuxTS")),
0042       firstSample_(conf.getParameter<int>("firstSample")),
0043       samplesToAdd_(conf.getParameter<int>("samplesToAdd")),
0044       tsFromDB_(conf.getParameter<bool>("tsFromDB")),
0045       useLeakCorrection_(conf.getParameter<bool>("useLeakCorrection")),
0046       dataOOTCorrectionName_(""),
0047       dataOOTCorrectionCategory_("Data"),
0048       mcOOTCorrectionName_(""),
0049       mcOOTCorrectionCategory_("MC"),
0050       setPileupCorrection_(nullptr) {
0051   // register for data access
0052   tok_ho_ = consumes<HODigiCollection>(inputLabel_);
0053   tok_hf_ = consumes<HFDigiCollection>(inputLabel_);
0054   tok_calib_ = consumes<HcalCalibDigiCollection>(inputLabel_);
0055 
0056   std::string subd = conf.getParameter<std::string>("Subdetector");
0057   //Set all FlagSetters to 0
0058   /* Important to do this!  Otherwise, if the setters are turned off,
0059      the "if (XSetter_) delete XSetter_;" commands can crash
0060   */
0061 
0062   recoParamsFromDB_ = conf.getParameter<bool>("recoParamsFromDB");
0063   //  recoParamsFromDB_ = false ; //  trun off for now.
0064 
0065   // std::cout<<"  HcalHitReconstructor   recoParamsFromDB_ "<<recoParamsFromDB_<<std::endl;
0066 
0067   if (conf.existsAs<bool>("setNegativeFlags"))
0068     setNegativeFlags_ = conf.getParameter<bool>("setNegativeFlags");
0069 
0070   hfdigibit_ = nullptr;
0071 
0072   hfS9S1_ = nullptr;
0073   hfS8S1_ = nullptr;
0074   hfPET_ = nullptr;
0075   saturationFlagSetter_ = nullptr;
0076   HFTimingTrustFlagSetter_ = nullptr;
0077   digiTimeFromDB_ = false;  // only need for HF
0078 
0079   if (setSaturationFlags_) {
0080     const edm::ParameterSet& pssat = conf.getParameter<edm::ParameterSet>("saturationParameters");
0081     saturationFlagSetter_ = new HcalADCSaturationFlag(pssat.getParameter<int>("maxADCvalue"));
0082   }
0083 
0084   if (!strcasecmp(subd.c_str(), "HO")) {
0085     subdet_ = HcalOuter;
0086     // setPileupCorrection_ = &HcalSimpleRecAlgo::setHOPileupCorrection;
0087     setPileupCorrection_ = nullptr;
0088     produces<HORecHitCollection>();
0089   } else if (!strcasecmp(subd.c_str(), "HF")) {
0090     subdet_ = HcalForward;
0091     // setPileupCorrection_ = &HcalSimpleRecAlgo::setHFPileupCorrection;
0092     setPileupCorrection_ = nullptr;
0093     digiTimeFromDB_ = conf.getParameter<bool>("digiTimeFromDB");
0094 
0095     if (setTimingTrustFlags_) {
0096       const edm::ParameterSet& pstrust = conf.getParameter<edm::ParameterSet>("hfTimingTrustParameters");
0097       HFTimingTrustFlagSetter_ = new HFTimingTrustFlag(pstrust.getParameter<int>("hfTimingTrustLevel1"),
0098                                                        pstrust.getParameter<int>("hfTimingTrustLevel2"));
0099     }
0100 
0101     if (setNoiseFlags_) {
0102       const edm::ParameterSet& psdigi = conf.getParameter<edm::ParameterSet>("digistat");
0103       const edm::ParameterSet& psTimeWin = conf.getParameter<edm::ParameterSet>("HFInWindowStat");
0104       hfdigibit_ = new HcalHFStatusBitFromDigis(psdigi, psTimeWin);
0105 
0106       const edm::ParameterSet& psS9S1 = conf.getParameter<edm::ParameterSet>("S9S1stat");
0107       hfS9S1_ = new HcalHF_S9S1algorithm(psS9S1.getParameter<std::vector<double> >("short_optimumSlope"),
0108                                          psS9S1.getParameter<std::vector<double> >("shortEnergyParams"),
0109                                          psS9S1.getParameter<std::vector<double> >("shortETParams"),
0110                                          psS9S1.getParameter<std::vector<double> >("long_optimumSlope"),
0111                                          psS9S1.getParameter<std::vector<double> >("longEnergyParams"),
0112                                          psS9S1.getParameter<std::vector<double> >("longETParams"),
0113                                          psS9S1.getParameter<int>("HcalAcceptSeverityLevel"),
0114                                          psS9S1.getParameter<bool>("isS8S1"));
0115 
0116       const edm::ParameterSet& psS8S1 = conf.getParameter<edm::ParameterSet>("S8S1stat");
0117       hfS8S1_ = new HcalHF_S9S1algorithm(psS8S1.getParameter<std::vector<double> >("short_optimumSlope"),
0118                                          psS8S1.getParameter<std::vector<double> >("shortEnergyParams"),
0119                                          psS8S1.getParameter<std::vector<double> >("shortETParams"),
0120                                          psS8S1.getParameter<std::vector<double> >("long_optimumSlope"),
0121                                          psS8S1.getParameter<std::vector<double> >("longEnergyParams"),
0122                                          psS8S1.getParameter<std::vector<double> >("longETParams"),
0123                                          psS8S1.getParameter<int>("HcalAcceptSeverityLevel"),
0124                                          psS8S1.getParameter<bool>("isS8S1"));
0125 
0126       const edm::ParameterSet& psPET = conf.getParameter<edm::ParameterSet>("PETstat");
0127       hfPET_ = new HcalHF_PETalgorithm(psPET.getParameter<std::vector<double> >("short_R"),
0128                                        psPET.getParameter<std::vector<double> >("shortEnergyParams"),
0129                                        psPET.getParameter<std::vector<double> >("shortETParams"),
0130                                        psPET.getParameter<std::vector<double> >("long_R"),
0131                                        psPET.getParameter<std::vector<double> >("longEnergyParams"),
0132                                        psPET.getParameter<std::vector<double> >("longETParams"),
0133                                        psPET.getParameter<int>("HcalAcceptSeverityLevel"),
0134                                        psPET.getParameter<std::vector<double> >("short_R_29"),
0135                                        psPET.getParameter<std::vector<double> >("long_R_29"));
0136     }
0137     produces<HFRecHitCollection>();
0138   } else if (!strcasecmp(subd.c_str(), "ZDC")) {
0139     det_ = DetId::Calo;
0140     subdet_ = HcalZDCDetId::SubdetectorId;
0141     produces<ZDCRecHitCollection>();
0142   } else if (!strcasecmp(subd.c_str(), "CALIB")) {
0143     subdet_ = HcalOther;
0144     subdetOther_ = HcalCalibration;
0145     produces<HcalCalibRecHitCollection>();
0146   } else {
0147     edm::LogWarning("Configuration") << "HcalHitReconstructor is not associated with a specific subdetector!"
0148                                      << std::endl;
0149   }
0150 
0151   // If no valid OOT pileup correction name specified,
0152   // disable the correction
0153   if (conf.existsAs<std::string>("dataOOTCorrectionName"))
0154     dataOOTCorrectionName_ = conf.getParameter<std::string>("dataOOTCorrectionName");
0155   if (conf.existsAs<std::string>("dataOOTCorrectionCategory"))
0156     dataOOTCorrectionCategory_ = conf.getParameter<std::string>("dataOOTCorrectionCategory");
0157   if (conf.existsAs<std::string>("mcOOTCorrectionName"))
0158     mcOOTCorrectionName_ = conf.getParameter<std::string>("mcOOTCorrectionName");
0159   if (conf.existsAs<std::string>("mcOOTCorrectionCategory"))
0160     mcOOTCorrectionCategory_ = conf.getParameter<std::string>("mcOOTCorrectionCategory");
0161   if (dataOOTCorrectionName_.empty() && mcOOTCorrectionName_.empty())
0162     setPileupCorrection_ = nullptr;
0163 
0164   // ES tokens
0165   htopoToken_ = esConsumes<HcalTopology, HcalRecNumberingRecord, edm::Transition::BeginRun>();
0166   if (tsFromDB_ || recoParamsFromDB_)
0167     paramsToken_ = esConsumes<HcalRecoParams, HcalRecoParamsRcd, edm::Transition::BeginRun>();
0168   if (digiTimeFromDB_)
0169     digiTimeToken_ = esConsumes<HcalFlagHFDigiTimeParams, HcalFlagHFDigiTimeParamsRcd, edm::Transition::BeginRun>();
0170   conditionsToken_ = esConsumes<HcalDbService, HcalDbRecord>();
0171   qualToken_ = esConsumes<HcalChannelQuality, HcalChannelQualityRcd>(edm::ESInputTag("", "withTopo"));
0172   sevToken_ = esConsumes<HcalSeverityLevelComputer, HcalSeverityLevelComputerRcd>();
0173 }
0174 
0175 HcalHitReconstructor::~HcalHitReconstructor() {
0176   delete hfdigibit_;
0177 
0178   delete hfS9S1_;
0179   delete hfS8S1_;
0180   delete hfPET_;
0181   delete saturationFlagSetter_;
0182   delete HFTimingTrustFlagSetter_;
0183 }
0184 
0185 void HcalHitReconstructor::beginRun(edm::Run const& r, edm::EventSetup const& es) {
0186   const HcalTopology& htopo = es.getData(htopoToken_);
0187 
0188   if (tsFromDB_ || recoParamsFromDB_) {
0189     const HcalRecoParams& p = es.getData(paramsToken_);
0190     paramTS_ = std::make_unique<HcalRecoParams>(p);
0191     paramTS_->setTopo(&htopo);
0192 
0193     // std::cout<<" skdump in HcalHitReconstructor::beginRun   dupm RecoParams "<<std::endl;
0194     // std::ofstream skfile("skdumpRecoParamsNewFormat.txt");
0195     // HcalDbASCIIIO::dumpObject(skfile, (*paramTS_) );
0196   }
0197 
0198   if (digiTimeFromDB_) {
0199     const HcalFlagHFDigiTimeParams& p = es.getData(digiTimeToken_);
0200     hFDigiTimeParams_ = std::make_unique<HcalFlagHFDigiTimeParams>(p);
0201     hFDigiTimeParams_->setTopo(&htopo);
0202   }
0203 
0204   reco_.beginRun(es);
0205 }
0206 
0207 void HcalHitReconstructor::endRun(edm::Run const& r, edm::EventSetup const& es) { reco_.endRun(); }
0208 
0209 void HcalHitReconstructor::produce(edm::Event& e, const edm::EventSetup& eventSetup) {
0210   // get conditions
0211   const HcalDbService* conditions = &eventSetup.getData(conditionsToken_);
0212   const HcalChannelQuality* myqual = &eventSetup.getData(qualToken_);
0213   const HcalSeverityLevelComputer* mySeverity = &eventSetup.getData(sevToken_);
0214 
0215   if (useLeakCorrection_)
0216     reco_.setLeakCorrection();
0217 
0218   // GET THE BEAM CROSSING INFO HERE, WHEN WE UNDERSTAND HOW THINGS WORK.
0219   // Then, call "setBXInfo" method of the reco_ object.
0220   // Also remember to call SetBXInfo in the negative energy flag setter.
0221 
0222   if (det_ == DetId::Hcal) {
0223     //  HO ------------------------------------------------------------------
0224     if (subdet_ == HcalOuter) {
0225       edm::Handle<HODigiCollection> digi;
0226       e.getByToken(tok_ho_, digi);
0227 
0228       // create empty output
0229       auto rec = std::make_unique<HORecHitCollection>();
0230       rec->reserve(digi->size());
0231       // run the algorithm
0232       HODigiCollection::const_iterator i;
0233 
0234       // Vote on majority TS0 CapId
0235       int favorite_capid = 0;
0236       if (correctTiming_) {
0237         long capid_votes[4] = {0, 0, 0, 0};
0238         for (i = digi->begin(); i != digi->end(); i++) {
0239           capid_votes[(*i)[0].capid()]++;
0240         }
0241         for (int k = 0; k < 4; k++)
0242           if (capid_votes[k] > capid_votes[favorite_capid])
0243             favorite_capid = k;
0244       }
0245 
0246       for (i = digi->begin(); i != digi->end(); i++) {
0247         HcalDetId cell = i->id();
0248         DetId detcell = (DetId)cell;
0249         // firstSample & samplesToAdd
0250         if (tsFromDB_ || recoParamsFromDB_) {
0251           const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0252           if (tsFromDB_) {
0253             firstSample_ = param_ts->firstSample();
0254             samplesToAdd_ = param_ts->samplesToAdd();
0255           }
0256           if (recoParamsFromDB_) {
0257             bool correctForTimeslew = param_ts->correctForTimeslew();
0258             bool correctForPhaseContainment = param_ts->correctForPhaseContainment();
0259             float phaseNS = param_ts->correctionPhaseNS();
0260             useLeakCorrection_ = param_ts->useLeakCorrection();
0261             correctTiming_ = param_ts->correctTiming();
0262             firstAuxTS_ = param_ts->firstAuxTS();
0263             int pileupCleaningID = param_ts->pileupCleaningID();
0264             reco_.setRecoParams(
0265                 correctForTimeslew, correctForPhaseContainment, useLeakCorrection_, pileupCleaningID, phaseNS);
0266           }
0267         }
0268 
0269         int first = firstSample_;
0270         int toadd = samplesToAdd_;
0271 
0272         if (first >= i->size() || first < 0)
0273           edm::LogWarning("Configuration")
0274               << "HcalHitReconstructor: illegal firstSample" << first << "  in subdet " << subdet_ << std::endl;
0275 
0276         // check on cells to be ignored and dropped: (rof,20.Feb.09)
0277         const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
0278         if (mySeverity->dropChannel(mydigistatus->getValue()))
0279           continue;
0280         if (dropZSmarkedPassed_)
0281           if (i->zsMarkAndPass())
0282             continue;
0283 
0284         const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0285         const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0286         const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0287         HcalCoderDb coder(*channelCoder, *shape);
0288 
0289         rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0290 
0291         // Set auxiliary flag
0292         int auxflag = 0;
0293         int fTS = firstAuxTS_;
0294         if (fTS < 0)
0295           fTS = 0;  //silly protection against negative time slice values
0296         for (int xx = fTS; xx < fTS + 4 && xx < i->size(); ++xx)
0297           auxflag += (i->sample(xx).adc())
0298                      << (7 *
0299                          (xx - fTS));  // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
0300         // bits 28 and 29 are reserved for capid of the first time slice saved in aux
0301         auxflag += ((i->sample(fTS).capid()) << 28);
0302         (rec->back()).setAux(auxflag);
0303         // (rec->back()).setFlags(0);  Don't want to do this because the algorithm
0304         //                             can already set some flags
0305         // Fill Presample ADC flag
0306         if (fTS > 0)
0307           (rec->back()).setFlagField((i->sample(fTS - 1).adc()), HcalCaloFlagLabels::PresampleADC, 7);
0308 
0309         if (setSaturationFlags_)
0310           saturationFlagSetter_->setSaturationFlag(rec->back(), *i);
0311         if (correctTiming_)
0312           HcalTimingCorrector::Correct(rec->back(), *i, favorite_capid);
0313       }
0314       // return result
0315       e.put(std::move(rec));
0316 
0317       // HF -------------------------------------------------------------------
0318     } else if (subdet_ == HcalForward) {
0319       edm::Handle<HFDigiCollection> digi;
0320       e.getByToken(tok_hf_, digi);
0321 
0322       ///////////////////////////////////////////////////////////////// HF
0323       // create empty output
0324       auto rec = std::make_unique<HFRecHitCollection>();
0325       rec->reserve(digi->size());
0326       // run the algorithm
0327       HFDigiCollection::const_iterator i;
0328 
0329       // Vote on majority TS0 CapId
0330       int favorite_capid = 0;
0331       if (correctTiming_) {
0332         long capid_votes[4] = {0, 0, 0, 0};
0333         for (i = digi->begin(); i != digi->end(); i++) {
0334           capid_votes[(*i)[0].capid()]++;
0335         }
0336         for (int k = 0; k < 4; k++)
0337           if (capid_votes[k] > capid_votes[favorite_capid])
0338             favorite_capid = k;
0339       }
0340 
0341       for (i = digi->begin(); i != digi->end(); i++) {
0342         HcalDetId cell = i->id();
0343         DetId detcell = (DetId)cell;
0344 
0345         if (tsFromDB_ || recoParamsFromDB_) {
0346           const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0347           if (tsFromDB_) {
0348             firstSample_ = param_ts->firstSample();
0349             samplesToAdd_ = param_ts->samplesToAdd();
0350           }
0351           if (recoParamsFromDB_) {
0352             bool correctForTimeslew = param_ts->correctForTimeslew();
0353             bool correctForPhaseContainment = param_ts->correctForPhaseContainment();
0354             float phaseNS = param_ts->correctionPhaseNS();
0355             useLeakCorrection_ = param_ts->useLeakCorrection();
0356             correctTiming_ = param_ts->correctTiming();
0357             firstAuxTS_ = param_ts->firstAuxTS();
0358             int pileupCleaningID = param_ts->pileupCleaningID();
0359             reco_.setRecoParams(
0360                 correctForTimeslew, correctForPhaseContainment, useLeakCorrection_, pileupCleaningID, phaseNS);
0361           }
0362         }
0363 
0364         int first = firstSample_;
0365         int toadd = samplesToAdd_;
0366 
0367         if (first >= i->size() || first < 0)
0368           edm::LogWarning("Configuration")
0369               << "HcalHitReconstructor: illegal firstSample" << first << "  in subdet " << subdet_ << std::endl;
0370 
0371         // check on cells to be ignored and dropped: (rof,20.Feb.09)
0372         const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
0373         if (mySeverity->dropChannel(mydigistatus->getValue()))
0374           continue;
0375         if (dropZSmarkedPassed_)
0376           if (i->zsMarkAndPass())
0377             continue;
0378 
0379         const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0380         const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0381         const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0382         HcalCoderDb coder(*channelCoder, *shape);
0383 
0384         // Set HFDigiTime flag values from digiTimeFromDB_
0385         if (digiTimeFromDB_ && hfdigibit_) {
0386           const HcalFlagHFDigiTimeParam* hfDTparam = hFDigiTimeParams_->getValues(detcell.rawId());
0387           hfdigibit_->resetParamsFromDB(hfDTparam->HFdigiflagFirstSample(),
0388                                         hfDTparam->HFdigiflagSamplesToAdd(),
0389                                         hfDTparam->HFdigiflagExpectedPeak(),
0390                                         hfDTparam->HFdigiflagMinEThreshold(),
0391                                         hfDTparam->HFdigiflagCoefficients());
0392         }
0393 
0394         //std::cout << "TOADDHF " << toadd << " " << first << " " << std::endl;
0395         rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0396 
0397         // Set auxiliary flag
0398         int auxflag = 0;
0399         int fTS = firstAuxTS_;
0400         if (fTS < 0)
0401           fTS = 0;  // silly protection against negative time slice values
0402         for (int xx = fTS; xx < fTS + 4 && xx < i->size(); ++xx)
0403           auxflag += (i->sample(xx).adc())
0404                      << (7 *
0405                          (xx - fTS));  // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
0406         // bits 28 and 29 are reserved for capid of the first time slice saved in aux
0407         auxflag += ((i->sample(fTS).capid()) << 28);
0408         (rec->back()).setAux(auxflag);
0409 
0410         // (rec->back()).setFlags(0);  Don't want to do this because the algorithm
0411         //                             can already set some flags
0412 
0413         // Fill Presample ADC flag
0414         if (fTS > 0)
0415           (rec->back()).setFlagField((i->sample(fTS - 1).adc()), HcalCaloFlagLabels::PresampleADC, 7);
0416 
0417         // This calls the code for setting the HF noise bit determined from digi shape
0418         if (setNoiseFlags_)
0419           hfdigibit_->hfSetFlagFromDigi(rec->back(), *i, coder, calibrations);
0420         if (setSaturationFlags_)
0421           saturationFlagSetter_->setSaturationFlag(rec->back(), *i);
0422         if (setTimingTrustFlags_)
0423           HFTimingTrustFlagSetter_->setHFTimingTrustFlag(rec->back(), *i);
0424         if (correctTiming_)
0425           HcalTimingCorrector::Correct(rec->back(), *i, favorite_capid);
0426       }  // for (i=digi->begin(); i!=digi->end(); i++) -- loop on all HF digis
0427 
0428       // The following flags require the full set of rechits
0429       // These need to be set consecutively, so an energy check should be the first
0430       // test performed on these hits (to minimize the loop time)
0431       if (setNoiseFlags_) {
0432         // Step 1:  Set PET flag  (short fibers of |ieta|==29)
0433         // Neighbor/partner channels that are flagged by Pulse Shape algorithm (HFDigiTime)
0434         // won't be considered in these calculations
0435         for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
0436           int depth = i->id().depth();
0437           int ieta = i->id().ieta();
0438           // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
0439           if (depth == 2 || abs(ieta) == 29)
0440             hfPET_->HFSetFlagFromPET(*i, *rec, myqual, mySeverity);
0441         }
0442 
0443         // Step 2:  Set S8S1 flag (short fibers or |ieta|==29)
0444         for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
0445           int depth = i->id().depth();
0446           int ieta = i->id().ieta();
0447           // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
0448           if (depth == 2 || abs(ieta) == 29)
0449             hfS8S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
0450         }
0451 
0452         // Set 3:  Set S9S1 flag (long fibers)
0453         for (HFRecHitCollection::iterator i = rec->begin(); i != rec->end(); ++i) {
0454           int depth = i->id().depth();
0455           int ieta = i->id().ieta();
0456           // Short fibers and all channels at |ieta|=29 use PET settings in Algo 3
0457           if (depth == 1 && abs(ieta) != 29)
0458             hfS9S1_->HFSetFlagFromS9S1(*i, *rec, myqual, mySeverity);
0459         }
0460       }
0461 
0462       // return result
0463       e.put(std::move(rec));
0464     } else if (subdet_ == HcalOther && subdetOther_ == HcalCalibration) {
0465       edm::Handle<HcalCalibDigiCollection> digi;
0466       e.getByToken(tok_calib_, digi);
0467 
0468       // create empty output
0469       auto rec = std::make_unique<HcalCalibRecHitCollection>();
0470       rec->reserve(digi->size());
0471       // run the algorithm
0472       int first = firstSample_;
0473       int toadd = samplesToAdd_;
0474 
0475       HcalCalibDigiCollection::const_iterator i;
0476       for (i = digi->begin(); i != digi->end(); i++) {
0477         HcalCalibDetId cell = i->id();
0478         //  HcalDetId cellh = i->id();
0479         DetId detcell = (DetId)cell;
0480         // check on cells to be ignored and dropped: (rof,20.Feb.09)
0481         const HcalChannelStatus* mydigistatus = myqual->getValues(detcell.rawId());
0482         if (mySeverity->dropChannel(mydigistatus->getValue()))
0483           continue;
0484         if (dropZSmarkedPassed_)
0485           if (i->zsMarkAndPass())
0486             continue;
0487 
0488         const HcalCalibrations& calibrations = conditions->getHcalCalibrations(cell);
0489         const HcalQIECoder* channelCoder = conditions->getHcalCoder(cell);
0490         const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
0491         HcalCoderDb coder(*channelCoder, *shape);
0492 
0493         // firstSample & samplesToAdd
0494         if (tsFromDB_) {
0495           const HcalRecoParam* param_ts = paramTS_->getValues(detcell.rawId());
0496           first = param_ts->firstSample();
0497           toadd = param_ts->samplesToAdd();
0498         }
0499         rec->push_back(reco_.reconstruct(*i, first, toadd, coder, calibrations));
0500 
0501         /*
0502       // Flag setting not available for calibration rechits
0503     // Set auxiliary flag
0504     int auxflag=0;
0505         int fTS = firstAuxTS_;
0506     for (int xx=fTS; xx<fTS+4 && xx<i->size();++xx)
0507       auxflag+=(i->sample(xx).adc())<<(7*(xx-fTS)); // store the time slices in the first 28 bits of aux, a set of 4 7-bit adc values
0508     // bits 28 and 29 are reserved for capid of the first time slice saved in aux
0509     auxflag+=((i->sample(fTS).capid())<<28);
0510     (rec->back()).setAux(auxflag);
0511 
0512     (rec->back()).setFlags(0); // Not yet implemented for HcalCalibRecHit
0513     */
0514       }
0515       // return result
0516       e.put(std::move(rec));
0517     }
0518   }
0519   //DL  delete myqual;
0520 }  // void HcalHitReconstructor::produce(...)