CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
ZdcSimpleRecAlgo_Run3 Class Reference

#include <ZdcSimpleRecAlgo_Run3.h>

Public Member Functions

void initCorrectionMethod (const int method, const int ZdcSection)
 
void initRatioSubtraction (const float ratio, const float frac, const int ZdcSection)
 
void initTemplateFit (const std::vector< unsigned int > &bxTs, const std::vector< double > &chargeRatios, const int nTs, const int ZdcSection)
 
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
 
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
 
 ZdcSimpleRecAlgo_Run3 (int recoMethod)
 

Private Attributes

std::map< int, int > correctionMethod_
 
int nTs_
 
std::map< int, float > ootpuFrac_
 
std::map< int, float > ootpuRatio_
 
int recoMethod_
 
std::map< int, bool > templateFitValid_
 
std::map< int, std::vector< double > > templateFitValues_
 

Detailed Description

This class reconstructs RecHits from Digis for ZDC by addition of selected time samples, pedestal subtraction, and gain application. The time of the hit is reconstructed using a weighted peak bin calculation supplemented by precise time lookup table. A consumer of this class also has the option of correcting the reconstructed time for energy-dependent time slew associated with the QIE.

A sencon method based on a based on a event by event substraction is also implelented. signal = (S4 + S5 - 2*(S1+S2+S3 + S7+S8+S9+S10))*(ft-Gev constant) where SN is the signal in the nth time slice

Author
E. Garcia CSU & J. Gomez UMD

Definition at line 39 of file ZdcSimpleRecAlgo_Run3.h.

Constructor & Destructor Documentation

◆ ZdcSimpleRecAlgo_Run3()

ZdcSimpleRecAlgo_Run3::ZdcSimpleRecAlgo_Run3 ( int  recoMethod)

Simple constructor for PMT-based detectors

Definition at line 11 of file ZdcSimpleRecAlgo_Run3.cc.

Member Function Documentation

◆ initCorrectionMethod()

void ZdcSimpleRecAlgo_Run3::initCorrectionMethod ( const int  method,
const int  ZdcSection 
)

◆ initRatioSubtraction()

void ZdcSimpleRecAlgo_Run3::initRatioSubtraction ( const float  ratio,
const float  frac,
const int  ZdcSection 
)

Definition at line 55 of file ZdcSimpleRecAlgo_Run3.cc.

References DivergingColor::frac, ootpuFrac_, and ootpuRatio_.

Referenced by ZdcHitReconstructor_Run3::ZdcHitReconstructor_Run3().

55  {
56  ootpuRatio_[ZdcSection] = ratio;
57  ootpuFrac_[ZdcSection] = frac;
58 }
std::map< int, float > ootpuRatio_
std::map< int, float > ootpuFrac_

◆ initTemplateFit()

void ZdcSimpleRecAlgo_Run3::initTemplateFit ( const std::vector< unsigned int > &  bxTs,
const std::vector< double > &  chargeRatios,
const int  nTs,
const int  ZdcSection 
)

Definition at line 19 of file ZdcSimpleRecAlgo_Run3.cc.

References a, b, mps_fire::i, dqmiolumiharvest::j, nTs_, push_back(), templateFitValid_, templateFitValues_, and heppy_batch::val.

Referenced by ZdcHitReconstructor_Run3::ZdcHitReconstructor_Run3().

22  {
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 }
std::map< int, std::vector< double > > templateFitValues_
std::map< int, bool > templateFitValid_
double b
Definition: hdecay.h:120
deadvectors [0] push_back({0.0175431, 0.538005, 6.80997, 13.29})
double a
Definition: hdecay.h:121

◆ reco0()

ZDCRecHit ZdcSimpleRecAlgo_Run3::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

Definition at line 84 of file ZdcSimpleRecAlgo_Run3.cc.

References HcalCoder::adc2fC(), correctionMethod_, HBHEDarkening_cff::energy, Hcal_Conditions_forGlobalTag_cff::gain, zdchelper::getNoiseOOTPURatio(), HcalPedestal::getValue(), HcalPedestal::getWidth(), QIE10DataFrame::id(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::iv, LogDebug, WZElectronSkims53X_cff::max, findQualityFiles::maxI, hgchebackDigitizer_cfi::noise, nTs_, ootpuFrac_, ootpuRatio_, QIE10Task_cfi::ped, HcalCalibrations::respcorrgain(), QIE10DataFrame::samples(), HcalZDCDetId::section(), ZDCRecHit::setEnergySOIp1(), CaloSamples::size(), zdchelper::subPedestal(), templateFitValid_, templateFitValues_, hcalRecHitTable_cff::time, and ApeEstimator_cff::width.

Referenced by reconstruct().

89  {
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)*ta;
235  }
236  }
237 
238  if (tmp_energy > 0)
239  chargeWeightedTime = (tmp_TSWeightedEnergy / tmp_energy) * 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  int le_tdc = digi[maxI].le_tdc();
246  // TDC error codes will be 60=-1, 61 = -2, 62 = -3, 63 = -4
247  if (le_tdc >= 60)
248  tmp_tdctime = -1 * (le_tdc - 59);
249  else
250  tmp_tdctime = maxI * 25. + (le_tdc / 2.0);
251  rh.setTDCtime(tmp_tdctime);
252  }
253 
254  rh.setChargeWeightedTime(chargeWeightedTime);
255  rh.setRatioSOIp1(ratioSOIp1);
256  return rh;
257 }
constexpr edm::DataFrame::id_type id() const
int size() const
get the size
Definition: CaloSamples.h:24
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
std::map< int, std::vector< double > > templateFitValues_
std::map< int, bool > templateFitValid_
std::map< int, int > correctionMethod_
constexpr double respcorrgain(int fCapId) const
get response corrected gain for capid=0..3
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
constexpr Section section() const
get the section
Definition: HcalZDCDetId.h:92
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const =0
std::map< int, float > ootpuRatio_
std::map< int, float > ootpuFrac_
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)
#define LogDebug(id)

◆ reconstruct()

ZDCRecHit ZdcSimpleRecAlgo_Run3::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

Definition at line 259 of file ZdcSimpleRecAlgo_Run3.cc.

References Exception, and reco0().

Referenced by ZdcHitReconstructor_Run3::produce().

264  {
265  return ZdcSimpleRecAlgo_Run3::reco0(digi, coder, calibs, effPeds, myNoiseTS, mySignalTS);
266 
267  edm::LogError("ZDCSimpleRecAlgoImpl::reconstruct, recoMethod was not declared");
268  throw cms::Exception("ZDCSimpleRecoAlgos::exiting process");
269 }
Log< level::Error, false > LogError
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

Member Data Documentation

◆ correctionMethod_

std::map<int, int> ZdcSimpleRecAlgo_Run3::correctionMethod_
private

Definition at line 71 of file ZdcSimpleRecAlgo_Run3.h.

Referenced by initCorrectionMethod(), and reco0().

◆ nTs_

int ZdcSimpleRecAlgo_Run3::nTs_
private

Definition at line 66 of file ZdcSimpleRecAlgo_Run3.h.

Referenced by initTemplateFit(), and reco0().

◆ ootpuFrac_

std::map<int, float> ZdcSimpleRecAlgo_Run3::ootpuFrac_
private

Definition at line 70 of file ZdcSimpleRecAlgo_Run3.h.

Referenced by initRatioSubtraction(), and reco0().

◆ ootpuRatio_

std::map<int, float> ZdcSimpleRecAlgo_Run3::ootpuRatio_
private

Definition at line 69 of file ZdcSimpleRecAlgo_Run3.h.

Referenced by initRatioSubtraction(), and reco0().

◆ recoMethod_

int ZdcSimpleRecAlgo_Run3::recoMethod_
private

Definition at line 65 of file ZdcSimpleRecAlgo_Run3.h.

◆ templateFitValid_

std::map<int, bool> ZdcSimpleRecAlgo_Run3::templateFitValid_
private

Definition at line 68 of file ZdcSimpleRecAlgo_Run3.h.

Referenced by initTemplateFit(), and reco0().

◆ templateFitValues_

std::map<int, std::vector<double> > ZdcSimpleRecAlgo_Run3::templateFitValues_
private

Definition at line 67 of file ZdcSimpleRecAlgo_Run3.h.

Referenced by initTemplateFit(), and reco0().