Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-05 23:51:50

0001 #ifndef SiPixelDigitizer_h
0002 #define SiPixelDigitizer_h
0003 
0004 /** \class SiPixelDigitizer
0005  *
0006  * SiPixelDigitizer produces digis from SimHits
0007  * The real algorithm is in SiPixelDigitizerAlgorithm
0008  *
0009  * \author Michele Pioppi-INFN Perugia
0010  *
0011  * \version   Sep 26 2005  
0012 
0013  *
0014  ************************************************************/
0015 
0016 #include <map>
0017 #include <memory>
0018 #include <string>
0019 #include <vector>
0020 
0021 #include "DataFormats/Common/interface/Handle.h"
0022 #include "DataFormats/Provenance/interface/EventID.h"
0023 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0024 #include "FWCore/Framework/interface/ProducesCollector.h"
0025 #include "FWCore/Framework/interface/FrameworkfwdMostUsed.h"
0026 #include "FWCore/Utilities/interface/ESGetToken.h"
0027 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0028 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0029 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0030 #include "SimGeneral/MixingModule/interface/DigiAccumulatorMixMod.h"
0031 
0032 class MagneticField;
0033 class PileUpEventPrincipal;
0034 class PixelGeomDetUnit;
0035 class PSimHit;
0036 class SiPixelDigitizerAlgorithm;
0037 class TrackerGeometry;
0038 
0039 namespace CLHEP {
0040   class HepRandomEngine;
0041 }
0042 
0043 namespace cms {
0044   class SiPixelDigitizer : public DigiAccumulatorMixMod {
0045   public:
0046     explicit SiPixelDigitizer(const edm::ParameterSet& conf, edm::ProducesCollector, edm::ConsumesCollector& iC);
0047 
0048     ~SiPixelDigitizer() override;
0049 
0050     void initializeEvent(edm::Event const& e, edm::EventSetup const& c) override;
0051     void accumulate(edm::Event const& e, edm::EventSetup const& c) override;
0052     void accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& c, edm::StreamID const&) override;
0053     void finalizeEvent(edm::Event& e, edm::EventSetup const& c) override;
0054 
0055     virtual void beginJob() {}
0056 
0057     void StorePileupInformation(std::vector<int>& numInteractionList,
0058                                 std::vector<int>& bunchCrossingList,
0059                                 std::vector<float>& TrueInteractionList,
0060                                 std::vector<edm::EventID>& eventInfoList,
0061                                 int bunchSpacing) override {
0062       PileupInfo_ = std::make_unique<PileupMixingContent>(
0063           numInteractionList, bunchCrossingList, TrueInteractionList, eventInfoList, bunchSpacing);
0064     }
0065 
0066     PileupMixingContent* getEventPileupInfo() override { return PileupInfo_.get(); }
0067 
0068   private:
0069     void accumulatePixelHits(edm::Handle<std::vector<PSimHit> >,
0070                              size_t globalSimHitIndex,
0071                              const unsigned int tofBin,
0072                              edm::EventSetup const& c);
0073 
0074     bool firstInitializeEvent_;
0075     bool firstFinalizeEvent_;
0076     bool applyLateReweighting_;
0077     bool usePixelExtraLiteFormat_;
0078     const bool store_SimHitEntryExitPoints_;
0079     const bool store_SimHitEntryExitPointsLite_;
0080     bool makeDigiSimLinks_;
0081     std::unique_ptr<SiPixelDigitizerAlgorithm> _pixeldigialgo;
0082     /** @brief Offset to add to the index of each sim hit to account for which crossing it's in.
0083 *
0084 * I need to know what each sim hit index will be when the hits from all crossing frames are merged into
0085 * one collection (assuming the MixingModule is configured to create the crossing frame for all sim hits).
0086 * To do this I'll record how many hits were in each crossing, and then add that on to the index for a given
0087 * hit in a given crossing. This assumes that the crossings are processed in the same order here as they are
0088 * put into the crossing frame, which I'm pretty sure is true.<br/>
0089 * The key is the name of the sim hit collection. */
0090     std::map<std::string, size_t> crossingSimHitIndexOffset_;
0091 
0092     typedef std::vector<std::string> vstring;
0093     const std::string hitsProducer;
0094     const vstring trackerContainers;
0095     const TrackerGeometry* pDD = nullptr;
0096     const MagneticField* pSetup = nullptr;
0097     std::map<unsigned int, PixelGeomDetUnit const*> detectorUnits;
0098     CLHEP::HepRandomEngine* randomEngine_ = nullptr;
0099 
0100     std::unique_ptr<PileupMixingContent> PileupInfo_;
0101 
0102     const bool pilotBlades;         // Default = false
0103     const int NumberOfEndcapDisks;  // Default = 2
0104 
0105     const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoToken_;
0106     const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> pDDToken_;
0107     const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> pSetupToken_;
0108 
0109     // infrastructure to reject dead pixels as defined in db (added by F.Blekman)
0110   };
0111 }  // namespace cms
0112 
0113 #endif