Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-01-21 00:19:54

0001 //
0002 //
0003 // Package:    Phase2TrackerDigitizer
0004 // Class:      Phase2TrackerDigitizer
0005 //
0006 // *\class SiPhase2TrackerDigitizer Phase2TrackerDigitizer.cc SimTracker/SiPhase2Digitizer/src/Phase2TrackerDigitizer.cc
0007 //
0008 // Author: Suchandra Dutta, Suvankar Roy Chowdhury, Subir Sarkar
0009 // Date: January 29, 2016
0010 // Description: <one line class summary>
0011 //
0012 // Implementation:
0013 //     <Notes on implementation>
0014 //
0015 #include <memory>
0016 #include <set>
0017 #include <iostream>
0018 
0019 #include "SimTracker/SiPhase2Digitizer/plugins/Phase2TrackerDigitizer.h"
0020 #include "SimTracker/SiPhase2Digitizer/plugins/Phase2TrackerDigitizerAlgorithm.h"
0021 #include "SimTracker/SiPhase2Digitizer/plugins/SSDigitizerAlgorithm.h"
0022 #include "SimTracker/SiPhase2Digitizer/plugins/PSSDigitizerAlgorithm.h"
0023 #include "SimTracker/SiPhase2Digitizer/plugins/PSPDigitizerAlgorithm.h"
0024 #include "SimTracker/SiPhase2Digitizer/plugins/PixelDigitizerAlgorithm.h"
0025 #include "SimTracker/SiPhase2Digitizer/plugins/Pixel3DDigitizerAlgorithm.h"
0026 #include "SimTracker/Common/interface/DigitizerUtility.h"
0027 
0028 #include "FWCore/Framework/interface/EventSetup.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/InputTag.h"
0031 #include "FWCore/ServiceRegistry/interface/Service.h"
0032 #include "FWCore/Framework/interface/ConsumesCollector.h"
0033 
0034 #include "FWCore/Framework/interface/Event.h"
0035 #include "FWCore/Framework/interface/LuminosityBlock.h"
0036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0037 
0038 #include "DataFormats/Common/interface/Handle.h"
0039 #include "DataFormats/SiPixelDigi/interface/PixelDigi.h"
0040 #include "DataFormats/SiPixelDigi/interface/PixelDigiCollection.h"
0041 #include "DataFormats/Phase2TrackerDigi/interface/Phase2TrackerDigi.h"
0042 #include "DataFormats/Common/interface/DetSet.h"
0043 #include "DataFormats/Common/interface/DetSetVector.h"
0044 
0045 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0046 #include "SimDataFormats/TrackerDigiSimLink/interface/PixelDigiSimLink.h"
0047 
0048 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0049 #include "Geometry/CommonDetUnit/interface/PixelGeomDetUnit.h"
0050 
0051 #include "MagneticField/Engine/interface/MagneticField.h"
0052 
0053 #include "SimGeneral/MixingModule/interface/PileUpEventPrincipal.h"
0054 
0055 // Random Number
0056 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0057 #include "FWCore/Utilities/interface/Exception.h"
0058 
0059 namespace cms {
0060 
0061   Phase2TrackerDigitizer::Phase2TrackerDigitizer(const edm::ParameterSet& iConfig,
0062                                                  edm::ProducesCollector producesCollector,
0063                                                  edm::ConsumesCollector& iC)
0064       : first_(true),
0065         hitsProducer_(iConfig.getParameter<std::string>("hitsProducer")),
0066         trackerContainers_(iConfig.getParameter<std::vector<std::string> >("ROUList")),
0067         pDDToken_(iC.esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("GeometryType")))),
0068         pSetupToken_(iC.esConsumes()),
0069         tTopoToken_(iC.esConsumes()),
0070         isOuterTrackerReadoutAnalog_(iConfig.getParameter<bool>("isOTreadoutAnalog")),
0071         premixStage1_(iConfig.getParameter<bool>("premixStage1")),
0072         makeDigiSimLinks_(
0073             iConfig.getParameter<edm::ParameterSet>("AlgorithmCommon").getUntrackedParameter<bool>("makeDigiSimLinks")) {
0074     const std::string alias1("simSiPixelDigis");
0075     producesCollector.produces<edm::DetSetVector<PixelDigi> >("Pixel").setBranchAlias(alias1);
0076     if (makeDigiSimLinks_)
0077       producesCollector.produces<edm::DetSetVector<PixelDigiSimLink> >("Pixel").setBranchAlias(alias1);
0078 
0079     if (!iConfig.getParameter<bool>("isOTreadoutAnalog")) {
0080       const std::string alias2("simSiTrackerDigis");
0081       if (premixStage1_) {
0082         // Premixing exploits the ADC field of PixelDigi to store the collected charge
0083         // But we still want everything else to be treated like for Phase2TrackerDigi
0084         producesCollector.produces<edm::DetSetVector<PixelDigi> >("Tracker").setBranchAlias(alias2);
0085       } else {
0086         producesCollector.produces<edm::DetSetVector<Phase2TrackerDigi> >("Tracker").setBranchAlias(alias2);
0087       }
0088       if (makeDigiSimLinks_)
0089         producesCollector.produces<edm::DetSetVector<PixelDigiSimLink> >("Tracker").setBranchAlias(alias2);
0090     }
0091     // creating algorithm objects and pushing them into the map
0092     algomap_[AlgorithmType::InnerPixel] = std::make_unique<PixelDigitizerAlgorithm>(iConfig, iC);
0093     algomap_[AlgorithmType::InnerPixel3D] = std::make_unique<Pixel3DDigitizerAlgorithm>(iConfig, iC);
0094     algomap_[AlgorithmType::PixelinPS] = std::make_unique<PSPDigitizerAlgorithm>(iConfig, iC);
0095     algomap_[AlgorithmType::StripinPS] = std::make_unique<PSSDigitizerAlgorithm>(iConfig, iC);
0096     algomap_[AlgorithmType::TwoStrip] = std::make_unique<SSDigitizerAlgorithm>(iConfig, iC);
0097   }
0098 
0099   Phase2TrackerDigitizer::~Phase2TrackerDigitizer() {}
0100   void Phase2TrackerDigitizer::accumulatePixelHits(edm::Handle<std::vector<PSimHit> > hSimHits,
0101                                                    size_t globalSimHitIndex,
0102                                                    const uint32_t tofBin) {
0103     if (hSimHits.isValid()) {
0104       std::set<uint32_t> detIds;
0105       auto const& simHits = *(hSimHits.product());
0106       for (auto it = std::begin(simHits), itEnd = std::end(simHits); it != itEnd; ++it, ++globalSimHitIndex) {
0107         uint32_t detId_raw = (*it).detUnitId();
0108         auto fiter = detectorUnits_.find(detId_raw);
0109         if (fiter == detectorUnits_.end())
0110           continue;
0111 
0112         if (detIds.insert(detId_raw).second) {
0113           // The insert succeeded, so this detector element has not yet been processed.
0114           const Phase2TrackerGeomDetUnit* phase2det = fiter->second;
0115 
0116           // access to magnetic field in global coordinates
0117           GlobalVector bfield = pSetup_->inTesla(phase2det->surface().position());
0118           LogDebug("PixelDigitizer") << "B-field(T) at " << phase2det->surface().position()
0119                                      << " (cm): " << pSetup_->inTesla(phase2det->surface().position());
0120 
0121           auto kiter = algomap_.find(getAlgoType(detId_raw));
0122           if (kiter != algomap_.end())
0123             kiter->second->accumulateSimHits(it, itEnd, globalSimHitIndex, tofBin, phase2det, bfield);
0124           else
0125             edm::LogInfo("Phase2TrackerDigitizer") << "Unsupported algorithm: ";
0126         }
0127       }
0128     }
0129   }
0130 
0131   void Phase2TrackerDigitizer::initializeEvent(edm::Event const& e, edm::EventSetup const& iSetup) {
0132     edm::Service<edm::RandomNumberGenerator> rng;
0133     if (!rng.isAvailable()) {
0134       throw cms::Exception("Configuration")
0135           << "Phase2TrackerDigitizer requires the RandomNumberGeneratorService\n"
0136              "which is not present in the configuration file.  You must add the service\n"
0137              "in the configuration file or remove the modules that require it.";
0138     }
0139 
0140     pSetup_ = &iSetup.getData(pSetupToken_);
0141     tTopo_ = &iSetup.getData(tTopoToken_);
0142 
0143     if (theTkDigiGeomWatcher_.check(iSetup)) {
0144       pDD_ = &iSetup.getData(pDDToken_);
0145 
0146       // reset cache
0147       ModuleTypeCache().swap(moduleTypeCache_);
0148       detectorUnits_.clear();
0149       for (auto const& det_u : pDD_->detUnits()) {
0150         uint32_t rawId = det_u->geographicalId().rawId();
0151         if (DetId(rawId).det() == DetId::Detector::Tracker) {
0152           const Phase2TrackerGeomDetUnit* pixdet = dynamic_cast<const Phase2TrackerGeomDetUnit*>(det_u);
0153           assert(pixdet);
0154           detectorUnits_.emplace(rawId, pixdet);
0155         }
0156       }
0157     }
0158 
0159     // Must initialize all the algorithms
0160     for (auto const& el : algomap_) {
0161       if (first_)
0162         el.second->init(iSetup);
0163 
0164       el.second->initializeEvent(rng->getEngine(e.streamID()));
0165     }
0166     first_ = false;
0167     // Make sure that the first crossing processed starts indexing the sim hits from zero.
0168     // This variable is used so that the sim hits from all crossing frames have sequential
0169     // indices used to create the digi-sim link (if configured to do so) rather than starting
0170     // from zero for each crossing.
0171     crossingSimHitIndexOffset_.clear();
0172   }
0173   void Phase2TrackerDigitizer::accumulate(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
0174     accumulate_local<edm::Event>(iEvent, iSetup);
0175   }
0176 
0177   void Phase2TrackerDigitizer::accumulate(PileUpEventPrincipal const& iEvent,
0178                                           edm::EventSetup const& iSetup,
0179                                           edm::StreamID const&) {
0180     accumulate_local<PileUpEventPrincipal>(iEvent, iSetup);
0181   }
0182 
0183   template <class T>
0184   void Phase2TrackerDigitizer::accumulate_local(T const& iEvent, edm::EventSetup const& iSetup) {
0185     for (auto const& v : trackerContainers_) {
0186       edm::Handle<std::vector<PSimHit> > simHits;
0187       edm::InputTag tag(hitsProducer_, v);
0188       iEvent.getByLabel(tag, simHits);
0189 
0190       //edm::EDGetTokenT< std::vector<PSimHit> > simHitToken_(consumes< std::vector<PSimHit>(tag));
0191       //iEvent.getByToken(simHitToken_, simHits);
0192 
0193       uint32_t tofBin = PixelDigiSimLink::LowTof;
0194       if (v.find(std::string("HighTof")) != std::string::npos)
0195         tofBin = PixelDigiSimLink::HighTof;
0196       accumulatePixelHits(simHits, crossingSimHitIndexOffset_[tag.encode()], tofBin);
0197       // Now that the hits have been processed, I'll add the amount of hits in this crossing on to
0198       // the global counter. Next time accumulateStripHits() is called it will count the sim hits
0199       // as though they were on the end of this collection.
0200       // Note that this is only used for creating digi-sim links (if configured to do so).
0201       if (simHits.isValid())
0202         crossingSimHitIndexOffset_[tag.encode()] += simHits->size();
0203     }
0204   }
0205 
0206   // For premixing
0207   void Phase2TrackerDigitizer::loadAccumulator(const std::map<uint32_t, std::map<int, float> >& accumulator) {
0208     for (const auto& detMap : accumulator) {
0209       AlgorithmType algoType = getAlgoType(detMap.first);
0210       auto& algo = *(algomap_.at(algoType));
0211       algo.loadAccumulator(detMap.first, detMap.second);
0212     }
0213   }
0214 
0215   void Phase2TrackerDigitizer::finalizeEvent(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0216     //Decide if we want analog readout for Outer Tracker.
0217     addPixelCollection(iEvent, iSetup, isOuterTrackerReadoutAnalog_);
0218     if (!isOuterTrackerReadoutAnalog_) {
0219       if (premixStage1_)
0220         addOuterTrackerCollection<PixelDigi>(iEvent, iSetup);
0221       else
0222         addOuterTrackerCollection<Phase2TrackerDigi>(iEvent, iSetup);
0223     }
0224   }
0225   Phase2TrackerDigitizer::AlgorithmType Phase2TrackerDigitizer::getAlgoType(uint32_t detId_raw) {
0226     // get mType either from the geometry or from our cache (faster)
0227     TrackerGeometry::ModuleType mType = TrackerGeometry::ModuleType::UNKNOWN;
0228     auto itr = moduleTypeCache_.find(detId_raw);
0229     if (itr != moduleTypeCache_.end()) {
0230       mType = itr->second;
0231     } else {
0232       mType = pDD_->getDetectorType(DetId(detId_raw));
0233       moduleTypeCache_.emplace(detId_raw, mType);
0234     }
0235 
0236     AlgorithmType algotype = AlgorithmType::Unknown;
0237     switch (mType) {
0238       case TrackerGeometry::ModuleType::Ph1PXB:
0239         algotype = AlgorithmType::InnerPixel;
0240         break;
0241       case TrackerGeometry::ModuleType::Ph1PXF:
0242         algotype = AlgorithmType::InnerPixel;
0243         break;
0244       case TrackerGeometry::ModuleType::Ph2PXB:
0245         algotype = AlgorithmType::InnerPixel;
0246         break;
0247       case TrackerGeometry::ModuleType::Ph2PXF:
0248         algotype = AlgorithmType::InnerPixel;
0249         break;
0250       case TrackerGeometry::ModuleType::Ph2PXB3D:
0251         algotype = AlgorithmType::InnerPixel3D;
0252         break;
0253       case TrackerGeometry::ModuleType::Ph2PXF3D:
0254         algotype = AlgorithmType::InnerPixel3D;
0255         break;
0256       case TrackerGeometry::ModuleType::Ph2PSP:
0257         algotype = AlgorithmType::PixelinPS;
0258         break;
0259       case TrackerGeometry::ModuleType::Ph2PSS:
0260         algotype = AlgorithmType::StripinPS;
0261         break;
0262       case TrackerGeometry::ModuleType::Ph2SS:
0263         algotype = AlgorithmType::TwoStrip;
0264         break;
0265       default:
0266         edm::LogError("Phase2TrackerDigitizer") << "ERROR - Wrong Detector Type, No Algorithm available ";
0267     }
0268 
0269     return algotype;
0270   }
0271   void Phase2TrackerDigitizer::addPixelCollection(edm::Event& iEvent,
0272                                                   const edm::EventSetup& iSetup,
0273                                                   const bool ot_analog) {
0274     std::vector<edm::DetSet<PixelDigi> > digiVector;
0275     std::vector<edm::DetSet<PixelDigiSimLink> > digiLinkVector;
0276     for (auto const& det_u : pDD_->detUnits()) {
0277       uint32_t rawId = det_u->geographicalId().rawId();
0278       auto algotype = getAlgoType(rawId);
0279       auto fiter = algomap_.find(algotype);
0280       if (fiter == algomap_.end())
0281         continue;
0282 
0283       // Decide if we want analog readout for Outer Tracker.
0284       if (!ot_analog && algotype != AlgorithmType::InnerPixel && algotype != AlgorithmType::InnerPixel3D)
0285         continue;
0286 
0287       std::map<int, digitizerUtility::DigiSimInfo> digi_map;
0288       fiter->second->digitize(dynamic_cast<const Phase2TrackerGeomDetUnit*>(det_u), digi_map, tTopo_);
0289 
0290       edm::DetSet<PixelDigi> collector(rawId);
0291       edm::DetSet<PixelDigiSimLink> linkcollector(rawId);
0292       for (auto const& digi_p : digi_map) {
0293         digitizerUtility::DigiSimInfo info = digi_p.second;
0294         const auto& ip = PixelDigi::channelToPixel(digi_p.first);
0295         collector.data.emplace_back(ip.first, ip.second, info.sig_tot);
0296         for (auto const& sim_p : info.simInfoList) {
0297           linkcollector.data.emplace_back(digi_p.first,
0298                                           sim_p.second->trackId(),
0299                                           sim_p.second->hitIndex(),
0300                                           sim_p.second->tofBin(),
0301                                           sim_p.second->eventId(),
0302                                           sim_p.first);
0303         }
0304       }
0305       if (!collector.data.empty())
0306         digiVector.push_back(std::move(collector));
0307       if (!linkcollector.data.empty())
0308         digiLinkVector.push_back(std::move(linkcollector));
0309     }
0310 
0311     // Step C: create collection with the cache vector of DetSet
0312     auto output = std::make_unique<edm::DetSetVector<PixelDigi> >(digiVector);
0313     auto outputlink = std::make_unique<edm::DetSetVector<PixelDigiSimLink> >(digiLinkVector);
0314 
0315     // Step D: write output to file
0316     iEvent.put(std::move(output), "Pixel");
0317     if (makeDigiSimLinks_)
0318       iEvent.put(std::move(outputlink), "Pixel");
0319   }
0320 }  // namespace cms
0321 namespace {
0322   void addToCollector(edm::DetSet<PixelDigi>& collector, const int channel, const digitizerUtility::DigiSimInfo& info) {
0323     // For premixing stage1 the channel must be decoded with PixelDigi
0324     // so that when the row and column are inserted to PixelDigi the
0325     // coded channel stays the same (so that it can then be decoded
0326     // with Phase2TrackerDigi in stage2).
0327     const auto& ip = PixelDigi::channelToPixel(channel);
0328     collector.data.emplace_back(ip.first, ip.second, info.sig_tot);
0329   }
0330   void addToCollector(edm::DetSet<Phase2TrackerDigi>& collector,
0331                       const int channel,
0332                       const digitizerUtility::DigiSimInfo& info) {
0333     const auto& ip = Phase2TrackerDigi::channelToPixel(channel);
0334     collector.data.emplace_back(ip.first, ip.second, info.ot_bit);
0335   }
0336 }  // namespace
0337 namespace cms {
0338   template <typename DigiType>
0339   void Phase2TrackerDigitizer::addOuterTrackerCollection(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0340     std::vector<edm::DetSet<DigiType> > digiVector;
0341     std::vector<edm::DetSet<PixelDigiSimLink> > digiLinkVector;
0342     for (auto const& det_u : pDD_->detUnits()) {
0343       uint32_t rawId = det_u->geographicalId().rawId();
0344       auto algotype = getAlgoType(rawId);
0345 
0346       auto fiter = algomap_.find(algotype);
0347       if (fiter == algomap_.end() || algotype == AlgorithmType::InnerPixel || algotype == AlgorithmType::InnerPixel3D)
0348         continue;
0349 
0350       std::map<int, digitizerUtility::DigiSimInfo> digi_map;
0351       fiter->second->digitize(dynamic_cast<const Phase2TrackerGeomDetUnit*>(det_u), digi_map, tTopo_);
0352 
0353       edm::DetSet<DigiType> collector(rawId);
0354       edm::DetSet<PixelDigiSimLink> linkcollector(rawId);
0355       for (auto const& digi_p : digi_map) {
0356         digitizerUtility::DigiSimInfo info = digi_p.second;
0357         addToCollector(collector, digi_p.first, info);
0358         for (auto const& sim_p : info.simInfoList) {
0359           linkcollector.data.emplace_back(digi_p.first,
0360                                           sim_p.second->trackId(),
0361                                           sim_p.second->hitIndex(),
0362                                           sim_p.second->tofBin(),
0363                                           sim_p.second->eventId(),
0364                                           sim_p.first);
0365         }
0366       }
0367 
0368       if (!collector.data.empty())
0369         digiVector.push_back(std::move(collector));
0370       if (!linkcollector.data.empty())
0371         digiLinkVector.push_back(std::move(linkcollector));
0372     }
0373 
0374     // Step C: create collection with the cache vector of DetSet
0375     auto output = std::make_unique<edm::DetSetVector<DigiType> >(digiVector);
0376     auto outputlink = std::make_unique<edm::DetSetVector<PixelDigiSimLink> >(digiLinkVector);
0377 
0378     // Step D: write output to file
0379     iEvent.put(std::move(output), "Tracker");
0380     if (makeDigiSimLinks_)
0381       iEvent.put(std::move(outputlink), "Tracker");
0382   }
0383 }  // namespace cms
0384 
0385 #include "FWCore/Framework/interface/MakerMacros.h"
0386 #include "SimGeneral/MixingModule/interface/DigiAccumulatorMixModFactory.h"
0387 
0388 using cms::Phase2TrackerDigitizer;
0389 DEFINE_DIGI_ACCUMULATOR(Phase2TrackerDigitizer);