Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "SimMuon/GEMDigitizer/interface/GEMDigiModule.h"
0002 
0003 #include "SimMuon/GEMDigitizer/interface/GEMSignalModel.h"
0004 #include "SimMuon/GEMDigitizer/interface/GEMBkgModel.h"
0005 #include "SimMuon/GEMDigitizer/interface/GEMNoiseModel.h"
0006 #include "Geometry/GEMGeometry/interface/GEMEtaPartition.h"
0007 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
0008 
0009 GEMDigiModule::GEMDigiModule(const edm::ParameterSet& config) {
0010   bool simulateBkgNoise_(config.getParameter<bool>("simulateBkgNoise"));
0011   bool simulateIntrinsicNoise_(config.getParameter<bool>("simulateIntrinsicNoise"));
0012   if (simulateIntrinsicNoise_) {
0013     models.push_back(std::make_unique<GEMNoiseModel>(config));
0014   }
0015   if (simulateBkgNoise_) {
0016     models.push_back(std::make_unique<GEMBkgModel>(config));
0017   }
0018   models.push_back(std::make_unique<GEMSignalModel>(config));
0019 }
0020 
0021 GEMDigiModule::~GEMDigiModule() = default;
0022 
0023 void GEMDigiModule::simulate(const GEMEtaPartition* roll,
0024                              const edm::PSimHitContainer& simHits,
0025                              CLHEP::HepRandomEngine* engine) {
0026   stripDigiSimLinks_.clear();
0027   detectorHitMap_.clear();
0028   stripDigiSimLinks_ = StripDigiSimLinks(roll->id().rawId());
0029   theGemDigiSimLinks_.clear();
0030   theGemDigiSimLinks_ = GEMDigiSimLinks(roll->id().rawId());
0031   for (auto&& model : models) {
0032     model->simulate(roll, simHits, engine, strips_, detectorHitMap_);
0033   }
0034   return;
0035 }
0036 
0037 void GEMDigiModule::fillDigis(int rollDetId, GEMDigiCollection& digis) {
0038   for (const auto& d : strips_) {
0039     if (d.second == -999)
0040       continue;
0041     // (strip, bx)
0042     GEMDigi digi(d.first, d.second);
0043     digis.insertDigi(GEMDetId(rollDetId), digi);
0044     addLinks(d.first, d.second);
0045     addLinksWithPartId(d.first, d.second);
0046   }
0047   strips_.clear();
0048 }
0049 
0050 void GEMDigiModule::addLinks(unsigned int strip, int bx) {
0051   std::pair<unsigned int, int> digi(strip, bx);
0052   auto channelHitItr = detectorHitMap_.equal_range(digi);
0053 
0054   // find the fraction contribution for each SimTrack
0055   std::map<int, float> simTrackChargeMap;
0056   std::map<int, EncodedEventId> eventIdMap;
0057   float totalCharge(0.);
0058   for (auto hitItr = channelHitItr.first; hitItr != channelHitItr.second; ++hitItr) {
0059     const PSimHit* hit(hitItr->second);
0060     // might be zero for unit tests and such
0061     if (hit == nullptr)
0062       continue;
0063 
0064     int simTrackId(hit->trackId());
0065     //float charge = hit->getCharge();
0066     const float charge(1.f);
0067     auto chargeItr = simTrackChargeMap.find(simTrackId);
0068     if (chargeItr == simTrackChargeMap.end()) {
0069       simTrackChargeMap[simTrackId] = charge;
0070       eventIdMap[simTrackId] = hit->eventId();
0071     } else {
0072       chargeItr->second += charge;
0073     }
0074     totalCharge += charge;
0075   }
0076 
0077   for (const auto& charge : simTrackChargeMap) {
0078     const int simTrackId(charge.first);
0079     auto link(StripDigiSimLink(strip, simTrackId, eventIdMap[simTrackId], charge.second / totalCharge));
0080     stripDigiSimLinks_.push_back(link);
0081   }
0082 }
0083 
0084 void GEMDigiModule::addLinksWithPartId(unsigned int strip, int bx) {
0085   std::pair<unsigned int, int> digi(strip, bx);
0086   std::pair<DetectorHitMap::iterator, DetectorHitMap::iterator> channelHitItr = detectorHitMap_.equal_range(digi);
0087 
0088   for (DetectorHitMap::iterator hitItr = channelHitItr.first; hitItr != channelHitItr.second; ++hitItr) {
0089     const PSimHit* hit = (hitItr->second);
0090     // might be zero for unit tests and such
0091     if (hit == nullptr)
0092       continue;
0093 
0094     theGemDigiSimLinks_.push_back(GEMDigiSimLink(strip, bx, hit->particleType(), hit->trackId(), hit->eventId()));
0095   }
0096 }
0097 
0098 void GEMDigiModule::setGeometry(const GEMGeometry* geom) {
0099   for (auto&& model : models) {
0100     model->setGeometry(geom);
0101   }
0102 }