Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef __SimTracker_SiPhase2Digitizer_Phase2TrackerDigitizer_h
0002 #define __SimTracker_SiPhase2Digitizer_Phase2TrackerDigitizer_h
0003 
0004 //-------------------------------------------------------------
0005 // class Phase2TrackerDigitizer
0006 //
0007 // Phase2TrackerDigitizer produces digis from SimHits
0008 // The real algorithm is in Phase2TrackerDigitizerAlgorithm
0009 //
0010 // Author: Suchandra Dutta, Suvankar Roy Chowdhury, Subir Sarkar
0011 //
0012 //--------------------------------------------------------------
0013 
0014 #include <map>
0015 #include <string>
0016 #include <vector>
0017 #include <unordered_map>
0018 
0019 #include "FWCore/Framework/interface/ESWatcher.h"
0020 #include "FWCore/Framework/interface/ProducesCollector.h"
0021 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0022 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0023 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0024 #include "DataFormats/DetId/interface/DetId.h"
0025 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0026 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0027 #include "SimGeneral/MixingModule/interface/DigiAccumulatorMixMod.h"
0028 #include "SimTracker/SiPhase2Digitizer/plugins/Phase2TrackerDigitizerFwd.h"
0029 
0030 // Forward declaration
0031 namespace CLHEP {
0032   class HepRandomEngine;
0033 }
0034 
0035 namespace edm {
0036   class Event;
0037   class EventSetup;
0038   class ParameterSet;
0039   template <typename T>
0040   class Handle;
0041   class ConsumesCollector;
0042 }  // namespace edm
0043 
0044 class MagneticField;
0045 class PileUpEventPrincipal;
0046 class PSimHit;
0047 class Phase2TrackerDigitizerAlgorithm;
0048 class TrackerDigiGeometryRecord;
0049 
0050 namespace cms {
0051   class Phase2TrackerDigitizer : public DigiAccumulatorMixMod {
0052   public:
0053     using ModuleTypeCache = std::unordered_map<uint32_t, TrackerGeometry::ModuleType>;
0054 
0055     explicit Phase2TrackerDigitizer(const edm::ParameterSet& iConfig,
0056                                     edm::ProducesCollector,
0057                                     edm::ConsumesCollector& iC);
0058     ~Phase2TrackerDigitizer() override;
0059 
0060     void initializeEvent(edm::Event const& e, edm::EventSetup const& c) override;
0061     void accumulate(edm::Event const& e, edm::EventSetup const& c) override;
0062     void accumulate(PileUpEventPrincipal const& e, edm::EventSetup const& c, edm::StreamID const&) override;
0063     void finalizeEvent(edm::Event& e, edm::EventSetup const& c) override;
0064     virtual void beginJob() {}
0065 
0066     template <class T>
0067     void accumulate_local(T const& iEvent, edm::EventSetup const& iSetup);
0068 
0069     // For premixing
0070     void loadAccumulator(const std::map<uint32_t, std::map<int, float> >& accumulator);
0071 
0072   private:
0073     using vstring = std::vector<std::string>;
0074 
0075     // constants of different algorithm types
0076     enum class AlgorithmType { InnerPixel, InnerPixel3D, PixelinPS, StripinPS, TwoStrip, Unknown };
0077     AlgorithmType getAlgoType(uint32_t idet);
0078 
0079     void accumulatePixelHits(edm::Handle<std::vector<PSimHit> >, size_t globalSimHitIndex, const uint32_t tofBin);
0080     void addPixelCollection(edm::Event& iEvent, const edm::EventSetup& iSetup, const bool ot_analog);
0081 
0082     // Templated for premixing
0083     template <typename DigiType>
0084     void addOuterTrackerCollection(edm::Event& iEvent, const edm::EventSetup& iSetup);
0085 
0086     bool first_;
0087 
0088     /** @brief Offset to add to the index of each sim hit to account for which crossing it's in.
0089      *
0090      * I need to know what each sim hit index will be when the hits from all crossing frames are merged into
0091      * one collection (assuming the MixingModule is configured to create the crossing frame for all sim hits).
0092      * To do this I'll record how many hits were in each crossing, and then add that on to the index for a given
0093      * hit in a given crossing. This assumes that the crossings are processed in the same order here as they are
0094      * put into the crossing frame, which I'm pretty sure is true.<br/>
0095      * The key is the name of the sim hit collection. */
0096     std::map<std::string, size_t> crossingSimHitIndexOffset_;
0097     std::map<AlgorithmType, std::unique_ptr<Phase2TrackerDigitizerAlgorithm> > algomap_;
0098     const std::string hitsProducer_;
0099     const vstring trackerContainers_;
0100     const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> pDDToken_;
0101     const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> pSetupToken_;
0102     const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> tTopoToken_;
0103     const TrackerGeometry* pDD_ = nullptr;
0104     const MagneticField* pSetup_ = nullptr;
0105     std::map<uint32_t, const Phase2TrackerGeomDetUnit*> detectorUnits_;
0106     const TrackerTopology* tTopo_ = nullptr;
0107     edm::ESWatcher<TrackerDigiGeometryRecord> theTkDigiGeomWatcher_;
0108     const bool isOuterTrackerReadoutAnalog_;
0109     const bool usePseudoPixel3DAlgo_;
0110     const bool premixStage1_;
0111     const bool makeDigiSimLinks_;
0112     // cache for detector types
0113     ModuleTypeCache moduleTypeCache_;
0114   };
0115 }  // namespace cms
0116 #endif