CMS 3D CMS Logo

Phase2TrackerDigitizer.h
Go to the documentation of this file.
1 #ifndef __SimTracker_SiPhase2Digitizer_Phase2TrackerDigitizer_h
2 #define __SimTracker_SiPhase2Digitizer_Phase2TrackerDigitizer_h
3 
4 //-------------------------------------------------------------
5 // class Phase2TrackerDigitizer
6 //
7 // Phase2TrackerDigitizer produces digis from SimHits
8 // The real algorithm is in Phase2TrackerDigitizerAlgorithm
9 //
10 // Author: Suchandra Dutta, Suvankar Roy Chowdhury, Subir Sarkar
11 //
12 //--------------------------------------------------------------
13 
14 #include <map>
15 #include <string>
16 #include <vector>
17 #include <unordered_map>
18 
29 
30 // Forward declaration
31 namespace CLHEP {
32  class HepRandomEngine;
33 }
34 
35 namespace edm {
36  class Event;
37  class EventSetup;
38  class ParameterSet;
39  template <typename T>
40  class Handle;
41  class ConsumesCollector;
42 } // namespace edm
43 
44 class MagneticField;
46 class PSimHit;
49 
50 namespace cms {
52  public:
53  using ModuleTypeCache = std::unordered_map<uint32_t, TrackerGeometry::ModuleType>;
54 
55  explicit Phase2TrackerDigitizer(const edm::ParameterSet& iConfig,
58  ~Phase2TrackerDigitizer() override;
59 
60  void initializeEvent(edm::Event const& e, edm::EventSetup const& c) override;
61  void accumulate(edm::Event const& e, edm::EventSetup const& c) override;
62  void accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& c, edm::StreamID const&) override;
63  void finalizeEvent(edm::Event& e, edm::EventSetup const& c) override;
64  virtual void beginJob() {}
65 
66  template <class T>
67  void accumulate_local(T const& iEvent, edm::EventSetup const& iSetup);
68 
69  // For premixing
70  void loadAccumulator(const std::map<uint32_t, std::map<int, float> >& accumulator);
71 
72  private:
73  using vstring = std::vector<std::string>;
74 
75  // constants of different algorithm types
77  AlgorithmType getAlgoType(uint32_t idet);
78 
79  void accumulatePixelHits(edm::Handle<std::vector<PSimHit> >, size_t globalSimHitIndex, const uint32_t tofBin);
80  void addPixelCollection(edm::Event& iEvent, const edm::EventSetup& iSetup, const bool ot_analog);
81 
82  // Templated for premixing
83  template <typename DigiType>
85 
86  bool first_;
87 
96  std::map<std::string, size_t> crossingSimHitIndexOffset_;
97  std::map<AlgorithmType, std::unique_ptr<Phase2TrackerDigitizerAlgorithm> > algomap_;
103  const TrackerGeometry* pDD_ = nullptr;
104  const MagneticField* pSetup_ = nullptr;
105  std::map<uint32_t, const Phase2TrackerGeomDetUnit*> detectorUnits_;
106  const TrackerTopology* tTopo_ = nullptr;
109  const bool premixStage1_;
110  const bool makeDigiSimLinks_;
111  // cache for detector types
113  };
114 } // namespace cms
115 #endif
edm::ESWatcher< TrackerDigiGeometryRecord > theTkDigiGeomWatcher_
Phase2TrackerDigitizer(const edm::ParameterSet &iConfig, edm::ProducesCollector, edm::ConsumesCollector &iC)
const TrackerTopology * tTopo_
void accumulate(edm::Event const &e, edm::EventSetup const &c) override
void accumulatePixelHits(edm::Handle< std::vector< PSimHit > >, size_t globalSimHitIndex, const uint32_t tofBin)
void finalizeEvent(edm::Event &e, edm::EventSetup const &c) override
AlgorithmType getAlgoType(uint32_t idet)
void initializeEvent(edm::Event const &e, edm::EventSetup const &c) override
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
void accumulate_local(T const &iEvent, edm::EventSetup const &iSetup)
int iEvent
Definition: GenABIO.cc:224
std::unordered_map< uint32_t, TrackerGeometry::ModuleType > ModuleTypeCache
void loadAccumulator(const std::map< uint32_t, std::map< int, float > > &accumulator)
Namespace of DDCMS conversion namespace.
std::map< uint32_t, const Phase2TrackerGeomDetUnit * > detectorUnits_
void addOuterTrackerCollection(edm::Event &iEvent, const edm::EventSetup &iSetup)
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > pSetupToken_
std::map< AlgorithmType, std::unique_ptr< Phase2TrackerDigitizerAlgorithm > > algomap_
HLT enums.
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > pDDToken_
void addPixelCollection(edm::Event &iEvent, const edm::EventSetup &iSetup, const bool ot_analog)
long double T
std::vector< std::string > vstring
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. ...