CMS 3D CMS Logo

HGCDigitizerBase.h
Go to the documentation of this file.
1 #ifndef SimCalorimetry_HGCSimProducers_hgcdigitizerbase
2 #define SimCalorimetry_HGCSimProducers_hgcdigitizerbase
3 
4 #include <array>
5 #include <iostream>
6 #include <vector>
7 #include <memory>
8 #include <unordered_map>
9 #include <unordered_set>
10 
15 
17 
19 
21 
22 namespace hgc = hgc_digi;
23 
24 template <class DFr>
26  public:
27 
28  typedef DFr DigiType;
29 
31 
39  void run(std::unique_ptr<DColl> &digiColl, hgc::HGCSimHitDataAccumulator &simData,
40  const CaloSubdetectorGeometry* theGeom, const std::unordered_set<DetId>& validIds,
41  uint32_t digitizationType,CLHEP::HepRandomEngine* engine);
42 
46  float keV2fC() const { return keV2fC_; }
48  float tdcOnset() const { return myFEelectronics_->getTDCOnset(); }
49  std::array<float,3> tdcForToAOnset() const { return myFEelectronics_->getTDCForToAOnset(); }
50 
54  void runSimple(std::unique_ptr<DColl> &coll, hgc::HGCSimHitDataAccumulator &simData,
55  const CaloSubdetectorGeometry* theGeom, const std::unordered_set<DetId>& validIds,
56  CLHEP::HepRandomEngine* engine);
57 
61  void updateOutput(std::unique_ptr<DColl> &coll, const DFr& rawDataFrame);
62 
66  virtual void runDigitizer(std::unique_ptr<DColl> &coll, hgc::HGCSimHitDataAccumulator &simData,
67  const CaloSubdetectorGeometry* theGeom, const std::unordered_set<DetId>& validIds,
68  uint32_t digitizerType, CLHEP::HepRandomEngine* engine)
69  {
70  throw cms::Exception("HGCDigitizerBaseException") << " Failed to find specialization of runDigitizer";
71  }
72 
76  virtual ~HGCDigitizerBase()
77  { };
78 
79 
80 
81  protected:
82 
83  //baseline configuration
85 
86  //1keV in fC
87  float keV2fC_;
88 
89  //noise level
90  std::vector<float> noise_fC_;
91 
92  //charge collection efficiency
93  std::vector<double> cce_;
94 
95  //front-end electronics model
96  std::unique_ptr<HGCFEElectronics<DFr> > myFEelectronics_;
97 
98  //bunch time
99  double bxTime_;
100 
101  //if true will put both in time and out-of-time samples in the event
103 
104 };
105 
106 #endif
std::vector< float > noise_fC_
void runSimple(std::unique_ptr< DColl > &coll, hgc::HGCSimHitDataAccumulator &simData, const CaloSubdetectorGeometry *theGeom, const std::unordered_set< DetId > &validIds, CLHEP::HepRandomEngine *engine)
a trivial digitization: sum energies and digitize without noise
void updateOutput(std::unique_ptr< DColl > &coll, const DFr &rawDataFrame)
prepares the output according to the number of time samples to produce
std::array< float, 3 > tdcForToAOnset() const
virtual ~HGCDigitizerBase()
DTOR.
edm::ParameterSet myCfg_
std::unordered_map< uint32_t, HGCCellInfo > HGCSimHitDataAccumulator
edm::SortedCollection< DFr > DColl
float tdcOnset() const
void run(std::unique_ptr< DColl > &digiColl, hgc::HGCSimHitDataAccumulator &simData, const CaloSubdetectorGeometry *theGeom, const std::unordered_set< DetId > &validIds, uint32_t digitizationType, CLHEP::HepRandomEngine *engine)
steer digitization mode
std::unique_ptr< HGCFEElectronics< DFr > > myFEelectronics_
HGCDigitizerBase(const edm::ParameterSet &ps)
CTOR.
float keV2fC() const
getters
std::vector< double > cce_
JetCorrectorParametersCollection coll
Definition: classes.h:10
virtual void runDigitizer(std::unique_ptr< DColl > &coll, hgc::HGCSimHitDataAccumulator &simData, const CaloSubdetectorGeometry *theGeom, const std::unordered_set< DetId > &validIds, uint32_t digitizerType, CLHEP::HepRandomEngine *engine)
to be specialized by top class
models the behavior of the front-end electronics
bool toaModeByEnergy() const