Hall C ROOT/C++ Analyzer (hcana)
THcRawAdcHit.cxx
Go to the documentation of this file.
1 
184 // TODO: Disallow using both SetData and SetDataTimePedestalPeak.
185 
186 
187 #include "THcRawAdcHit.h"
188 #include <stdexcept>
189 #include "TString.h"
190  const Double_t THcRawAdcHit::fNAdcChan = 4096.0; // Number of FADC channels in units of ADC channels
191  const Double_t THcRawAdcHit::fAdcRange = 1.0; // Dynamic range of FADCs in units of V, // TO-DO: Get fAdcRange from pre-start event
192  const Double_t THcRawAdcHit::fAdcImpedence = 50.0; // FADC input impedence in units of Ohms
193  const Double_t THcRawAdcHit::fAdcTimeSample = 4000.0; // Length of FADC time sample in units of ps
194  const Double_t THcRawAdcHit::fAdcTimeRes = 0.0625; // FADC time resolution in units of ns
195 
197  TObject(),
198  fNPedestalSamples(4), fNPeakSamples(9),
199  fPeakPedestalRatio(1.0*fNPeakSamples/fNPedestalSamples),
200  fSubsampleToTimeFactor(0.0625),
201  fPed(0), fPulseInt(), fPulseAmp(), fPulseTime(), fSample(),
202  fRefTime(0), fHasMulti(kFALSE), fHasRefTime(kFALSE), fNPulses(0), fNSamples(0)
203 {}
204 
206  TObject::operator=(right);
207 
208  if (this != &right) {
209  fPed = right.fPed;
210  for (UInt_t i=0; i<fMaxNPulses; ++i) {
211  fPulseInt[i] = right.fPulseInt[i];
212  fPulseAmp[i] = right.fPulseAmp[i];
213  fPulseTime[i] = right.fPulseTime[i];
214  }
215  for (UInt_t i=0; i<fMaxNSamples; ++i) {
216  fSample[i] = right.fSample[i];
217  }
218  fHasMulti = right.fHasMulti;
219  fNPulses = right.fNPulses;
220  fNSamples = right.fNSamples;
221  fRefTime = right.fRefTime;
222  fHasRefTime = right.fHasRefTime;
223  }
224 
225  return *this;
226 }
227 
229 
231  TObject::Clear(opt);
232 
233  fPed = 0;
234  for (UInt_t i=0; i<fNPulses; ++i) {
235  fPulseInt[i] = 0;
236  fPulseAmp[i] = 0;
237  fPulseTime[i] = 0;
238  }
239  for (UInt_t i=0; i<fNSamples; ++i) {
240  fSample[i] = 0 ;
241  }
242  fHasMulti = kFALSE;
243  fNPulses = 0;
244  fNSamples = 0;
245  fRefTime = 0;
247 }
248 
250  if (fNPulses >= fMaxNPulses) {
251  throw std::out_of_range(
252  "`THcRawAdcHit::SetData`: too many pulses!"
253  );
254  }
255  fPulseInt[fNPulses] = data;
256  ++fNPulses;
257 }
258 
260  fRefTime = refTime;
261  fHasRefTime = kTRUE;
262 }
263 
265  fRefDiffTime = refDiffTime;
266 }
267 
269  if (fNSamples >= fMaxNSamples) {
270  throw std::out_of_range(
271  "`THcRawAdcHit::SetSample`: too many samples!"
272  );
273  }
274  fSample[fNSamples] = data;
275  ++fNSamples;
276 }
277 
279  Int_t data, Int_t time, Int_t pedestal, Int_t peak
280 ) {
281  if (fNPulses >= fMaxNPulses) {
282  throw std::out_of_range(
283  "`THcRawAdcHit::SetDataTimePedestalPeak`: too many pulses!"
284  );
285  }
286  fPulseInt[fNPulses] = data;
287  fPulseTime[fNPulses] = time;
288  fPed = pedestal;
289  fPulseAmp[fNPulses] = peak;
290  fHasMulti = kTRUE;
291  ++fNPulses;
292 }
293 
295  if (iPulse >= fNPulses && iPulse != 0) {
296  TString msg = TString::Format(
297  "`THcRawAdcHit::GetRawData`: requested pulse %d where only %d pulses available!",
298  iPulse, fNPulses
299  );
300  throw std::out_of_range(msg.Data());
301  }
302  else if (iPulse >= fNPulses && iPulse == 0) {
303  return 0;
304  }
305  else {
306  return fPulseInt[iPulse];
307  }
308 }
309 
310 Double_t THcRawAdcHit::GetAverage(UInt_t iSampleLow, UInt_t iSampleHigh) const {
311  if (iSampleHigh >= fNSamples || iSampleLow >= fNSamples) {
312  TString msg = TString::Format(
313  "`THcRawAdcHit::GetAverage`: not this many samples available!"
314  );
315  throw std::out_of_range(msg.Data());
316  }
317  else {
318  Double_t average = 0.0;
319  for (UInt_t i=iSampleLow; i<=iSampleHigh; ++i) {
320  average += fSample[i];
321  }
322  return average / (iSampleHigh - iSampleLow + 1);
323  }
324 }
325 
326 
327 Int_t THcRawAdcHit::GetIntegral(UInt_t iSampleLow, UInt_t iSampleHigh) const {
328  if (iSampleHigh >= fNSamples || iSampleLow >= fNSamples) {
329  TString msg = TString::Format(
330  "`THcRawAdcHit::GetAverage`: not this many samples available!"
331  );
332  throw std::out_of_range(msg.Data());
333  }
334  else {
335  Int_t integral = 0;
336  for (UInt_t i=iSampleLow; i<=iSampleHigh; ++i) {
337  integral += fSample[i];
338  }
339  return integral;
340  }
341 }
342 
344  UInt_t iPedLow, UInt_t iPedHigh, UInt_t iIntLow, UInt_t iIntHigh
345 ) const {
346  return
347  GetIntegral(iIntLow, iIntHigh)
348  - GetAverage(iPedHigh, iPedLow) * (iIntHigh - iIntLow + 1);
349 }
350 
352  return fNPulses;
353 }
354 
356  return fNSamples;
357 }
358 
360  return fHasMulti;
361 }
362 
364  return fPed;
365 }
366 
368  if (iPulse < fNPulses) {
369  return fPulseInt[iPulse];
370  }
371  else if (iPulse == 0) {
372  return 0;
373  }
374  else {
375  TString msg = TString::Format(
376  "`THcRawAdcHit::GetPulseIntRaw`: Trying to get pulse %d where only %d pulses available!",
377  iPulse, fNPulses
378  );
379  throw std::out_of_range(msg.Data());
380  }
381 }
382 
384  if (iPulse < fNPulses) {
385  return fPulseAmp[iPulse];
386  }
387  else if (iPulse == 0) {
388  return 0;
389  }
390  else {
391  TString msg = TString::Format(
392  "`THcRawAdcHit::GetPulseAmpRaw`: Trying to get pulse %d where only %d pulses available!",
393  iPulse, fNPulses
394  );
395  throw std::out_of_range(msg.Data());
396  }
397 }
398 
400  if (iPulse < fNPulses) {
401  return fPulseTime[iPulse];
402  }
403  else if (iPulse == 0) {
404  return 0;
405  }
406  else {
407  TString msg = TString::Format(
408  "`THcRawAdcHit::GetPulseTimeRaw`: Trying to get pulse %d where only %d pulses available!",
409  iPulse, fNPulses
410  );
411  throw std::out_of_range(msg.Data());
412  }
413 }
414 
416  if (iSample < fNSamples) {
417  return fSample[iSample];
418  }
419  else {
420  TString msg = TString::Format(
421  "`THcRawAdcHit::GetSampleRaw`: Trying to get sample %d where only %d samples available!",
422  iSample, fNSamples
423  );
424  throw std::out_of_range(msg.Data());
425  }
426 }
427 
429  return (static_cast<Double_t>(fPed)/static_cast<Double_t>(fNPedestalSamples))*GetAdcTomV();
430 }
431 
433  return (static_cast<Double_t>(fPulseInt[iPulse]) - static_cast<Double_t>(fPed)*fPeakPedestalRatio)*GetAdcTopC();
434 }
435 
437  return (static_cast<Double_t>(fPulseAmp[iPulse]) - static_cast<Double_t>(fPed)/static_cast<Double_t>(fNPedestalSamples))*GetAdcTomV();
438 }
439 
441  Int_t rawtime = fPulseTime[iPulse];
442  if (fHasRefTime) {
443  rawtime -= fRefTime;
444  }
445  return (static_cast<Double_t>(rawtime)*GetAdcTons());
446 }
447 
449  Int_t integral = 0;
450 
451  for (UInt_t iSample=0; iSample<fNSamples; ++iSample) {
452  integral += fSample[iSample];
453  }
454 
455  return integral;
456 }
457 
459  return static_cast<Double_t>(GetSampleIntRaw()) - GetPed()*static_cast<Double_t>(fNSamples);
460 }
461 
463  if (NSA < 0 || NSB < 0 || NPED < 0) {
464  TString msg = TString::Format(
465  "`THcRawAdcHit::SetF250Params`: One of the params is negative! NSA = %d NSB = %d NPED = %d",
466  NSA, NSB, NPED
467  );
468  throw std::invalid_argument(msg.Data());
469  }
470  fNPedestalSamples = NPED;
471  fNPeakSamples = NSA + NSB;
473 }
474 
475 // FADC conversion factors
476 // Convert pedestal and amplitude to mV
478  // 1000 mV / 4096 ADC channels
479  return (fAdcRange*1000. / fNAdcChan);
480 }
481 
482 // Convert integral to pC
484  // (1 V / 4096 adc channels) * (4000 ps time sample / 50 ohms input resistance) = 0.020 pc/channel
486 }
487 
488 // Convert time sub samples to ns
490  return fAdcTimeRes;
491 }
492 
494  if (fHasRefTime) {
495  return fRefTime;
496  }
497  else {
498  TString msg = TString::Format(
499  "`THcRawAdcHit::GetRefTime`: Reference time not available!"
500  );
501  throw std::runtime_error(msg.Data());
502  }
503 }
504 
506  if (fHasRefTime) {
507  return fRefDiffTime;
508  }
509  else {
510  TString msg = TString::Format(
511  "`THcRawAdcHit::GetRefTime`: Reference time not available!"
512  );
513  throw std::runtime_error(msg.Data());
514  }
515 }
516 
518  return fHasRefTime;
519 }
520 
Double_t GetData(UInt_t iPedLow, UInt_t iPedHigh, UInt_t iIntLow, UInt_t iIntHigh) const
Gets pedestal subtracted integral of samples. In channels.
virtual void Clear(Option_t *="")
static const Double_t fNAdcChan
Definition: THcRawAdcHit.h:66
Int_t fRefTime
Definition: THcRawAdcHit.h:82
static const Double_t fAdcTimeSample
Definition: THcRawAdcHit.h:69
Double_t GetSampleInt() const
Gets pedestal subtracted integral of samples. In channels.
UInt_t GetNSamples() const
Gets number of set samples.
Int_t fNPeakSamples
Definition: THcRawAdcHit.h:73
void SetSample(Int_t data)
Sets raw signal sample.
const char Option_t
Bool_t fHasRefTime
Definition: THcRawAdcHit.h:86
Int_t GetRawData(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
Int_t GetPulseTimeRaw(UInt_t iPulse=0) const
Gets raw pulse time. In subsamples.
int Int_t
bool Bool_t
Double_t GetAverage(UInt_t iSampleLow, UInt_t iSampleHigh) const
Gets average of raw samples. In channels.
Int_t GetSampleRaw(UInt_t iSample=0) const
Gets raw sample. In channels.
static const Double_t fAdcRange
Definition: THcRawAdcHit.h:67
void SetRefTime(Int_t refTime)
Sets reference time. In channels.
const char * Data() const
static TString Format(const char *fmt,...)
Int_t fPulseAmp[fMaxNPulses]
Definition: THcRawAdcHit.h:79
UInt_t fNPulses
Definition: THcRawAdcHit.h:87
Int_t fPulseTime[fMaxNPulses]
Definition: THcRawAdcHit.h:80
Bool_t fHasMulti
Definition: THcRawAdcHit.h:85
static const UInt_t fMaxNPulses
Definition: THcRawAdcHit.h:62
TObject & operator=(const TObject &rhs)
Double_t GetAdcTons() const
Double_t GetPed() const
Gets sample pedestal. In channels.
void SetF250Params(Int_t NSA, Int_t NSB, Int_t NPED)
Sets F250 parameters used for pedestal subtraction.
Bool_t HasMulti() const
Queries whether data is from flash 250 module.
Int_t GetPulseAmpRaw(UInt_t iPulse=0) const
Gets raw pulse amplitude. In channels.
Int_t fSample[fMaxNSamples]
Definition: THcRawAdcHit.h:81
Double_t GetPulseAmp(UInt_t iPulse=0) const
Gets pedestal subtracted pulse amplitude. In channels.
Double_t GetPulseInt(UInt_t iPulse=0) const
Gets pedestal subtracted pulse integral. In channels.
void SetRefDiffTime(Int_t refDiffTime)
unsigned int UInt_t
Int_t GetSampleIntRaw() const
Gets raw integral of sTimeFacamples. In channels.
static const UInt_t fMaxNSamples
Definition: THcRawAdcHit.h:63
static const Double_t fAdcImpedence
Definition: THcRawAdcHit.h:68
Int_t fPulseInt[fMaxNPulses]
Definition: THcRawAdcHit.h:78
void SetDataTimePedestalPeak(Int_t data, Int_t time, Int_t pedestal, Int_t peak)
Sets various bits of ADC data.
static const Double_t fAdcTimeRes
Definition: THcRawAdcHit.h:70
const Bool_t kFALSE
THcRawAdcHit()
Constructor.
virtual void Clear(Option_t *opt="")
Clears variables before next event.
Int_t GetPulseIntRaw(UInt_t iPulse=0) const
Gets raw pulse integral. In channels.
Int_t GetRefTime() const
double Double_t
Bool_t HasRefTime() const
UInt_t GetNPulses() const
Gets number of set pulses.
void SetData(Int_t data)
Sets raw ADC value.
Int_t fRefDiffTime
Definition: THcRawAdcHit.h:83
Int_t GetIntegral(UInt_t iSampleLow, UInt_t iSampleHigh) const
Gets integral of raw samples. In channels.
ClassImp(THcDCLookupTTDConv) THcDCLookupTTDConv
Double_t fPeakPedestalRatio
Definition: THcRawAdcHit.h:74
Int_t GetPedRaw() const
Gets raw signal pedestal. In channels.
THcRawAdcHit & operator=(const THcRawAdcHit &right)
Assignment operator.
Double_t GetPulseTime(UInt_t iPulse=0) const
Int_t fNPedestalSamples
Definition: THcRawAdcHit.h:72
Int_t GetRefDiffTime() const
virtual ~THcRawAdcHit()
Destructor.
const Bool_t kTRUE
Double_t GetAdcTopC() const
Class representing a single raw ADC hit.
Definition: THcRawAdcHit.h:7
Double_t GetAdcTomV() const
UInt_t fNSamples
Definition: THcRawAdcHit.h:88