Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // File: DataMixingModule.cc
0002 // Description:  see DataMixingModule.h
0003 // Author:  Mike Hildreth, University of Notre Dame
0004 //
0005 //--------------------------------------------
0006 
0007 #include <functional>
0008 #include <iostream>
0009 #include <map>
0010 #include <memory>
0011 
0012 #include "DataFormats/Common/interface/Handle.h"
0013 #include "FWCore/Framework/interface/ConsumesCollector.h"
0014 #include "FWCore/Framework/interface/ModuleContextSentry.h"
0015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0016 #include "FWCore/ServiceRegistry/interface/InternalContext.h"
0017 #include "FWCore/ServiceRegistry/interface/ModuleCallingContext.h"
0018 #include "FWCore/ServiceRegistry/interface/ParentContext.h"
0019 //
0020 //
0021 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0022 #include "DataMixingModule.h"
0023 #include "SimDataFormats/CrossingFrame/interface/CrossingFramePlaybackInfoNew.h"
0024 #include "SimGeneral/MixingModule/interface/PileUpEventPrincipal.h"
0025 
0026 using namespace std;
0027 
0028 namespace edm {
0029 
0030   // Constructor
0031   DataMixingModule::DataMixingModule(const edm::ParameterSet &ps, MixingCache::Config const *globalConf)
0032       : BMixingModule(ps, globalConf),
0033         EBPileInputTag_(ps.getParameter<edm::InputTag>("EBPileInputTag")),
0034         EEPileInputTag_(ps.getParameter<edm::InputTag>("EEPileInputTag")),
0035         ESPileInputTag_(ps.getParameter<edm::InputTag>("ESPileInputTag")),
0036         HBHEPileInputTag_(ps.getParameter<edm::InputTag>("HBHEPileInputTag")),
0037         HOPileInputTag_(ps.getParameter<edm::InputTag>("HOPileInputTag")),
0038         HFPileInputTag_(ps.getParameter<edm::InputTag>("HFPileInputTag")),
0039         ZDCPileInputTag_(ps.getParameter<edm::InputTag>("ZDCPileInputTag")),
0040         QIE10PileInputTag_(ps.getParameter<edm::InputTag>("QIE10PileInputTag")),
0041         QIE11PileInputTag_(ps.getParameter<edm::InputTag>("QIE11PileInputTag")),
0042         tokDB_(esConsumes<HcalDbService, HcalDbRecord>()),
0043         label_(ps.getParameter<std::string>("Label")) {
0044     // prepare for data access in DataMixingHcalDigiWorkerProd
0045     tok_hbhe_ = consumes<HBHEDigitizerTraits::DigiCollection>(HBHEPileInputTag_);
0046     tok_ho_ = consumes<HODigitizerTraits::DigiCollection>(HOPileInputTag_);
0047     tok_hf_ = consumes<HFDigitizerTraits::DigiCollection>(HFPileInputTag_);
0048     tok_zdc_ = consumes<ZDCDigitizerTraits::DigiCollection>(ZDCPileInputTag_);
0049     tok_qie10_ = consumes<HcalQIE10DigitizerTraits::DigiCollection>(QIE10PileInputTag_);
0050     tok_qie11_ = consumes<HcalQIE11DigitizerTraits::DigiCollection>(QIE11PileInputTag_);
0051 
0052     // get the subdetector names
0053     this->getSubdetectorNames();  // something like this may be useful to check
0054                                   // what we are supposed to do...
0055 
0056     // For now, list all of them here.  Later, make this selectable with input
0057     // parameters
0058     //
0059 
0060     // Check to see if we are working in Full or Fast Simulation
0061 
0062     MergeTrackerDigis_ = (ps.getParameter<std::string>("TrackerMergeType")) == "Digis";
0063     MergeEMDigis_ = (ps.getParameter<std::string>("EcalMergeType")) == "Digis";
0064     MergeHcalDigis_ = (ps.getParameter<std::string>("HcalMergeType")) == "Digis";
0065     if (MergeHcalDigis_)
0066       MergeHcalDigisProd_ = (ps.getParameter<std::string>("HcalDigiMerge") == "FullProd");
0067 
0068     // Put Fast Sim Sequences here for Simplification: Fewer options!
0069 
0070     if (MergeEMDigis_) {
0071       // cout<<"EM Digis TRUE!!!"<<endl;
0072 
0073       EBDigiCollectionDM_ = ps.getParameter<std::string>("EBDigiCollectionDM");
0074       EEDigiCollectionDM_ = ps.getParameter<std::string>("EEDigiCollectionDM");
0075       ESDigiCollectionDM_ = ps.getParameter<std::string>("ESDigiCollectionDM");
0076       //   nMaxPrintout_            =
0077       //   ps.getUntrackedParameter<int>("nMaxPrintout",10);
0078 
0079       produces<EBDigiCollection>(EBDigiCollectionDM_);
0080       produces<EEDigiCollection>(EEDigiCollectionDM_);
0081       produces<ESDigiCollection>(ESDigiCollectionDM_);
0082 
0083       EMDigiWorker_ = new DataMixingEMDigiWorker(ps, consumesCollector());
0084     } else {  // merge RecHits
0085       EBRecHitCollectionDM_ = ps.getParameter<std::string>("EBRecHitCollectionDM");
0086       EERecHitCollectionDM_ = ps.getParameter<std::string>("EERecHitCollectionDM");
0087       ESRecHitCollectionDM_ = ps.getParameter<std::string>("ESRecHitCollectionDM");
0088       //   nMaxPrintout_            =
0089       //   ps.getUntrackedParameter<int>("nMaxPrintout",10);
0090 
0091       produces<EBRecHitCollection>(EBRecHitCollectionDM_);
0092       produces<EERecHitCollection>(EERecHitCollectionDM_);
0093       produces<ESRecHitCollection>(ESRecHitCollectionDM_);
0094 
0095       EMWorker_ = new DataMixingEMWorker(ps, consumesCollector());
0096     }
0097     // Hcal next
0098 
0099     if (MergeHcalDigis_) {
0100       HBHEDigiCollectionDM_ = ps.getParameter<std::string>("HBHEDigiCollectionDM");
0101       HODigiCollectionDM_ = ps.getParameter<std::string>("HODigiCollectionDM");
0102       HFDigiCollectionDM_ = ps.getParameter<std::string>("HFDigiCollectionDM");
0103       ZDCDigiCollectionDM_ = ps.getParameter<std::string>("ZDCDigiCollectionDM");
0104       QIE10DigiCollectionDM_ = ps.getParameter<std::string>("QIE10DigiCollectionDM");
0105       QIE11DigiCollectionDM_ = ps.getParameter<std::string>("QIE11DigiCollectionDM");
0106 
0107       produces<HBHEDigiCollection>();
0108       produces<HODigiCollection>();
0109       produces<HFDigiCollection>();
0110       produces<ZDCDigiCollection>();
0111 
0112       produces<QIE10DigiCollection>("HFQIE10DigiCollection");
0113       produces<QIE11DigiCollection>("HBHEQIE11DigiCollection");
0114 
0115       if (ps.getParameter<bool>("debugCaloSamples")) {
0116         produces<CaloSamplesCollection>("HcalSamples");
0117       }
0118       if (ps.getParameter<bool>("injectTestHits")) {
0119         produces<edm::PCaloHitContainer>("HcalHits");
0120       }
0121 
0122       if (MergeHcalDigisProd_)
0123         HcalDigiWorkerProd_ = new DataMixingHcalDigiWorkerProd(ps, consumesCollector(), tokDB_);
0124       else
0125         HcalDigiWorker_ = new DataMixingHcalDigiWorker(ps, consumesCollector());
0126 
0127     } else {
0128       HBHERecHitCollectionDM_ = ps.getParameter<std::string>("HBHERecHitCollectionDM");
0129       HORecHitCollectionDM_ = ps.getParameter<std::string>("HORecHitCollectionDM");
0130       HFRecHitCollectionDM_ = ps.getParameter<std::string>("HFRecHitCollectionDM");
0131       ZDCRecHitCollectionDM_ = ps.getParameter<std::string>("ZDCRecHitCollectionDM");
0132 
0133       produces<HBHERecHitCollection>(HBHERecHitCollectionDM_);
0134       produces<HORecHitCollection>(HORecHitCollectionDM_);
0135       produces<HFRecHitCollection>(HFRecHitCollectionDM_);
0136       produces<ZDCRecHitCollection>(ZDCRecHitCollectionDM_);
0137 
0138       HcalWorker_ = new DataMixingHcalWorker(ps, consumesCollector());
0139     }
0140 
0141     // Muons
0142 
0143     DTDigiCollectionDM_ = ps.getParameter<std::string>("DTDigiCollectionDM");
0144     RPCDigiCollectionDM_ = ps.getParameter<std::string>("RPCDigiCollectionDM");
0145     CSCStripDigiCollectionDM_ = ps.getParameter<std::string>("CSCStripDigiCollectionDM");
0146     CSCWireDigiCollectionDM_ = ps.getParameter<std::string>("CSCWireDigiCollectionDM");
0147     CSCComparatorDigiCollectionDM_ = ps.getParameter<std::string>("CSCComparatorDigiCollectionDM");
0148 
0149     produces<DTDigiCollection>();
0150     produces<RPCDigiCollection>();
0151     produces<CSCStripDigiCollection>(CSCStripDigiCollectionDM_);
0152     produces<CSCWireDigiCollection>(CSCWireDigiCollectionDM_);
0153     produces<CSCComparatorDigiCollection>(CSCComparatorDigiCollectionDM_);
0154 
0155     MuonWorker_ = new DataMixingMuonWorker(ps, consumesCollector());
0156 
0157     // Si-Strips
0158 
0159     useSiStripRawDigi_ = ps.exists("SiStripRawDigiSource")
0160                              ? ps.getParameter<std::string>("SiStripRawDigiSource") == "PILEUP" ||
0161                                    ps.getParameter<std::string>("SiStripRawDigiSource") == "SIGNAL"
0162                              : false;
0163 
0164     SiStripDigiCollectionDM_ = ps.getParameter<std::string>("SiStripDigiCollectionDM");
0165 
0166     if (useSiStripRawDigi_) {
0167       produces<edm::DetSetVector<SiStripRawDigi>>(SiStripDigiCollectionDM_);
0168       SiStripRawWorker_ = new DataMixingSiStripRawWorker(ps, consumesCollector());
0169 
0170     } else {
0171       produces<edm::DetSetVector<SiStripDigi>>(SiStripDigiCollectionDM_);
0172 
0173       SiStripWorker_ = new DataMixingSiStripWorker(ps, consumesCollector());
0174     }
0175 
0176     // Pixels
0177 
0178     PixelDigiCollectionDM_ = ps.getParameter<std::string>("PixelDigiCollectionDM");
0179 
0180     produces<edm::DetSetVector<PixelDigi>>(PixelDigiCollectionDM_);
0181 
0182     SiPixelWorker_ = new DataMixingSiPixelWorker(ps, consumesCollector());
0183 
0184     // Pileup Information: if doing pre-mixing, we have to save the pileup
0185     // information from the Secondary stream
0186 
0187     MergePileup_ = ps.getParameter<bool>("MergePileupInfo");
0188 
0189     if (MergePileup_) {
0190       produces<std::vector<PileupSummaryInfo>>();
0191       produces<int>("bunchSpacing");
0192       produces<CrossingFramePlaybackInfoNew>();
0193 
0194       std::vector<edm::InputTag> GenPUProtonsInputTags;
0195       GenPUProtonsInputTags = ps.getParameter<std::vector<edm::InputTag>>("GenPUProtonsInputTags");
0196       for (std::vector<edm::InputTag>::const_iterator it_InputTag = GenPUProtonsInputTags.begin();
0197            it_InputTag != GenPUProtonsInputTags.end();
0198            ++it_InputTag)
0199         produces<std::vector<reco::GenParticle>>(it_InputTag->label());
0200 
0201       PUWorker_ = new DataMixingPileupCopy(ps, consumesCollector());
0202     }
0203   }
0204 
0205   void DataMixingModule::getSubdetectorNames() {
0206     // get subdetector names
0207     // edm::Service<edm::ConstProductRegistry> reg;
0208     // Loop over provenance of products in registry.
0209     // for (edm::ProductRegistry::ProductList::const_iterator it =
0210     // reg->productList().begin(); it != reg->productList().end(); ++it) {
0211 
0212     //  **** Check this out.... ****
0213 
0214     // See FWCore/Framework/interface/BranchDescription.h
0215     // BranchDescription contains all the information for the product.
0216 
0217     // This section not very backwards-compatible in terms of digi-merging.  Need
0218     // to be able to specify here which data format to look at...
0219 
0220     //      edm::BranchDescription desc = it->second;
0221     // if (!desc.friendlyClassName_.compare(0,9,"EBRecHitC")) {
0222     //  Subdetectors_.push_back(desc.productInstanceName_);
0223     // LogInfo("DataMixingModule") <<"Adding container
0224     // "<<desc.productInstanceName_ <<" for pileup treatment";
0225     //}
0226     // else if (!desc.friendlyClassName_.compare(0,9,"EERecHitC")) {
0227     //      else if (!desc.friendlyClassName_.compare(0,9,"EErechitC") &&
0228     //      desc.productInstanceName_.compare(0,11,"TrackerHits")) {
0229     //  Subdetectors_.push_back(desc.productInstanceName_);
0230     // LogInfo("DataMixingModule") <<"Adding container
0231     // "<<desc.productInstanceName_ <<" for pileup treatment";
0232     //}
0233     // else if (!desc.friendlyClassName_.compare(0,9,"HBRecHitC")) {
0234     //  Subdetectors_.push_back(desc.productInstanceName_);
0235     // LogInfo("DataMixingModule") <<"Adding container
0236     // "<<desc.productInstanceName_ <<" for pileup treatment";
0237     //}
0238     // else if (!desc.friendlyClassName_.compare(0,9,"HERecHitC")) {
0239     //  Subdetectors_.push_back(desc.productInstanceName_);
0240     // LogInfo("DataMixingModule") <<"Adding container
0241     // "<<desc.productInstanceName_ <<" for pileup treatment";
0242     // }
0243     // and so on with other detector types...
0244     // }
0245   }
0246 
0247   void DataMixingModule::initializeEvent(const edm::Event &e, const edm::EventSetup &ES) {}
0248 
0249   void DataMixingModule::beginRun(edm::Run const &run, const edm::EventSetup &ES) { BMixingModule::beginRun(run, ES); }
0250 
0251   void DataMixingModule::endRun(edm::Run const &run, const edm::EventSetup &ES) { BMixingModule::endRun(run, ES); }
0252 
0253   // Virtual destructor needed.
0254   DataMixingModule::~DataMixingModule() {
0255     if (MergeEMDigis_) {
0256       delete EMDigiWorker_;
0257     } else {
0258       delete EMWorker_;
0259     }
0260     if (MergeHcalDigis_) {
0261       if (MergeHcalDigisProd_) {
0262         delete HcalDigiWorkerProd_;
0263       } else {
0264         delete HcalDigiWorker_;
0265       }
0266     } else {
0267       delete HcalWorker_;
0268     }
0269     if (MuonWorker_)
0270       delete MuonWorker_;
0271     if (useSiStripRawDigi_)
0272       delete SiStripRawWorker_;
0273     else
0274       delete SiStripWorker_;
0275     delete SiPixelWorker_;
0276     if (MergePileup_) {
0277       delete PUWorker_;
0278     }
0279   }
0280 
0281   void DataMixingModule::addSignals(const edm::Event &e, const edm::EventSetup &ES) {
0282     // fill in maps of hits
0283 
0284     LogDebug("DataMixingModule") << "===============> adding MC signals for " << e.id();
0285 
0286     // Ecal
0287     if (MergeEMDigis_) {
0288       EMDigiWorker_->addEMSignals(e, ES);
0289     } else {
0290       EMWorker_->addEMSignals(e);
0291     }
0292 
0293     // Hcal
0294     if (MergeHcalDigis_) {
0295       if (MergeHcalDigisProd_) {
0296         HcalDigiWorkerProd_->addHcalSignals(e, ES);
0297       } else {
0298         HcalDigiWorker_->addHcalSignals(e, ES);
0299       }
0300     } else {
0301       HcalWorker_->addHcalSignals(e);
0302     }
0303 
0304     // Muon
0305     MuonWorker_->addMuonSignals(e);
0306 
0307     // SiStrips
0308     if (useSiStripRawDigi_)
0309       SiStripRawWorker_->addSiStripSignals(e);
0310     else
0311       SiStripWorker_->addSiStripSignals(e);
0312 
0313     // SiPixels
0314     SiPixelWorker_->addSiPixelSignals(e);
0315     AddedPileup_ = false;
0316   }  // end of addSignals
0317 
0318   bool DataMixingModule::pileWorker(
0319       const EventPrincipal &ep, int bcr, int eventNr, const edm::EventSetup &ES, edm::ModuleCallingContext const *mcc) {
0320     InternalContext internalContext(ep.id(), mcc);
0321     ParentContext parentContext(&internalContext);
0322     ModuleCallingContext moduleCallingContext(&moduleDescription());
0323     ModuleContextSentry moduleContextSentry(&moduleCallingContext, parentContext);
0324 
0325     LogDebug("DataMixingModule") << "\n===============> adding pileups from event  " << ep.id() << " for bunchcrossing "
0326                                  << bcr;
0327 
0328     // Note:  setupPileUpEvent may modify the run and lumi numbers of the
0329     // EventPrincipal to match that of the primary event.
0330     setupPileUpEvent(ES);
0331 
0332     // check and see if we need to copy the pileup information from
0333     // secondary stream to the output stream
0334     // We only have the pileup event here, so pick the first time and store the
0335     // info
0336 
0337     if (MergePileup_ && !AddedPileup_) {
0338       PUWorker_->addPileupInfo(&ep, eventNr, &moduleCallingContext);
0339 
0340       AddedPileup_ = true;
0341     }
0342 
0343     // fill in maps of hits; same code as addSignals, except now applied to the
0344     // pileup events
0345 
0346     // Ecal
0347     if (MergeEMDigis_) {
0348       EMDigiWorker_->addEMPileups(bcr, &ep, eventNr, ES, &moduleCallingContext);
0349     } else {
0350       EMWorker_->addEMPileups(bcr, &ep, eventNr, &moduleCallingContext);
0351     }
0352 
0353     // Hcal
0354     if (MergeHcalDigis_) {
0355       if (MergeHcalDigisProd_) {
0356         HcalDigiWorkerProd_->addHcalPileups(bcr, &ep, eventNr, ES, &moduleCallingContext);
0357       } else {
0358         HcalDigiWorker_->addHcalPileups(bcr, &ep, eventNr, ES, &moduleCallingContext);
0359       }
0360     } else {
0361       HcalWorker_->addHcalPileups(bcr, &ep, eventNr, &moduleCallingContext);
0362     }
0363 
0364     // Muon
0365     MuonWorker_->addMuonPileups(bcr, &ep, eventNr, &moduleCallingContext);
0366 
0367     // SiStrips
0368     if (useSiStripRawDigi_)
0369       SiStripRawWorker_->addSiStripPileups(bcr, &ep, eventNr, &moduleCallingContext);
0370     else
0371       SiStripWorker_->addSiStripPileups(bcr, &ep, eventNr, &moduleCallingContext);
0372 
0373     // SiPixels
0374     // whoops this should be for the MC worker ?????
0375     // SiPixelWorker_->setPileupInfo(ps,bunchSpacing);
0376     SiPixelWorker_->addSiPixelPileups(bcr, &ep, eventNr, &moduleCallingContext);
0377 
0378     return true;
0379   }
0380 
0381   void DataMixingModule::doPileUp(edm::Event &e, const edm::EventSetup &ES) {
0382     using namespace std::placeholders;
0383 
0384     std::vector<edm::SecondaryEventIDAndFileInfo> recordEventID;
0385     std::vector<int> PileupList;
0386     PileupList.clear();
0387     TrueNumInteractions_.clear();
0388 
0389     ModuleCallingContext const *mcc = e.moduleCallingContext();
0390 
0391     for (int bunchCrossing = minBunch_; bunchCrossing <= maxBunch_; ++bunchCrossing) {
0392       for (unsigned int isource = 0; isource < maxNbSources_; ++isource) {
0393         std::shared_ptr<PileUp> source = inputSources_[isource];
0394         if (!source || !(source->doPileUp(bunchCrossing)))
0395           continue;
0396 
0397         if (isource == 0)
0398           source->CalculatePileup(minBunch_, maxBunch_, PileupList, TrueNumInteractions_, e.streamID());
0399 
0400         int NumPU_Events = 0;
0401         if (isource == 0) {
0402           NumPU_Events = PileupList[bunchCrossing - minBunch_];
0403         } else {
0404           // non-minbias pileup only gets one event for now. Fix later if desired.
0405           NumPU_Events = 1;
0406         }
0407 
0408         source->readPileUp(
0409             e.id(),
0410             recordEventID,
0411             std::bind(&DataMixingModule::pileWorker, std::ref(*this), _1, bunchCrossing, _2, std::cref(ES), mcc),
0412             NumPU_Events,
0413             e.streamID());
0414       }
0415     }
0416   }
0417 
0418   void DataMixingModule::put(edm::Event &e, const edm::EventSetup &ES) {
0419     // individual workers...
0420 
0421     // move pileup first so we have access to the information for the put step
0422 
0423     std::vector<PileupSummaryInfo> ps;
0424     int bunchSpacing = 10000;
0425 
0426     if (MergePileup_) {
0427       PUWorker_->getPileupInfo(ps, bunchSpacing);
0428       PUWorker_->putPileupInfo(e);
0429     }
0430 
0431     // Ecal
0432     if (MergeEMDigis_) {
0433       EMDigiWorker_->putEM(e, ES);
0434     } else {
0435       EMWorker_->putEM(e);
0436     }
0437 
0438     // Hcal
0439     if (MergeHcalDigis_) {
0440       if (MergeHcalDigisProd_) {
0441         HcalDigiWorkerProd_->putHcal(e, ES);
0442       } else {
0443         HcalDigiWorker_->putHcal(e, ES);
0444       }
0445     } else {
0446       HcalWorker_->putHcal(e);
0447     }
0448 
0449     // Muon
0450     MuonWorker_->putMuon(e);
0451 
0452     // SiStrips
0453     if (useSiStripRawDigi_)
0454       SiStripRawWorker_->putSiStrip(e);
0455     else
0456       SiStripWorker_->putSiStrip(e);
0457 
0458     // SiPixels
0459     SiPixelWorker_->putSiPixel(e);
0460   }
0461 
0462   void DataMixingModule::beginLuminosityBlock(LuminosityBlock const &l1, EventSetup const &c) {
0463     BMixingModule::beginLuminosityBlock(l1, c);
0464   }
0465 
0466   void DataMixingModule::endLuminosityBlock(LuminosityBlock const &l1, EventSetup const &c) {
0467     BMixingModule::endLuminosityBlock(l1, c);
0468   }
0469 
0470 }  // namespace edm