CMS 3D CMS Logo

ZdcSimpleRecAlgo_Run3.cc
Go to the documentation of this file.
3 #include <algorithm> // for "max"
4 #include <iostream>
5 #include <cmath>
6 
7 #include <Eigen/Dense> // for TemplateFit Method
8 
9 // #include "DataFormats/METReco/interface/HcalCaloFlagLabels.h"
10 
12 
13 void ZdcSimpleRecAlgo_Run3::initCorrectionMethod(const int method, const int ZdcSection) {
14  correctionMethod_[ZdcSection] = method;
15 };
16 
17 // Template fit is a linear combination of timeslices in a digi assuming there is a potential charge peak at each Bx
18 // and the charges of the Ts follwing a peak are consistent with the chargeRatios
19 void ZdcSimpleRecAlgo_Run3::initTemplateFit(const std::vector<unsigned int>& bxTs,
20  const std::vector<double>& chargeRatios,
21  const int nTs,
22  const int ZdcSection) {
23  int nRatios = chargeRatios.size();
24  int nBx = bxTs.size();
25  int nCols = nBx + 1;
26  double val = 0;
27  int index = 0;
28  int timeslice = 0;
29  nTs_ = nTs;
30  Eigen::MatrixXf a(nTs, nCols);
31  for (int j = 0; j < nBx; j++) {
32  timeslice = bxTs.at(j);
33  for (int i = 0; i < nTs; i++) {
34  val = 0;
35  index = i - timeslice;
36  if (index >= 0 && index < nRatios)
37  val = chargeRatios.at(index);
38  a(i, j) = val;
39  }
40  }
41  for (int i = 0; i < nTs; i++)
42  a(i, nBx) = 1;
43  Eigen::MatrixXf b = a.transpose() * a;
44  if (std::fabs(b.determinant()) < 1E-8) {
45  templateFitValid_[ZdcSection] = false;
46  return;
47  }
48  templateFitValid_[ZdcSection] = true;
49  Eigen::MatrixXf tfMatrix;
50  tfMatrix = b.inverse() * a.transpose();
51  for (int i = 0; i < nTs; i++)
52  templateFitValues_[ZdcSection].push_back(tfMatrix.coeff(1, i));
53 }
54 
55 void ZdcSimpleRecAlgo_Run3::initRatioSubtraction(const float ratio, const float frac, const int ZdcSection) {
56  ootpuRatio_[ZdcSection] = ratio;
57  ootpuFrac_[ZdcSection] = frac;
58 }
59 
60 // helper functions for pedestal subtraction and noise calculations
61 namespace zdchelper {
62 
63  inline double subPedestal(const float charge, const float ped, const float width) {
64  if (charge - ped > width)
65  return (charge - ped);
66  else
67  return (0);
68  }
69 
70  // gets the noise with a check if the ratio of energy0 / energy1 > ootpuRatio
71  // energy1 is noiseTS , energy0 is noiseTs-1
72  inline double getNoiseOOTPURatio(const float energy0,
73  const float energy1,
74  const float ootpuRatio,
75  const float ootpuFrac) {
76  if (energy0 >= ootpuRatio * energy1 || ootpuRatio < 0)
77  return (ootpuFrac * energy1);
78  else
79  return (energy1);
80  }
81 
82 } // namespace zdchelper
83 
85  const HcalCoder& coder,
86  const HcalCalibrations& calibs,
87  const HcalPedestal& effPeds,
88  const std::vector<unsigned int>& myNoiseTS,
89  const std::vector<unsigned int>& mySignalTS) const {
90  CaloSamples tool;
91  coder.adc2fC(digi, tool);
92  // Reads noiseTS and signalTS from database
93  int ifirst = mySignalTS[0];
94  double ampl = 0;
95  int maxI = -1;
96  double maxA = -1e10;
97  double ta = 0;
98  double energySOIp1 = 0;
99  double ratioSOIp1 = -1.0;
100  double chargeWeightedTime = -99.0;
101 
102  double Allnoise = 0;
103  int noiseslices = 0;
104  int CurrentTS = 0;
105  double noise = 0;
106  int digi_size = digi.samples();
107  HcalZDCDetId cell = digi.id();
108  int zdcsection = cell.section();
109 
110  // determining noise
111  for (unsigned int iv = 0; iv < myNoiseTS.size(); ++iv) {
112  CurrentTS = myNoiseTS[iv];
113  int capid = digi[CurrentTS].capid();
114  float ped = effPeds.getValue(capid);
115  float width = effPeds.getWidth(capid);
116  float gain = calibs.respcorrgain(capid);
117  if (CurrentTS >= digi_size)
118  continue;
119  float energy1 = zdchelper::subPedestal(tool[CurrentTS], ped, width) * gain;
120  float energy0 = 0;
121  if (CurrentTS > 0) {
122  capid = digi[CurrentTS - 1].capid();
123  ped = effPeds.getValue(capid);
124  width = effPeds.getWidth(capid);
125  gain = calibs.respcorrgain(capid);
126  energy0 = zdchelper::subPedestal(tool[CurrentTS - 1], ped, width) * gain;
127  }
128  Allnoise += zdchelper::getNoiseOOTPURatio(energy0, energy1, ootpuRatio_.at(zdcsection), ootpuFrac_.at(zdcsection));
129  noiseslices++;
130  }
131  if (noiseslices != 0) {
132  noise = (Allnoise) / double(noiseslices);
133  }
134 
135  // determining signal energy and max Ts
136  for (unsigned int ivs = 0; ivs < mySignalTS.size(); ++ivs) {
137  CurrentTS = mySignalTS[ivs];
138  if (CurrentTS >= digi_size)
139  continue;
140  float energy1 = -1;
141  int capid = digi[CurrentTS].capid();
142  // float ped = calibs.pedestal(capid);
143  float ped = effPeds.getValue(capid);
144  float width = effPeds.getWidth(capid);
145  float gain = calibs.respcorrgain(capid);
146  float energy0 = std::max(0.0, zdchelper::subPedestal(tool[CurrentTS], ped, width)) * gain;
147  if (CurrentTS < digi_size - 1) {
148  capid = digi[CurrentTS].capid();
149  ped = effPeds.getValue(capid);
150  width = effPeds.getWidth(capid);
151  gain = calibs.respcorrgain(capid);
152  energy1 = std::max(0.0, zdchelper::subPedestal(tool[CurrentTS + 1], ped, width)) * gain;
153  }
154  ta = energy0 - noise;
155  if (ta > 0)
156  ampl += ta;
157 
158  if (ta > maxA) {
159  ratioSOIp1 = (energy0 > 0 && energy1 > 0) ? energy0 / energy1 : -1.0;
160  maxA = ta;
161  maxI = CurrentTS;
162  }
163  }
164 
165  // determine energy if using Template Fit method
166  if (correctionMethod_.at(zdcsection) == 1 && templateFitValid_.at(zdcsection)) {
167  double energy = 0;
168  for (int iv = 0; iv < nTs_; iv++) {
169  int capid = digi[iv].capid();
170  float ped = effPeds.getValue(capid);
171  float width = effPeds.getWidth(capid);
172  float gain = calibs.respcorrgain(capid);
173  if (iv >= digi_size)
174  continue;
175  energy += zdchelper::subPedestal(tool[iv], ped, width) * (templateFitValues_.at(zdcsection)).at(iv) * gain;
176  }
177  ampl = std::max(0.0, energy);
178  }
179 
180  double time = -9999;
181  // Time based on regular energy
183  if (maxI == 0 || maxI == (tool.size() - 1)) {
184  LogDebug("HCAL Pulse") << "ZdcSimpleRecAlgo::reco2 :"
185  << " Invalid max amplitude position, "
186  << " max Amplitude: " << maxI << " first: " << ifirst << " last: " << (tool.size() - 1)
187  << std::endl;
188  } else {
189  int capid = digi[maxI - 1].capid();
190  double Energy0 = std::max(0.0, ((tool[maxI - 1]) * calibs.respcorrgain(capid)));
191 
192  capid = digi[maxI].capid();
193  double Energy1 = std::max(0.0, ((tool[maxI]) * calibs.respcorrgain(capid)));
194  capid = digi[maxI + 1].capid();
195  double Energy2 = std::max(0.0, ((tool[maxI + 1]) * calibs.respcorrgain(capid)));
196 
197  double TSWeightEnergy = ((maxI - 1) * Energy0 + maxI * Energy1 + (maxI + 1) * Energy2);
198  double EnergySum = Energy0 + Energy1 + Energy2;
199  double AvgTSPos = 0.;
200  if (EnergySum != 0)
201  AvgTSPos = TSWeightEnergy / EnergySum;
202  // If time is zero, set it to the "nonsensical" -99
203  // Time should be between 75ns and 175ns (Timeslices 3-7)
204  if (AvgTSPos == 0) {
205  time = -99;
206  } else {
207  time = (AvgTSPos * 25.0);
208  }
209  }
210 
211  // find energy for signal TS + 1
212  for (unsigned int ivs = 0; ivs < mySignalTS.size(); ++ivs) {
213  CurrentTS = mySignalTS[ivs] + 1;
214  if (CurrentTS >= digi_size)
215  continue;
216  int capid = digi[CurrentTS].capid();
217  // ta = tool[CurrentTS] - noise;
218  ta = tool[CurrentTS];
219  ta *= calibs.respcorrgain(capid); // fC --> GeV
220  if (ta > 0)
221  energySOIp1 += ta;
222  }
223 
224  double tmp_energy = 0;
225  double tmp_TSWeightedEnergy = 0;
226  for (int ts = 0; ts < digi_size; ++ts) {
227  int capid = digi[ts].capid();
228 
229  // max sure there are no negative values in time calculation
230  ta = std::max(0.0, tool[ts]);
231  ta *= calibs.respcorrgain(capid); // fC --> GeV
232  if (ta > 0) {
233  tmp_energy += ta;
234  tmp_TSWeightedEnergy += (ts + 1) * ta;
235  }
236  }
237 
238  if (tmp_energy > 0)
239  chargeWeightedTime = (tmp_TSWeightedEnergy / tmp_energy - 1) * 25.0;
240  auto rh = ZDCRecHit(digi.id(), ampl, time, -99);
241  rh.setEnergySOIp1(energySOIp1);
242 
243  if (maxI >= 0 && maxI < tool.size()) {
244  float tmp_tdctime = 0;
245  // TDC error codes will be 60=-1, 61 = -2, 62 = -3, 63 = -4
246  if (digi[maxI].le_tdc() >= 60)
247  tmp_tdctime = -1 * (digi[maxI].le_tdc() - 59);
248  else
249  tmp_tdctime = maxI * 25. + (digi[maxI].le_tdc() / 2);
250  rh.setTDCtime(tmp_tdctime);
251  }
252 
253  rh.setChargeWeightedTime(chargeWeightedTime);
254  rh.setRatioSOIp1(ratioSOIp1);
255  return rh;
256 }
257 
259  const std::vector<unsigned int>& myNoiseTS,
260  const std::vector<unsigned int>& mySignalTS,
261  const HcalCoder& coder,
262  const HcalCalibrations& calibs,
263  const HcalPedestal& effPeds) const {
264  return ZdcSimpleRecAlgo_Run3::reco0(digi, coder, calibs, effPeds, myNoiseTS, mySignalTS);
265 
266  edm::LogError("ZDCSimpleRecAlgoImpl::reconstruct, recoMethod was not declared");
267  throw cms::Exception("ZDCSimpleRecoAlgos::exiting process");
268 }
constexpr edm::DataFrame::id_type id() const
int size() const
get the size
Definition: CaloSamples.h:24
ZDCRecHit reconstruct(const QIE10DataFrame &digi, const std::vector< unsigned int > &myNoiseTS, const std::vector< unsigned int > &mySignalTS, const HcalCoder &coder, const HcalCalibrations &calibs, const HcalPedestal &effPeds) const
double subPedestal(const float charge, const float ped, const float width)
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalPedestal.h:20
Log< level::Error, false > LogError
std::map< int, std::vector< double > > templateFitValues_
void initRatioSubtraction(const float ratio, const float frac, const int ZdcSection)
std::map< int, bool > templateFitValid_
std::map< int, int > correctionMethod_
constexpr double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
ZdcSimpleRecAlgo_Run3(int recoMethod)
float getWidth(int fCapId) const
get width for capId = 0..3
Definition: HcalPedestal.h:25
constexpr void setEnergySOIp1(const float en)
Definition: ZDCRecHit.h:22
double b
Definition: hdecay.h:120
deadvectors [0] push_back({0.0175431, 0.538005, 6.80997, 13.29})
constexpr Section section() const
get the section
Definition: HcalZDCDetId.h:92
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
void initTemplateFit(const std::vector< unsigned int > &bxTs, const std::vector< double > &chargeRatios, const int nTs, const int ZdcSection)
double a
Definition: hdecay.h:121
std::map< int, float > ootpuRatio_
std::map< int, float > ootpuFrac_
void initCorrectionMethod(const int method, const int ZdcSection)
constexpr int samples() const
total number of samples in the digi
double getNoiseOOTPURatio(const float energy0, const float energy1, const float ootpuRatio, const float ootpuFrac)
ZDCRecHit reco0(const QIE10DataFrame &digi, const HcalCoder &coder, const HcalCalibrations &calibs, const HcalPedestal &effPeds, const std::vector< unsigned int > &myNoiseTS, const std::vector< unsigned int > &mySignalTS) const
#define LogDebug(id)