CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelDigitizer.h
Go to the documentation of this file.
1 #ifndef SiPixelDigitizer_h
2 #define SiPixelDigitizer_h
3 
16 #include <map>
17 #include <memory>
18 #include <string>
19 #include <vector>
20 
31 
32 class MagneticField;
34 class PixelGeomDetUnit;
35 class PSimHit;
37 class TrackerGeometry;
38 
39 namespace CLHEP {
40  class HepRandomEngine;
41 }
42 
43 namespace cms {
45  public:
47 
48  ~SiPixelDigitizer() override;
49 
50  void initializeEvent(edm::Event const& e, edm::EventSetup const& c) override;
51  void accumulate(edm::Event const& e, edm::EventSetup const& c) override;
52  void accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& c, edm::StreamID const&) override;
53  void finalizeEvent(edm::Event& e, edm::EventSetup const& c) override;
54 
55  virtual void beginJob() {}
56 
57  void StorePileupInformation(std::vector<int>& numInteractionList,
58  std::vector<int>& bunchCrossingList,
59  std::vector<float>& TrueInteractionList,
60  std::vector<edm::EventID>& eventInfoList,
61  int bunchSpacing) override {
62  PileupInfo_ = std::make_unique<PileupMixingContent>(
63  numInteractionList, bunchCrossingList, TrueInteractionList, eventInfoList, bunchSpacing);
64  }
65 
66  PileupMixingContent* getEventPileupInfo() override { return PileupInfo_.get(); }
67 
68  private:
69  void accumulatePixelHits(edm::Handle<std::vector<PSimHit> >,
70  size_t globalSimHitIndex,
71  const unsigned int tofBin,
72  edm::EventSetup const& c);
73 
76  std::unique_ptr<SiPixelDigitizerAlgorithm> _pixeldigialgo;
85  std::map<std::string, size_t> crossingSimHitIndexOffset_;
86 
87  typedef std::vector<std::string> vstring;
90  const TrackerGeometry* pDD = nullptr;
91  const MagneticField* pSetup = nullptr;
92  std::map<unsigned int, PixelGeomDetUnit const*> detectorUnits;
93  CLHEP::HepRandomEngine* randomEngine_ = nullptr;
94 
95  std::unique_ptr<PileupMixingContent> PileupInfo_;
96 
97  const bool pilotBlades; // Default = false
98  const int NumberOfEndcapDisks; // Default = 2
99 
103 
104  // infrastructure to reject dead pixels as defined in db (added by F.Blekman)
105  };
106 } // namespace cms
107 
108 #endif
void finalizeEvent(edm::Event &e, edm::EventSetup const &c) override
virtual void beginJob()
const edm::EventSetup & c
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
const TrackerGeometry * pDD
const MagneticField * pSetup
std::vector< std::string > vstring
std::unique_ptr< PileupMixingContent > PileupInfo_
const std::string hitsProducer
SiPixelDigitizer(const edm::ParameterSet &conf, edm::ProducesCollector, edm::ConsumesCollector &iC)
std::map< unsigned int, PixelGeomDetUnit const * > detectorUnits
void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
void accumulatePixelHits(edm::Handle< std::vector< PSimHit > >, size_t globalSimHitIndex, const unsigned int tofBin, edm::EventSetup const &c)
std::unique_ptr< SiPixelDigitizerAlgorithm > _pixeldigialgo
void accumulate(edm::Event const &e, edm::EventSetup const &c) override
CLHEP::HepRandomEngine * randomEngine_
void StorePileupInformation(std::vector< int > &numInteractionList, std::vector< int > &bunchCrossingList, std::vector< float > &TrueInteractionList, std::vector< edm::EventID > &eventInfoList, int bunchSpacing) override
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > pSetupToken_
const vstring trackerContainers
PileupMixingContent * getEventPileupInfo() override
std::map< std::string, size_t > crossingSimHitIndexOffset_
Offset to add to the index of each sim hit to account for which crossing it&#39;s in. ...