Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-10 02:59:12

0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "FWCore/ServiceRegistry/interface/Service.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "FWCore/Framework/interface/ConsumesCollector.h"
0007 #include "FWCore/Framework/interface/ProducesCollector.h"
0008 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0009 #include "SimGeneral/MixingModule/interface/PileUpEventPrincipal.h"
0010 
0011 //Data Formats
0012 #include "DataFormats/Common/interface/DetSetVector.h"
0013 #include "DataFormats/Common/interface/DetSet.h"
0014 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0015 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"  // not really needed
0016 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
0017 
0018 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0019 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0020 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0021 
0022 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0023 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0024 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0025 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0026 #include "Geometry/CommonDetUnit/interface/PixelGeomDetType.h"
0027 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0028 #include "CondFormats/SiPixelObjects/interface/PixelIndices.h"
0029 
0030 #include "CLHEP/Random/RandFlat.h"
0031 
0032 #include "SimGeneral/PreMixingModule/interface/PreMixingWorker.h"
0033 #include "SimGeneral/PreMixingModule/interface/PreMixingWorkerFactory.h"
0034 
0035 #include "SiPixelDigitizerAlgorithm.h"
0036 
0037 #include <map>
0038 #include <memory>
0039 
0040 class PreMixingSiPixelWorker : public PreMixingWorker {
0041 public:
0042   PreMixingSiPixelWorker(const edm::ParameterSet& ps, edm::ProducesCollector, edm::ConsumesCollector&& iC);
0043   ~PreMixingSiPixelWorker() override = default;
0044 
0045   void initializeEvent(edm::Event const& e, edm::EventSetup const& c) override;
0046   void addSignals(edm::Event const& e, edm::EventSetup const& es) override;
0047   void addPileups(PileUpEventPrincipal const&, edm::EventSetup const& es) override;
0048   void put(edm::Event& e, edm::EventSetup const& iSetup, std::vector<PileupSummaryInfo> const& ps, int bs) override;
0049 
0050 private:
0051   edm::InputTag pixeldigi_collectionSig_;   // secondary name given to collection of SiPixel digis
0052   edm::InputTag pixeldigi_collectionPile_;  // secondary name given to collection of SiPixel digis
0053   edm::InputTag pixeldigi_extraInfo_;       // secondary name given to collection of SiPixel digis
0054   std::string PixelDigiCollectionDM_;       // secondary name to be given to new SiPixel digis
0055 
0056   edm::EDGetTokenT<edm::DetSetVector<PixelDigi>> PixelDigiToken_;   // Token to retrieve information
0057   edm::EDGetTokenT<edm::DetSetVector<PixelDigi>> PixelDigiPToken_;  // Token to retrieve information
0058   edm::EDGetTokenT<edm::DetSetVector<PixelSimHitExtraInfo>> PixelDigiPExtraToken_;
0059   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoToken_;
0060   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> pDDToken_;
0061 
0062   SiPixelDigitizerAlgorithm digitizer_;
0063 
0064   //
0065   // Internal typedefs
0066 
0067   typedef int Amplitude;
0068   typedef std::map<int, Amplitude, std::less<int>> signal_map_type;  // from Digi.Skel.
0069   typedef std::map<uint32_t, signal_map_type> signalMaps;
0070 
0071   typedef std::multimap<int, PixelDigi>
0072       OneDetectorMap;  // maps by pixel ID for later combination - can have duplicate pixels
0073   typedef std::map<uint32_t, OneDetectorMap> SiGlobalIndex;  // map to all data for each detector ID
0074   typedef std::multimap<int, PixelSimHitExtraInfo> OneExtraInfoMap;
0075   typedef std::map<uint32_t, OneExtraInfoMap> SiPixelExtraInfo;
0076 
0077   SiGlobalIndex SiHitStorage_;
0078   SiPixelExtraInfo SiHitExtraStorage_;
0079 
0080   bool firstInitializeEvent_ = true;
0081   bool firstFinalizeEvent_ = true;
0082   bool applyLateReweighting_;
0083 };
0084 
0085 // Constructor
0086 PreMixingSiPixelWorker::PreMixingSiPixelWorker(const edm::ParameterSet& ps,
0087                                                edm::ProducesCollector producesCollector,
0088                                                edm::ConsumesCollector&& iC)
0089     : tTopoToken_(iC.esConsumes()),
0090       pDDToken_(iC.esConsumes(edm::ESInputTag("", ps.getParameter<std::string>("PixGeometryType")))),
0091       digitizer_(ps, iC) {
0092   // declare the products to produce
0093 
0094   pixeldigi_collectionSig_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionSig");
0095   pixeldigi_collectionPile_ = ps.getParameter<edm::InputTag>("pixeldigiCollectionPile");
0096   pixeldigi_extraInfo_ = ps.getParameter<edm::InputTag>("pixeldigiExtraCollectionPile");
0097   PixelDigiCollectionDM_ = ps.getParameter<std::string>("PixelDigiCollectionDM");
0098   applyLateReweighting_ = ps.getParameter<bool>("applyLateReweighting");
0099   LogDebug("PreMixingSiPixelWorker") << "applyLateReweighting_ in PreMixingSiPixelWorker  " << applyLateReweighting_;
0100 
0101   PixelDigiToken_ = iC.consumes<edm::DetSetVector<PixelDigi>>(pixeldigi_collectionSig_);
0102   PixelDigiPToken_ = iC.consumes<edm::DetSetVector<PixelDigi>>(pixeldigi_collectionPile_);
0103   PixelDigiPExtraToken_ = iC.consumes<edm::DetSetVector<PixelSimHitExtraInfo>>(pixeldigi_extraInfo_);
0104 
0105   producesCollector.produces<edm::DetSetVector<PixelDigi>>(PixelDigiCollectionDM_);
0106   producesCollector.produces<PixelFEDChannelCollection>(PixelDigiCollectionDM_);
0107 
0108   // clear local storage for this event
0109   SiHitStorage_.clear();
0110   SiHitExtraStorage_.clear();
0111 }
0112 
0113 // Need an event initialization
0114 
0115 void PreMixingSiPixelWorker::initializeEvent(edm::Event const& e, edm::EventSetup const& iSetup) {
0116   if (firstInitializeEvent_) {
0117     digitizer_.init(iSetup);
0118     firstInitializeEvent_ = false;
0119   }
0120   digitizer_.initializeEvent();
0121 }
0122 
0123 void PreMixingSiPixelWorker::addSignals(edm::Event const& e, edm::EventSetup const& es) {
0124   // fill in maps of hits
0125 
0126   LogDebug("PreMixingSiPixelWorker") << "===============> adding MC signals for " << e.id();
0127 
0128   edm::Handle<edm::DetSetVector<PixelDigi>> input;
0129 
0130   if (e.getByToken(PixelDigiToken_, input)) {
0131     //loop on all detsets (detectorIDs) inside the input collection
0132     edm::DetSetVector<PixelDigi>::const_iterator DSViter = input->begin();
0133     for (; DSViter != input->end(); DSViter++) {
0134 #ifdef DEBUG
0135       LogDebug("PreMixingSiPixelWorker") << "Processing DetID " << DSViter->id;
0136 #endif
0137 
0138       uint32_t detID = DSViter->id;
0139       edm::DetSet<PixelDigi>::const_iterator begin = (DSViter->data).begin();
0140       edm::DetSet<PixelDigi>::const_iterator end = (DSViter->data).end();
0141       edm::DetSet<PixelDigi>::const_iterator icopy;
0142 
0143       OneDetectorMap LocalMap;
0144 
0145       for (icopy = begin; icopy != end; icopy++) {
0146         LocalMap.insert(OneDetectorMap::value_type((icopy->channel()), *icopy));
0147       }
0148 
0149       SiHitStorage_.insert(SiGlobalIndex::value_type(detID, LocalMap));
0150     }
0151   }
0152 }  // end of addSiPixelSignals
0153 
0154 void PreMixingSiPixelWorker::addPileups(PileUpEventPrincipal const& pep, edm::EventSetup const& es) {
0155   LogDebug("PreMixingSiPixelWorker") << "\n===============> adding pileups from event  " << pep.principal().id()
0156                                      << " for bunchcrossing " << pep.bunchCrossing();
0157 
0158   // fill in maps of hits; same code as addSignals, except now applied to the pileup events
0159 
0160   edm::Handle<edm::DetSetVector<PixelDigi>> inputHandle;
0161   pep.getByLabel(pixeldigi_collectionPile_, inputHandle);
0162 
0163   // added for the Late CR
0164   edm::Handle<edm::DetSetVector<PixelSimHitExtraInfo>> pixelAddInfo;
0165   pep.getByLabel(pixeldigi_extraInfo_, pixelAddInfo);
0166   const TrackerTopology* tTopo = &es.getData(tTopoToken_);
0167   auto const& pDD = es.getData(pDDToken_);
0168   edm::Service<edm::RandomNumberGenerator> rng;
0169   CLHEP::HepRandomEngine* engine = &rng->getEngine(pep.principal().streamID());
0170 
0171   if (inputHandle.isValid()) {
0172     const auto& input = *inputHandle;
0173 
0174     bool loadExtraInformation = false;
0175 
0176     if (pixelAddInfo.isValid() && applyLateReweighting_) {
0177       // access the extra information
0178       loadExtraInformation = true;
0179       // Iterate on detector units
0180       edm::DetSetVector<PixelSimHitExtraInfo>::const_iterator detIdIter;
0181       for (detIdIter = pixelAddInfo->begin(); detIdIter != pixelAddInfo->end(); detIdIter++) {
0182         uint32_t detid = detIdIter->id;  // = rawid
0183         OneExtraInfoMap LocalExtraMap;
0184         edm::DetSet<PixelSimHitExtraInfo>::const_iterator di;
0185         for (di = detIdIter->data.begin(); di != detIdIter->data.end(); di++) {
0186           LocalExtraMap.insert(OneExtraInfoMap::value_type((di->hitIndex()), *di));
0187         }
0188         SiHitExtraStorage_.insert(SiPixelExtraInfo::value_type(detid, LocalExtraMap));
0189       }  // end loop on detIdIter
0190     }  // end if applyLateReweighting_
0191     else if (!pixelAddInfo.isValid() && applyLateReweighting_) {
0192       edm::LogError("PreMixingSiPixelWorker") << " Problem in accessing the Extra Pixel SimHit Collection  !!!! ";
0193       edm::LogError("PreMixingSiPixelWorker") << " The Late Charge Reweighting can not be applied ";
0194       throw cms::Exception("PreMixingSiPixelWorker")
0195           << " Problem in accessing the Extra Pixel SimHit Collection for Late Charge Reweighting \n";
0196     }
0197 
0198     //loop on all detsets (detectorIDs) inside the input collection
0199     edm::DetSetVector<PixelDigi>::const_iterator DSViter = input.begin();
0200     for (; DSViter != input.end(); DSViter++) {
0201 #ifdef DEBUG
0202       LogDebug("PreMixingSiPixelWorker") << "Pileups: Processing DetID " << DSViter->id;
0203 #endif
0204 
0205       uint32_t detID = DSViter->id;
0206       edm::DetSet<PixelDigi>::const_iterator begin = (DSViter->data).begin();
0207       edm::DetSet<PixelDigi>::const_iterator end = (DSViter->data).end();
0208       edm::DetSet<PixelDigi>::const_iterator icopy;
0209 
0210       // find correct local map (or new one) for this detector ID
0211 
0212       SiGlobalIndex::const_iterator itest;
0213 
0214       itest = SiHitStorage_.find(detID);
0215 
0216       std::vector<PixelDigi> TempDigis;
0217       for (icopy = begin; icopy != end; icopy++) {
0218         TempDigis.push_back(*icopy);
0219       }
0220       if (loadExtraInformation) {
0221         // apply the Late Charge Reweighthing on Pile-up digi
0222         SiPixelExtraInfo::const_iterator jtest;
0223         jtest = SiHitExtraStorage_.find(detID);
0224         OneExtraInfoMap LocalSimHitExtraMap = jtest->second;
0225         std::vector<PixelSimHitExtraInfo> TempSimExtra;
0226         for (auto& iLocal : LocalSimHitExtraMap) {
0227           TempSimExtra.push_back(iLocal.second);
0228         }
0229 
0230         for (const auto& iu : pDD.detUnits()) {
0231           if (iu->type().isTrackerPixel()) {
0232             uint32_t detIDinLoop = iu->geographicalId().rawId();
0233             if (detIDinLoop == detID) {
0234               digitizer_.lateSignalReweight(
0235                   dynamic_cast<const PixelGeomDetUnit*>(iu), TempDigis, TempSimExtra, tTopo, engine);
0236               break;
0237             }
0238           }
0239         }
0240       }
0241 
0242       if (itest != SiHitStorage_.end()) {  // this detID already has hits, add to existing map
0243 
0244         OneDetectorMap LocalMap = itest->second;
0245 
0246         // fill in local map with extra channels
0247         for (unsigned int ij = 0; ij < TempDigis.size(); ij++) {
0248           LocalMap.insert(OneDetectorMap::value_type((TempDigis[ij].channel()), TempDigis[ij]));
0249         }
0250 
0251         SiHitStorage_[detID] = LocalMap;
0252 
0253       } else {  // fill local storage with this information, put in global collection
0254 
0255         OneDetectorMap LocalMap;
0256 
0257         for (unsigned int ij = 0; ij < TempDigis.size(); ij++) {
0258           LocalMap.insert(OneDetectorMap::value_type((TempDigis[ij].channel()), TempDigis[ij]));
0259         }
0260 
0261         SiHitStorage_.insert(SiGlobalIndex::value_type(detID, LocalMap));
0262       }
0263     }
0264   }
0265 }
0266 
0267 void PreMixingSiPixelWorker::put(edm::Event& e,
0268                                  edm::EventSetup const& iSetup,
0269                                  std::vector<PileupSummaryInfo> const& ps,
0270                                  int bs) {
0271   // collection of Digis to put in the event
0272 
0273   std::vector<edm::DetSet<PixelDigi>> vPixelDigi;
0274 
0275   // loop through our collection of detectors, merging hits and putting new ones in the output
0276   signalMaps signal;
0277 
0278   // big loop over Detector IDs:
0279 
0280   for (SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin(); IDet != SiHitStorage_.end(); IDet++) {
0281     uint32_t detID = IDet->first;
0282 
0283     OneDetectorMap LocalMap = IDet->second;
0284 
0285     signal_map_type Signals;
0286     Signals.clear();
0287 
0288     //counter variables
0289     int formerPixel = -1;
0290     int currentPixel;
0291     int ADCSum = 0;
0292 
0293     OneDetectorMap::const_iterator iLocalchk;
0294 
0295     for (OneDetectorMap::const_iterator iLocal = LocalMap.begin(); iLocal != LocalMap.end(); ++iLocal) {
0296       currentPixel = iLocal->first;
0297 
0298       if (currentPixel == formerPixel) {  // we have to add these digis together
0299         ADCSum += (iLocal->second).adc();
0300       } else {
0301         if (formerPixel != -1) {  // ADC info stolen from SiStrips...
0302           if (ADCSum > 511)
0303             ADCSum = 255;
0304           else if (ADCSum > 253 && ADCSum < 512)
0305             ADCSum = 254;
0306 
0307           Signals.insert(std::make_pair(formerPixel, ADCSum));
0308         }
0309         // save pointers for next iteration
0310         formerPixel = currentPixel;
0311         ADCSum = (iLocal->second).adc();
0312       }
0313 
0314       iLocalchk = iLocal;
0315       if ((++iLocalchk) == LocalMap.end()) {  //make sure not to lose the last one
0316         if (ADCSum > 511)
0317           ADCSum = 255;
0318         else if (ADCSum > 253 && ADCSum < 512)
0319           ADCSum = 254;
0320         Signals.insert(std::make_pair(formerPixel, ADCSum));
0321       }
0322 
0323     }  // end of loop over one detector
0324 
0325     // stick this into the global vector of detector info
0326     signal.insert(std::make_pair(detID, Signals));
0327 
0328   }  // end of big loop over all detector IDs
0329 
0330   // put the collection of digis in the event
0331   edm::LogInfo("PreMixingSiPixelWorker") << "total # Merged Pixels: " << signal.size();
0332 
0333   std::vector<edm::DetSet<PixelDigi>> theDigiVector;
0334 
0335   // Load inefficiency constants (1st pass), set pileup information.
0336   if (firstFinalizeEvent_) {
0337     digitizer_.init_DynIneffDB(iSetup);
0338     firstFinalizeEvent_ = false;
0339   }
0340 
0341   digitizer_.calculateInstlumiFactor(ps, bs);
0342   digitizer_.setSimAccumulator(signal);
0343 
0344   edm::Service<edm::RandomNumberGenerator> rng;
0345   CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
0346 
0347   auto const& pDD = iSetup.getData(pDDToken_);
0348   const TrackerTopology* tTopo = &iSetup.getData(tTopoToken_);
0349 
0350   if (digitizer_.killBadFEDChannels()) {
0351     std::unique_ptr<PixelFEDChannelCollection> PixelFEDChannelCollection_ = digitizer_.chooseScenario(ps, engine);
0352     if (PixelFEDChannelCollection_ == nullptr) {
0353       throw cms::Exception("NullPointerError") << "PixelFEDChannelCollection not set in chooseScenario function.\n";
0354     }
0355     e.put(std::move(PixelFEDChannelCollection_), PixelDigiCollectionDM_);
0356   }
0357 
0358   for (const auto& iu : pDD.detUnits()) {
0359     if (iu->type().isTrackerPixel()) {
0360       edm::DetSet<PixelDigi> collector(iu->geographicalId().rawId());
0361       edm::DetSet<PixelDigiSimLink> linkcollector(
0362           iu->geographicalId().rawId());                // ignored as DigiSimLinks are combined separately
0363       std::vector<PixelDigiAddTempInfo> tempcollector;  // to be ignored ?
0364 
0365       digitizer_.digitize(
0366           dynamic_cast<const PixelGeomDetUnit*>(iu), collector.data, linkcollector.data, tempcollector, tTopo, engine);
0367       if (!collector.data.empty()) {
0368         theDigiVector.push_back(std::move(collector));
0369       }
0370     }
0371   }
0372 
0373   e.put(std::make_unique<edm::DetSetVector<PixelDigi>>(theDigiVector), PixelDigiCollectionDM_);
0374 
0375   // clear local storage for this event
0376   SiHitStorage_.clear();
0377   SiHitExtraStorage_.clear();
0378 }
0379 
0380 DEFINE_PREMIXING_WORKER(PreMixingSiPixelWorker);