CMS 3D CMS Logo

HGCDigitizerBase.cc
Go to the documentation of this file.
4 
5 using namespace hgc_digi;
6 
7 namespace {
8  void addCellMetadata(HGCCellInfo& info,
9  const HcalGeometry* geom,
10  const DetId& detid ) {
11  //base time samples for each DetId, initialized to 0
12  info.size = 1.0;
13  info.thickness = 1.0;
14  }
15 
16  void addCellMetadata(HGCCellInfo& info,
17  const HGCalGeometry* geom,
18  const DetId& detid ) {
19  const auto& topo = geom->topology();
20  const auto& dddConst = topo.dddConstants();
21  uint32_t id(detid.rawId());
22  int waferTypeL = 0;
23  bool isHalf = false;
24  HGCalDetId hid(id);
25  int wafer = HGCalDetId(id).wafer();
26  waferTypeL = dddConst.waferTypeL(wafer);
27  isHalf = dddConst.isHalfCell(wafer,hid.cell());
28  //base time samples for each DetId, initialized to 0
29  info.size = (isHalf ? 0.5 : 1.0);
30  info.thickness = waferTypeL;
31  }
32 
33  void addCellMetadata(HGCCellInfo& info,
34  const CaloSubdetectorGeometry* geom,
35  const DetId& detid ) {
36  if( DetId::Hcal == detid.det() ) {
37  const HcalGeometry* hc = static_cast<const HcalGeometry*>(geom);
38  addCellMetadata(info,hc,detid);
39  } else {
40  const HGCalGeometry* hg = static_cast<const HGCalGeometry*>(geom);
41  addCellMetadata(info,hg,detid);
42  }
43  }
44 
45 }
46 
47 
48 template<class DFr>
50  bxTime_ = ps.getParameter<double>("bxTime");
51  myCfg_ = ps.getParameter<edm::ParameterSet>("digiCfg");
52  doTimeSamples_ = myCfg_.getParameter< bool >("doTimeSamples");
53  if(myCfg_.exists("keV2fC")) keV2fC_ = myCfg_.getParameter<double>("keV2fC");
54  else keV2fC_ = 1.0;
55 
56  if( myCfg_.existsAs<std::vector<double> >( "chargeCollectionEfficiencies" ) ) {
57  cce_ = myCfg_.getParameter<std::vector<double> >("chargeCollectionEfficiencies");
58  } else {
59  std::vector<double>().swap(cce_);
60  }
61 
62  if(myCfg_.existsAs<double>("noise_fC")) {
63  noise_fC_.resize(1);
64  noise_fC_[0] = myCfg_.getParameter<double>("noise_fC");
65  } else if ( myCfg_.existsAs<std::vector<double> >("noise_fC") ) {
66  const auto& noises = myCfg_.getParameter<std::vector<double> >("noise_fC");
67  noise_fC_.resize(0);
68  noise_fC_.reserve(noises.size());
69  for( auto noise : noises ) { noise_fC_.push_back( noise ); }
70  } else {
71  noise_fC_.resize(1);
72  noise_fC_[0] = 1.f;
73  }
74  edm::ParameterSet feCfg = myCfg_.getParameter<edm::ParameterSet>("feCfg");
75  myFEelectronics_ = std::unique_ptr<HGCFEElectronics<DFr> >( new HGCFEElectronics<DFr>(feCfg) );
76 }
77 
78 template<class DFr>
79 void HGCDigitizerBase<DFr>::run( std::unique_ptr<HGCDigitizerBase::DColl> &digiColl,
80  HGCSimHitDataAccumulator &simData,
81  const CaloSubdetectorGeometry* theGeom,
82  const std::unordered_set<DetId>& validIds,
83  uint32_t digitizationType,
84  CLHEP::HepRandomEngine* engine) {
85  if(digitizationType==0) runSimple(digiColl,simData,theGeom,validIds,engine);
86  else runDigitizer(digiColl,simData,theGeom,validIds,digitizationType,engine);
87 }
88 
89 template<class DFr>
90 void HGCDigitizerBase<DFr>::runSimple(std::unique_ptr<HGCDigitizerBase::DColl> &coll,
91  HGCSimHitDataAccumulator &simData,
92  const CaloSubdetectorGeometry* theGeom,
93  const std::unordered_set<DetId>& validIds,
94  CLHEP::HepRandomEngine* engine) {
95  HGCSimHitData chargeColl,toa;
96 
97  // this represents a cell with no signal charge
98  HGCCellInfo zeroData;
99  zeroData.hit_info[0].fill(0.f); //accumulated energy
100  zeroData.hit_info[1].fill(0.f); //time-of-flight
101 
102  for( const auto& id : validIds ) {
103  chargeColl.fill(0.f);
104  toa.fill(0.f);
105  HGCSimHitDataAccumulator::iterator it = simData.find(id);
106  HGCCellInfo& cell = ( simData.end() == it ? zeroData : it->second );
107  addCellMetadata(cell,theGeom,id);
108 
109  for(size_t i=0; i<cell.hit_info[0].size(); i++) {
110  double rawCharge(cell.hit_info[0][i]);
111 
112  //time of arrival
113  toa[i]=cell.hit_info[1][i];
114  if(myFEelectronics_->toaMode()==HGCFEElectronics<DFr>::WEIGHTEDBYE && rawCharge>0)
115  toa[i]=cell.hit_info[1][i]/rawCharge;
116 
117  //convert total energy in GeV to charge (fC)
118  //double totalEn=rawEn*1e6*keV2fC_;
119  float totalCharge=rawCharge;
120 
121  //add noise (in fC)
122  //we assume it's randomly distributed and won't impact ToA measurement
123  totalCharge += std::max( (float)CLHEP::RandGaussQ::shoot(engine,0.0,cell.size*noise_fC_[cell.thickness-1]) , 0.f );
124  if(totalCharge<0.f) totalCharge=0.f;
125 
126  chargeColl[i]= totalCharge;
127  }
128 
129  //run the shaper to create a new data frame
130  DFr rawDataFrame( id );
131  if( !cce_.empty() )
132  myFEelectronics_->runShaper(rawDataFrame, chargeColl, toa, cell.thickness, engine, cce_[cell.thickness-1]);
133  else
134  myFEelectronics_->runShaper(rawDataFrame, chargeColl, toa, cell.thickness, engine);
135 
136  //update the output according to the final shape
137  updateOutput(coll,rawDataFrame);
138  }
139 }
140 
141 template<class DFr>
142 void HGCDigitizerBase<DFr>::updateOutput(std::unique_ptr<HGCDigitizerBase::DColl> &coll,
143  const DFr& rawDataFrame) {
144  int itIdx(9);
145  if(rawDataFrame.size()<=itIdx+2) return;
146 
147  DFr dataFrame( rawDataFrame.id() );
148  dataFrame.resize(5);
149  bool putInEvent(false);
150  for(int it=0;it<5; it++) {
151  dataFrame.setSample(it, rawDataFrame[itIdx-2+it]);
152  if(it==2) putInEvent = rawDataFrame[itIdx-2+it].threshold();
153  }
154 
155  if(putInEvent) {
156  coll->push_back(dataFrame);
157  }
158 }
159 
160 // cause the compiler to generate the appropriate code
162 template class HGCDigitizerBase<HGCEEDataFrame>;
163 template class HGCDigitizerBase<HGCBHDataFrame>;
T getParameter(std::string const &) const
static const TGPicture * info(bool iBackgroundIsBlack)
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< HGCSimHitData, 2 > hit_info
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
std::array< HGCSimData_t, nSamples > HGCSimHitData
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
std::unordered_map< uint32_t, HGCCellInfo > HGCSimHitDataAccumulator
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
const HGCalTopology & topology() const
Definition: HGCalGeometry.h:96
double f[11][100]
int wafer() const
get the wafer #
Definition: HGCalDetId.h:42
HGCDigitizerBase(const edm::ParameterSet &ps)
CTOR.
Definition: DetId.h:18
JetCorrectorParametersCollection coll
Definition: classes.h:10
const HGCalDDDConstants & dddConstants() const
susybsm::HSCParticleCollection hc
Definition: classes.h:25
models the behavior of the front-end electronics
Detector det() const
get the detector field from this detid
Definition: DetId.h:35