Hall C ROOT/C++ Analyzer (hcana)
Loading...
Searching...
No Matches
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 "TString.h"
189#include <iostream>
190
191using namespace std;
192
194 : fNPedestalSamples(4)
195 , fNPeakSamples(9)
196 , fPeakPedestalRatio(1.0*fNPeakSamples/fNPedestalSamples)
197 , fSubsampleToTimeFactor(0.0625)
198 , fSampThreshold(10.)
199 , fNSA(-1)
200 , fNSAT(-1)
201 , fNSB(-1)
202 , fPed(0)
203 , fSampPed(0)
204 , fPulseInt()
205 , fPulseAmp()
206 , fPulseTime()
207 , fSampPulseInt()
208 , fSampPulseAmp()
209 , fSampPulseTime()
210 , fSample()
211 , fRefTime(0)
212 , fRefDiffTime(0)
213 , fHasMulti(kFALSE)
214 , fHasRefTime(kFALSE)
215 , fNPulses(0)
216 , fNSampPulses(0)
217 , fNSamples(0)
218{}
219
221 TObject::operator=(right);
222
223 if (this != &right) {
224 fPed = right.fPed;
225 fSampPed = right.fSampPed;
226 for (UInt_t i=0; i<fMaxNPulses; ++i) {
227 fPulseInt[i] = right.fPulseInt[i];
228 fPulseAmp[i] = right.fPulseAmp[i];
229 fPulseTime[i] = right.fPulseTime[i];
230 fSampPulseInt[i] = right.fSampPulseInt[i];
231 fSampPulseAmp[i] = right.fSampPulseAmp[i];
232 fSampPulseTime[i] = right.fSampPulseTime[i];
233 }
234 for (UInt_t i=0; i<fMaxNSamples; ++i) {
235 fSample[i] = right.fSample[i];
236 }
237 fHasMulti = right.fHasMulti;
238 fNPulses = right.fNPulses;
240 fNSamples = right.fNSamples;
241 fRefTime = right.fRefTime;
242 fHasRefTime = right.fHasRefTime;
243 }
244
245 return *this;
246}
247
249
251 TObject::Clear(opt);
252
253 fPed = 0;
254 for (UInt_t i=0; i<fNPulses; ++i) {
255 fPulseInt[i] = 0;
256 fPulseAmp[i] = 0;
257 fPulseTime[i] = 0;
258 fSampPulseInt[i] = 0;
259 fSampPulseAmp[i] = 0;
260 fSampPulseTime[i] = 0;
261 }
262 for (UInt_t i=0; i<fNSamples; ++i) {
263 fSample[i] = 0 ;
264 }
266 fNPulses = 0;
267 fNSampPulses = 0;
268 fNSamples = 0;
269 fRefTime = 0;
271}
272
274 if (fNPulses >= fMaxNPulses) {
275 throw std::out_of_range(
276 "`THcRawAdcHit::SetData`: too many pulses!"
277 );
278 }
280 ++fNPulses;
281}
282
284 fRefTime = refTime;
286}
287
289 fRefDiffTime = refDiffTime;
290}
291
295
296
298 fNSAT = nsat;
299}
300
301
303 if (fNSamples >= fMaxNSamples) {
304 throw std::out_of_range(
305 "`THcRawAdcHit::SetSample`: too many samples!"
306 );
307 }
309 ++fNSamples;
310}
311
312
314 /*
315 if ( fNPulses ==0) {
316 std::cout << " in SetSampIntTimePedestalPeak " << fNPedestalSamples << " " << fNSamples<< std::endl;
317 if (fNPulses>0) {
318 for (UInt_t np=0;np<fNPulses;np++) {
319 std::cout << " Npulse = " << np+1 << " " << fPulseInt[np] << " " << GetPulseInt(np) << " " << fPulseAmp[np] << " " << GetPulseAmp(np) << " " << fPulseTime[np] << " " << GetPulseTime(np)<< std::endl;
320 }
321 } else {
322 std::cout << " No pulses found " << std::endl;
323 }
324 }
325 */
328 //
329 Int_t NS = fNPedestalSamples - 1;
330 fNSampPulses = 0;
331 if( fSampThreshold == 0 ) { // just want to get the pedestal for determining FADC threshold
333 TMath::Min((NS + fNSA - 1), int(fNSamples - 1)));
334 fNPeakSamples = TMath::Min((NS + fNSA - 1), int(fNSamples - 1)) - TMath::Max(NS - fNSB, 0) + 1;
337 fSampPulseTime[fNSampPulses] = 64 * NS;
338 fNSampPulses = 1;
339
340 } else {
341 //
342 Bool_t CheckSampBelowThres = kFALSE;
343 Int_t LastFourSampPed = GetIntegral(fNSamples - fNPedestalSamples, fNSamples - 1);
344 if ( LastFourSampPed < fSampPed) {
345 fSampPed= LastFourSampPed;
346 CheckSampBelowThres = kTRUE;
347 }
348 while( NS < int(fNSamples) && fNSampPulses < fMaxNPulses ) {
349 // GetSample(NS) Pedestal subtracted sample value (mV)
350 if( CheckSampBelowThres ) { // find sample below threshold to start search
351 if( GetSample(NS) < fSampThreshold ) CheckSampBelowThres = kFALSE;
352 } else {
353 Int_t ns_found = 0;
354 for( Int_t nt = NS; nt < TMath::Min((NS + fNSAT), int(fNSamples)); nt++ ) {
355 if( GetSample(nt) > fSampThreshold ) ns_found++;
356 }
357 if( ns_found == fNSAT ) { // NS is The TC bin
359 TMath::Min((NS + fNSA - 1), int(fNSamples - 1)));
360 fNPeakSamples = TMath::Min((NS + fNSA - 1), int(fNSamples - 1)) - TMath::Max(NS - fNSB, 0) + 1;
364 Int_t PeakBin = 0;
365 Double_t PeakVal = GetSample(NS);
366 for( Int_t nt = NS + 1; nt < TMath::Min((NS + fNSA), int(fNSamples)); nt++ ) {
367 if( GetSample(nt) < PeakVal && PeakBin == 0 ) {
368 PeakBin = nt - 1;
369 } else if( PeakBin == 0 ) {
370 PeakVal = GetSample(nt);
371 }
372 }
373 if( PeakBin > 0 ) {
375 Int_t Time = NS * 64; // Raw time in ADC bin number * 64 = "subsamples"
376 Double_t VMid = (GetSample(PeakBin)) / 2.;
377 for( Int_t nt = TMath::Max(NS - fNSB, 0); nt < TMath::Min(PeakBin, int(fNSamples - 1)); nt++ ) {
378 if( VMid >= GetSample(nt) && VMid < GetSample(nt + 1) ) {
379 Time = 64 * nt + int(64 * (VMid - GetSample(nt)) /
380 (GetSample(nt + 1) - GetSample(nt)));
381 }
382 }
384 } else {
386 fSampPulseTime[fNSampPulses] = 64 * NS;
387 }
388 // if (fNPulses ==0) std::cout << " NsampPulse = " << fNSampPulses+1 << " " << fSampPulseInt[fNSampPulses] << " " << GetSampPulseInt(fNSampPulses) << " npeak = " <<PeakBin << " " << fSampPulseAmp[fNSampPulses] << " " << GetSampPulseAmp(fNSampPulses) << " " << fSampPulseTime[fNSampPulses]<< " " << GetSampPulseTime(fNSampPulses) << " Vmid = " << VMid << " TC = " << NS << " TC+NSA-1 =" << NS+fNSA << std::endl;
389 fNSampPulses++;
390 NS = NS + fNSA;
391 CheckSampBelowThres = kTRUE;
392 }
393 }
394 NS++;
395 }
396 } // match else
397 /*
398 if ( fNSampPulses==0 && fNPulses==0) {
399 std::cout << " No sample pulses found" << std::endl;
400 for (Int_t nt=0;nt<fNSamples;nt++) {
401 cout << GetSample(nt) << " " ;
402 if (nt%10 ==0 && nt !=0) cout << endl;
403 }
404 cout << endl;
405 }
406 */
407}
408//
410 Int_t pedestal, Int_t peak ) {
411 if( fNPulses >= fMaxNPulses ) {
412 throw std::out_of_range(
413 "`THcRawAdcHit::SetDataTimePedestalPeak`: too many pulses!"
414 );
415 }
418 fPed = pedestal;
421 ++fNPulses;
422}
423
425 if( NSA < 0 || NSB < 0 || NPED < 0 ) {
426 TString msg = TString::Format(
427 "`THcRawAdcHit::SetF250Params`: One of the params is negative! NSA = %d NSB = %d NPED = %d",
428 NSA, NSB, NPED
429 );
430 throw std::invalid_argument(msg.Data());
431 }
432 fNPedestalSamples = NPED;
435 fNSA = NSA;
436 fNSB = NSB;
437 fNSAT = 2;
438}
439
440void THcRawAdcHit::throw_bad_index( const char* loc, UInt_t i, UInt_t N )
441{
442 throw std::out_of_range(Form("%s: Bad index %u %u", loc, i, N));
443}
444
445void THcRawAdcHit::throw_bad_range( const char* loc, UInt_t lo, UInt_t hi,
446 UInt_t N )
447{
448 throw std::out_of_range(Form("%s: bad range %u %u %u", loc, lo, hi, N));
449
450}
451
int Int_t
unsigned int UInt_t
uint32_t NSB
std::vector< uint32_t > peak
uint32_t NSA
uint32_t pedestal
uint32_t time
bool Bool_t
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
const char Option_t
#define N
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
ClassImp(VDC::AnalyticTTDConv) using namespace std
#define hi
char * Form(const char *fmt,...)
Class representing a single raw ADC hit.
Definition THcRawAdcHit.h:7
Bool_t fHasRefTime
void SetData(Int_t data)
Sets raw ADC value.
static void throw_bad_index(const char *loc, UInt_t i, UInt_t N)
Int_t fPulseTime[fMaxNPulses]
virtual ~THcRawAdcHit()
Destructor.
void SetSampNSAT(Int_t nsat)
void SetRefDiffTime(Int_t refDiffTime)
void SetDataTimePedestalPeak(Int_t data, Int_t time, Int_t pedestal, Int_t peak)
Sets various bits of ADC data.
THcRawAdcHit & operator=(const THcRawAdcHit &right)
Assignment operator.
void SetSampThreshold(Double_t thres)
Double_t GetSample(UInt_t iSample=0) const
static constexpr UInt_t fMaxNSamples
void SetSampIntTimePedestalPeak()
Int_t fSampPulseInt[fMaxNPulses]
Int_t fSample[fMaxNSamples]
Double_t fPeakPedestalRatio
static void throw_bad_range(const char *loc, UInt_t lo, UInt_t hi, UInt_t N)
Int_t fSampPulseAmp[fMaxNPulses]
static constexpr UInt_t fMaxNPulses
Int_t GetSampleRaw(UInt_t iSample=0) const
Gets raw sample. In channels.
Int_t fPulseInt[fMaxNPulses]
Double_t fSampThreshold
Int_t fPulseAmp[fMaxNPulses]
void SetSample(Int_t data)
Sets raw signal sample.
virtual void Clear(Option_t *opt="")
Clears variables before next event.
Int_t fNPeakSamples
Int_t GetIntegral(UInt_t iSampleLow, UInt_t iSampleHigh) const
Gets integral of raw samples. In channels.
THcRawAdcHit()
Constructor.
UInt_t fNSampPulses
void SetF250Params(Int_t NSA, Int_t NSB, Int_t NPED)
Sets F250 parameters used for pedestal subtraction.
Int_t fNPedestalSamples
void SetRefTime(Int_t refTime)
Sets reference time. In channels.
Int_t fSampPulseTime[fMaxNPulses]
virtual void Clear(Option_t *="")
TObject & operator=(const TObject &rhs)
const char * Data() const
Double_t Min(Double_t a, Double_t b)
Double_t Max(Double_t a, Double_t b)
STL namespace.