Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:32

0001 // File: DataMixingHcalDigiWorker.cc
0002 // Description:  see DataMixingHcalDigiWorker.h
0003 // Author:  Mike Hildreth, University of Notre Dame
0004 //
0005 //--------------------------------------------
0006 
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "FWCore/Framework/interface/ESHandle.h"
0009 #include "FWCore/Framework/interface/EventSetup.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include <map>
0012 #include <memory>
0013 // calibration headers, for future reference
0014 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
0015 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
0016 //
0017 //
0018 #include "DataMixingHcalDigiWorker.h"
0019 
0020 using namespace std;
0021 
0022 namespace {
0023 
0024   typedef std::multimap<DetId, CaloSamples> HcalDigiMap;
0025 
0026   template <class DIGI>
0027   void convertFc2adc(const CaloSamples &fC, const HcalDbService &conditions, DIGI &digi, int capIdOffset = 0) {
0028     HcalDetId id = fC.id();
0029     const HcalQIECoder *channelCoder = conditions.getHcalCoder(id);
0030     const HcalQIEShape *shape = conditions.getHcalShape(channelCoder);
0031     HcalCoderDb coder(*channelCoder, *shape);
0032     coder.fC2adc(fC, digi, capIdOffset);
0033   }
0034 
0035   template <class DIGI>
0036   void convertAdc2fChelper(const DIGI &digi, const HcalDbService &conditions, CaloSamples &fC, const HcalDetId &id) {
0037     const HcalCalibrations &calib = conditions.getHcalCalibrations(id);
0038     for (int i = 0; i < digi.size(); ++i) {
0039       int capId(digi.sample(i).capid());
0040       fC[i] -= calib.pedestal(capId);
0041     }
0042   }
0043 
0044   template <>
0045   void convertAdc2fChelper<QIE10DataFrame>(const QIE10DataFrame &digi,
0046                                            const HcalDbService &conditions,
0047                                            CaloSamples &fC,
0048                                            const HcalDetId &id) {
0049     const HcalCalibrations &calib = conditions.getHcalCalibrations(id);
0050     for (int i = 0; i < digi.samples(); ++i) {
0051       int capId(digi[i].capid());
0052       fC[i] -= calib.pedestal(capId);
0053     }
0054   }
0055 
0056   template <>
0057   void convertAdc2fChelper<QIE11DataFrame>(const QIE11DataFrame &digi,
0058                                            const HcalDbService &conditions,
0059                                            CaloSamples &fC,
0060                                            const HcalDetId &id) {
0061     const HcalCalibrations &calib = conditions.getHcalCalibrations(id);
0062     for (int i = 0; i < digi.samples(); ++i) {
0063       int capId(digi[i].capid());
0064       fC[i] -= calib.pedestal(capId);
0065     }
0066   }
0067 
0068   template <class DIGI>
0069   void convertAdc2fC(const DIGI &digi, const HcalDbService &conditions, bool keepPedestals, CaloSamples &fC) {
0070     HcalDetId id(digi.id());
0071     const HcalQIECoder *channelCoder = conditions.getHcalCoder(id);
0072     const HcalQIEShape *shape = conditions.getHcalShape(channelCoder);
0073     HcalCoderDb coder(*channelCoder, *shape);
0074     coder.adc2fC(digi, fC);
0075     if (!keepPedestals) {  // subtract pedestals
0076       convertAdc2fChelper(digi, conditions, fC, id);
0077     }
0078   }
0079 
0080   template <class DIGIS>
0081   void convertHcalDigis(const DIGIS &digis, const HcalDbService &conditions, bool keepPedestals, HcalDigiMap &map) {
0082     for (auto digi = digis.begin(); digi != digis.end(); ++digi) {
0083       CaloSamples fC;
0084       convertAdc2fC(*digi, conditions, keepPedestals, fC);
0085       if (!keepPedestals && map.find(digi->id()) == map.end()) {
0086         edm::LogWarning("DataMixingHcalDigiWorker")
0087             << "No signal hits found for HCAL cell " << digi->id() << " Pedestals may be lost for mixed hit";
0088       }
0089       map.insert(HcalDigiMap::value_type(digi->id(), fC));
0090     }
0091   }
0092 
0093   template <>
0094   void convertHcalDigis<QIE10DigiCollection>(const QIE10DigiCollection &digis,
0095                                              const HcalDbService &conditions,
0096                                              bool keepPedestals,
0097                                              HcalDigiMap &map) {
0098     for (auto digiItr = digis.begin(); digiItr != digis.end(); ++digiItr) {
0099       QIE10DataFrame digi(*digiItr);
0100       CaloSamples fC;
0101       convertAdc2fC(digi, conditions, keepPedestals, fC);
0102       if (!keepPedestals && map.find(digi.id()) == map.end()) {
0103         edm::LogWarning("DataMixingHcalDigiWorker")
0104             << "No signal hits found for HCAL cell " << digi.id() << " Pedestals may be lost for mixed hit";
0105       }
0106       map.insert(HcalDigiMap::value_type(digi.id(), fC));
0107     }
0108   }
0109 
0110   template <>
0111   void convertHcalDigis<QIE11DigiCollection>(const QIE11DigiCollection &digis,
0112                                              const HcalDbService &conditions,
0113                                              bool keepPedestals,
0114                                              HcalDigiMap &map) {
0115     for (auto digiItr = digis.begin(); digiItr != digis.end(); ++digiItr) {
0116       QIE11DataFrame digi(*digiItr);
0117       CaloSamples fC;
0118       convertAdc2fC(digi, conditions, keepPedestals, fC);
0119       if (!keepPedestals && map.find(digi.id()) == map.end()) {
0120         edm::LogWarning("DataMixingHcalDigiWorker")
0121             << "No signal hits found for HCAL cell " << digi.id() << " Pedestals may be lost for mixed hit";
0122       }
0123       map.insert(HcalDigiMap::value_type(digi.id(), fC));
0124     }
0125   }
0126 
0127   template <class DIGIS>
0128   bool convertSignalHcalDigis(const edm::Event &e,
0129                               const edm::EDGetTokenT<DIGIS> &token,
0130                               const HcalDbService &conditions,
0131                               HcalDigiMap &map) {
0132     edm::Handle<DIGIS> digis;
0133     if (!e.getByToken(token, digis))
0134       return false;
0135     convertHcalDigis(*digis, conditions, true, map);  // keep pedestals
0136     return true;
0137   }
0138 
0139   template <class DIGIS>
0140   bool convertPileupHcalDigis(const edm::EventPrincipal &ep,
0141                               const edm::InputTag &tag,
0142                               const edm::ModuleCallingContext *mcc,
0143                               const HcalDbService &conditions,
0144                               HcalDigiMap &map) {
0145     auto digis = edm::getProductByTag<DIGIS>(ep, tag, mcc);
0146     if (!digis)
0147       return false;
0148     convertHcalDigis(*(digis->product()), conditions, false,
0149                      map);  // subtract pedestals
0150     return true;
0151   }
0152 
0153   template <class DIGIS>
0154   void buildHcalDigisHelper(DIGIS &digis,
0155                             const DetId &formerID,
0156                             CaloSamples &resultSample,
0157                             const HcalDbService &conditions) {
0158     // make new digi
0159     digis.push_back(typename DIGIS::value_type(formerID));
0160     convertFc2adc(resultSample, conditions, digis.back(),
0161                   0);  // FR guess: simulation starts with 0
0162   }
0163 
0164   template <>
0165   void buildHcalDigisHelper<QIE10DigiCollection>(QIE10DigiCollection &digis,
0166                                                  const DetId &formerID,
0167                                                  CaloSamples &resultSample,
0168                                                  const HcalDbService &conditions) {
0169     // make new digi
0170     digis.push_back(formerID.rawId());
0171     QIE10DataFrame digi(digis.back());
0172     convertFc2adc(resultSample, conditions, digi,
0173                   0);  // FR guess: simulation starts with 0
0174   }
0175 
0176   template <>
0177   void buildHcalDigisHelper<QIE11DigiCollection>(QIE11DigiCollection &digis,
0178                                                  const DetId &formerID,
0179                                                  CaloSamples &resultSample,
0180                                                  const HcalDbService &conditions) {
0181     // make new digi
0182     digis.push_back(formerID.rawId());
0183     QIE11DataFrame digi(digis.back());
0184     convertFc2adc(resultSample, conditions, digi,
0185                   0);  // FR guess: simulation starts with 0
0186   }
0187 
0188   template <class DIGIS>
0189   std::unique_ptr<DIGIS> buildHcalDigis(const HcalDigiMap &map, const HcalDbService &conditions) {
0190     std::unique_ptr<DIGIS> digis(new DIGIS);
0191     // loop over the maps we have, re-making individual hits or digis if
0192     // necessary.
0193     DetId formerID = 0;
0194     CaloSamples resultSample;
0195 
0196     for (auto hit = map.begin(); hit != map.end(); ++hit) {
0197       DetId currentID = hit->first;
0198       const CaloSamples &hitSample = hit->second;
0199 
0200       if (currentID == formerID) {  // accumulating hits
0201         // loop over digi samples in each CaloSample
0202         unsigned int sizenew = (hitSample).size();
0203         unsigned int sizeold = resultSample.size();
0204         if (sizenew > sizeold) {  // extend sample
0205           resultSample.setSize(sizenew);
0206         }
0207         for (unsigned int isamp = 0; isamp < sizenew; isamp++) {  // add new values
0208           resultSample[isamp] += hitSample[isamp];                // for debugging below
0209         }
0210       }
0211       auto hit1 = hit;
0212       bool lastEntry = (++hit1 == map.end());
0213       if (currentID != formerID || lastEntry) {  // store current digi
0214         if (formerID > 0 || lastEntry) {
0215           // make new digi
0216           buildHcalDigisHelper(*digis, formerID, resultSample, conditions);
0217         }
0218         // reset pointers for next iteration
0219         formerID = currentID;
0220         resultSample = hitSample;
0221       }
0222     }
0223     return digis;
0224   }
0225 
0226 }  // namespace
0227 
0228 namespace edm {
0229 
0230   // Virtual constructor
0231 
0232   DataMixingHcalDigiWorker::DataMixingHcalDigiWorker() {}
0233 
0234   // Constructor
0235   DataMixingHcalDigiWorker::DataMixingHcalDigiWorker(const edm::ParameterSet &ps, edm::ConsumesCollector &&iC)
0236       : label_(ps.getParameter<std::string>("Label")) {
0237     // get the subdetector names
0238     //    this->getSubdetectorNames();  //something like this may be useful to
0239     //    check what we are supposed to do...
0240 
0241     // declare the products to produce
0242 
0243     // Hcal
0244 
0245     HBHEdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HBHEdigiCollectionSig");
0246     HOdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HOdigiCollectionSig");
0247     HFdigiCollectionSig_ = ps.getParameter<edm::InputTag>("HFdigiCollectionSig");
0248     ZDCdigiCollectionSig_ = ps.getParameter<edm::InputTag>("ZDCdigiCollectionSig");
0249     QIE10digiCollectionSig_ = ps.getParameter<edm::InputTag>("QIE10digiCollectionSig");
0250     QIE11digiCollectionSig_ = ps.getParameter<edm::InputTag>("QIE11digiCollectionSig");
0251 
0252     HBHEPileInputTag_ = ps.getParameter<edm::InputTag>("HBHEPileInputTag");
0253     HOPileInputTag_ = ps.getParameter<edm::InputTag>("HOPileInputTag");
0254     HFPileInputTag_ = ps.getParameter<edm::InputTag>("HFPileInputTag");
0255     ZDCPileInputTag_ = ps.getParameter<edm::InputTag>("ZDCPileInputTag");
0256     QIE10PileInputTag_ = ps.getParameter<edm::InputTag>("QIE10PileInputTag");
0257     QIE11PileInputTag_ = ps.getParameter<edm::InputTag>("QIE11PileInputTag");
0258 
0259     HBHEDigiToken_ = iC.consumes<HBHEDigiCollection>(HBHEdigiCollectionSig_);
0260     HODigiToken_ = iC.consumes<HODigiCollection>(HOdigiCollectionSig_);
0261     HFDigiToken_ = iC.consumes<HFDigiCollection>(HFdigiCollectionSig_);
0262     QIE10DigiToken_ = iC.consumes<QIE10DigiCollection>(QIE10digiCollectionSig_);
0263     QIE11DigiToken_ = iC.consumes<QIE11DigiCollection>(QIE11digiCollectionSig_);
0264 
0265     HBHEDigiPToken_ = iC.consumes<HBHEDigiCollection>(HBHEPileInputTag_);
0266     HODigiPToken_ = iC.consumes<HODigiCollection>(HOPileInputTag_);
0267     HFDigiPToken_ = iC.consumes<HFDigiCollection>(HFPileInputTag_);
0268     QIE10DigiPToken_ = iC.consumes<QIE10DigiCollection>(QIE10PileInputTag_);
0269     QIE11DigiPToken_ = iC.consumes<QIE11DigiCollection>(QIE11PileInputTag_);
0270 
0271     DoZDC_ = false;
0272     if (!ZDCPileInputTag_.label().empty())
0273       DoZDC_ = true;
0274 
0275     if (DoZDC_) {
0276       ZDCDigiToken_ = iC.consumes<ZDCDigiCollection>(ZDCdigiCollectionSig_);
0277       ZDCDigiPToken_ = iC.consumes<ZDCDigiCollection>(ZDCPileInputTag_);
0278     }
0279 
0280     HBHEDigiCollectionDM_ = ps.getParameter<std::string>("HBHEDigiCollectionDM");
0281     HODigiCollectionDM_ = ps.getParameter<std::string>("HODigiCollectionDM");
0282     HFDigiCollectionDM_ = ps.getParameter<std::string>("HFDigiCollectionDM");
0283     ZDCDigiCollectionDM_ = ps.getParameter<std::string>("ZDCDigiCollectionDM");
0284     QIE10DigiCollectionDM_ = ps.getParameter<std::string>("QIE10DigiCollectionDM");
0285     QIE11DigiCollectionDM_ = ps.getParameter<std::string>("QIE11DigiCollectionDM");
0286 
0287     dbToken_ = iC.esConsumes<HcalDbService, HcalDbRecord>();
0288   }
0289 
0290   // Virtual destructor needed.
0291   DataMixingHcalDigiWorker::~DataMixingHcalDigiWorker() {}
0292 
0293   void DataMixingHcalDigiWorker::addHcalSignals(const edm::Event &e, const edm::EventSetup &ES) {
0294     // Calibration stuff will look like this:
0295     // get conditions
0296     const auto &conditions = ES.getHandle(dbToken_);
0297 
0298     // fill in maps of hits
0299 
0300     LogInfo("DataMixingHcalDigiWorker") << "===============> adding MC signals for " << e.id();
0301     convertSignalHcalDigis<HBHEDigiCollection>(e, HBHEDigiToken_, *conditions, HBHEDigiStorage_);
0302     convertSignalHcalDigis<HODigiCollection>(e, HODigiToken_, *conditions, HODigiStorage_);
0303     convertSignalHcalDigis<HFDigiCollection>(e, HFDigiToken_, *conditions, HFDigiStorage_);
0304     convertSignalHcalDigis<QIE10DigiCollection>(e, QIE10DigiToken_, *conditions, QIE10DigiStorage_);
0305     convertSignalHcalDigis<QIE11DigiCollection>(e, QIE11DigiToken_, *conditions, QIE11DigiStorage_);
0306 
0307     // ZDC next
0308 
0309     if (DoZDC_) {
0310       Handle<ZDCDigiCollection> pZDCDigis;
0311 
0312       const ZDCDigiCollection *ZDCDigis = nullptr;
0313 
0314       if (e.getByToken(ZDCDigiToken_, pZDCDigis)) {
0315         ZDCDigis = pZDCDigis.product();  // get a ptr to the product
0316 #ifdef DEBUG
0317         LogDebug("DataMixingHcalDigiWorker") << "total # ZDC digis: " << ZDCDigis->size();
0318 #endif
0319       }
0320 
0321       if (ZDCDigis) {
0322         // loop over digis, storing them in a map so we can add pileup later
0323         for (ZDCDigiCollection::const_iterator it = ZDCDigis->begin(); it != ZDCDigis->end(); ++it) {
0324           // calibration, for future reference:  (same block for all Hcal types)
0325           // ZDC is different
0326           HcalZDCDetId cell = it->id();
0327           //         const HcalCalibrations&
0328           //         calibrations=conditions->getHcalCalibrations(cell);
0329           const HcalQIECoder *channelCoder = conditions->getHcalCoder(cell);
0330           const HcalQIEShape *shape = conditions->getHcalShape(channelCoder);  // this one is generic
0331           HcalCoderDb coder(*channelCoder, *shape);
0332 
0333           CaloSamples tool;
0334           coder.adc2fC((*it), tool);
0335 
0336           ZDCDigiStorage_.insert(ZDCDigiMap::value_type((it->id()), tool));
0337 
0338 #ifdef DEBUG
0339           // Commented out because this does not compile anymore
0340           // LogDebug("DataMixingHcalDigiWorker") << "processed ZDCDigi with
0341           // rawId: "
0342           //                                      << it->id() << "\n"
0343           //                                      << " digi energy: " <<
0344           //                                      it->energy();
0345 #endif
0346         }
0347       }
0348     }
0349 
0350   }  // end of addHCalSignals
0351 
0352   void DataMixingHcalDigiWorker::addHcalPileups(const int bcr,
0353                                                 const EventPrincipal *ep,
0354                                                 unsigned int eventNr,
0355                                                 const edm::EventSetup &ES,
0356                                                 ModuleCallingContext const *mcc) {
0357     LogDebug("DataMixingHcalDigiWorker") << "\n===============> adding pileups from event  " << ep->id()
0358                                          << " for bunchcrossing " << bcr;
0359 
0360     // get conditions
0361     const auto &conditions = ES.getHandle(dbToken_);
0362 
0363     convertPileupHcalDigis<HBHEDigiCollection>(*ep, HBHEPileInputTag_, mcc, *conditions, HBHEDigiStorage_);
0364     convertPileupHcalDigis<HODigiCollection>(*ep, HOPileInputTag_, mcc, *conditions, HODigiStorage_);
0365     convertPileupHcalDigis<HFDigiCollection>(*ep, HFPileInputTag_, mcc, *conditions, HFDigiStorage_);
0366     convertPileupHcalDigis<QIE10DigiCollection>(*ep, QIE10PileInputTag_, mcc, *conditions, QIE10DigiStorage_);
0367     convertPileupHcalDigis<QIE11DigiCollection>(*ep, QIE11PileInputTag_, mcc, *conditions, QIE11DigiStorage_);
0368 
0369     // ZDC Next
0370 
0371     if (DoZDC_) {
0372       std::shared_ptr<Wrapper<ZDCDigiCollection> const> ZDCDigisPTR =
0373           getProductByTag<ZDCDigiCollection>(*ep, ZDCPileInputTag_, mcc);
0374 
0375       if (ZDCDigisPTR) {
0376         const ZDCDigiCollection *ZDCDigis = const_cast<ZDCDigiCollection *>(ZDCDigisPTR->product());
0377 
0378         LogDebug("DataMixingHcalDigiWorker") << "total # ZDC digis: " << ZDCDigis->size();
0379 
0380         // loop over digis, adding these to the existing maps
0381         for (ZDCDigiCollection::const_iterator it = ZDCDigis->begin(); it != ZDCDigis->end(); ++it) {
0382           // calibration, for future reference:  (same block for all Hcal types)
0383           // ZDC is different
0384           HcalZDCDetId cell = it->id();
0385           const HcalQIECoder *channelCoder = conditions->getHcalCoder(cell);
0386           const HcalQIEShape *shape = conditions->getHcalShape(channelCoder);  // this one is generic
0387           HcalCoderDb coder(*channelCoder, *shape);
0388 
0389           CaloSamples tool;
0390           coder.adc2fC((*it), tool);
0391 
0392           ZDCDigiStorage_.insert(ZDCDigiMap::value_type((it->id()), tool));
0393 
0394 #ifdef DEBUG
0395           // Commented out because this does not compile anymore
0396           // LogDebug("DataMixingHcalDigiWorker") << "processed ZDCDigi with
0397           // rawId: "
0398           //                                      << it->id() << "\n"
0399           //                                      << " digi energy: " <<
0400           //                                      it->energy();
0401 #endif
0402         }
0403       }
0404     }
0405   }
0406 
0407   void DataMixingHcalDigiWorker::putHcal(edm::Event &e, const edm::EventSetup &ES) {
0408     const auto &conditions = ES.getHandle(dbToken_);
0409 
0410     // collection of digis to put in the event
0411     std::unique_ptr<HBHEDigiCollection> HBHEdigis = buildHcalDigis<HBHEDigiCollection>(HBHEDigiStorage_, *conditions);
0412     std::unique_ptr<HODigiCollection> HOdigis = buildHcalDigis<HODigiCollection>(HODigiStorage_, *conditions);
0413     std::unique_ptr<HFDigiCollection> HFdigis = buildHcalDigis<HFDigiCollection>(HFDigiStorage_, *conditions);
0414     std::unique_ptr<QIE10DigiCollection> QIE10digis =
0415         buildHcalDigis<QIE10DigiCollection>(QIE10DigiStorage_, *conditions);
0416     std::unique_ptr<QIE11DigiCollection> QIE11digis =
0417         buildHcalDigis<QIE11DigiCollection>(QIE11DigiStorage_, *conditions);
0418     std::unique_ptr<ZDCDigiCollection> ZDCdigis(new ZDCDigiCollection);
0419 
0420     // loop over the maps we have, re-making individual hits or digis if
0421     // necessary.
0422     DetId formerID = 0;
0423     DetId currentID;
0424 
0425     double fC_new;
0426     double fC_old;
0427     double fC_sum;
0428 
0429     // ZDC next...
0430 
0431     // loop over the maps we have, re-making individual hits or digis if
0432     // necessary.
0433     formerID = 0;
0434     CaloSamples ZDC_old;
0435 
0436     ZDCDigiMap::const_iterator iZDCchk;
0437 
0438     for (ZDCDigiMap::const_iterator iZDC = ZDCDigiStorage_.begin(); iZDC != ZDCDigiStorage_.end(); ++iZDC) {
0439       currentID = iZDC->first;
0440 
0441       if (currentID == formerID) {  // we have to add these digis together
0442 
0443         // loop over digi samples in each CaloSample
0444         unsigned int sizenew = (iZDC->second).size();
0445         unsigned int sizeold = ZDC_old.size();
0446 
0447         unsigned int max_samp = std::max(sizenew, sizeold);
0448 
0449         CaloSamples ZDC_bigger(currentID, max_samp);
0450 
0451         bool usenew = false;
0452 
0453         if (sizenew > sizeold)
0454           usenew = true;
0455 
0456         // samples from different events can be of different lengths - sum all
0457         // that overlap.
0458 
0459         for (unsigned int isamp = 0; isamp < max_samp; isamp++) {
0460           if (isamp < sizenew) {
0461             fC_new = (iZDC->second)[isamp];
0462           } else {
0463             fC_new = 0;
0464           }
0465 
0466           if (isamp < sizeold) {
0467             fC_old = ZDC_old[isamp];
0468           } else {
0469             fC_old = 0;
0470           }
0471 
0472           // add values
0473           fC_sum = fC_new + fC_old;
0474 
0475           if (usenew) {
0476             ZDC_bigger[isamp] = fC_sum;
0477           } else {
0478             ZDC_old[isamp] = fC_sum;
0479           }  // overwrite old sample, adding new info
0480         }
0481         if (usenew)
0482           ZDC_old = ZDC_bigger;  // save new, larger sized sample in "old" slot
0483 
0484       } else {
0485         if (formerID > 0) {
0486           // make new digi
0487           ZDCdigis->push_back(ZDCDataFrame(formerID));
0488 
0489           // set up information to convert back
0490 
0491           HcalZDCDetId cell = ZDC_old.id();
0492           const HcalQIECoder *channelCoder = conditions->getHcalCoder(cell);
0493           const HcalQIEShape *shape = conditions->getHcalShape(channelCoder);  // this one is generic
0494           HcalCoderDb coder(*channelCoder, *shape);
0495 
0496           unsigned int sizeold = ZDC_old.size();
0497           for (unsigned int isamp = 0; isamp < sizeold; isamp++) {
0498             coder.fC2adc(ZDC_old, (ZDCdigis->back()),
0499                          1);  // as per simulation, capid=1
0500           }
0501         }
0502         // save pointers for next iteration
0503         formerID = currentID;
0504         ZDC_old = iZDC->second;
0505       }
0506 
0507       iZDCchk = iZDC;
0508       if ((++iZDCchk) == ZDCDigiStorage_.end()) {  // make sure not to lose the last one
0509         // make new digi
0510         ZDCdigis->push_back(ZDCDataFrame(currentID));
0511 
0512         // set up information to convert back
0513 
0514         HcalZDCDetId cell = (iZDC->second).id();
0515         const HcalQIECoder *channelCoder = conditions->getHcalCoder(cell);
0516         const HcalQIEShape *shape = conditions->getHcalShape(channelCoder);  // this one is generic
0517         HcalCoderDb coder(*channelCoder, *shape);
0518 
0519         unsigned int sizeold = (iZDC->second).size();
0520         for (unsigned int isamp = 0; isamp < sizeold; isamp++) {
0521           coder.fC2adc(ZDC_old, (ZDCdigis->back()),
0522                        1);  // as per simulation, capid=1
0523         }
0524       }
0525     }
0526 
0527     // done merging
0528 
0529     // put the collection of recunstructed hits in the event
0530     LogInfo("DataMixingHcalDigiWorker") << "total # HBHE Merged digis: " << HBHEdigis->size();
0531     LogInfo("DataMixingHcalDigiWorker") << "total # HO Merged digis: " << HOdigis->size();
0532     LogInfo("DataMixingHcalDigiWorker") << "total # HF Merged digis: " << HFdigis->size();
0533     LogInfo("DataMixingHcalDigiWorker") << "total # QIE10 Merged digis: " << QIE10digis->size();
0534     LogInfo("DataMixingHcalDigiWorker") << "total # QIE11 Merged digis: " << QIE11digis->size();
0535     LogInfo("DataMixingHcalDigiWorker") << "total # ZDC Merged digis: " << ZDCdigis->size();
0536 
0537     e.put(std::move(HBHEdigis), HBHEDigiCollectionDM_);
0538     e.put(std::move(HOdigis), HODigiCollectionDM_);
0539     e.put(std::move(HFdigis), HFDigiCollectionDM_);
0540     e.put(std::move(QIE10digis), QIE10DigiCollectionDM_);
0541     e.put(std::move(QIE11digis), QIE11DigiCollectionDM_);
0542     e.put(std::move(ZDCdigis), ZDCDigiCollectionDM_);
0543 
0544     // clear local storage after this event
0545     HBHEDigiStorage_.clear();
0546     HODigiStorage_.clear();
0547     HFDigiStorage_.clear();
0548     QIE10DigiStorage_.clear();
0549     QIE11DigiStorage_.clear();
0550     ZDCDigiStorage_.clear();
0551   }
0552 
0553 }  // namespace edm