Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:17:52

0001 #include "MbdSig.h"
0002 #include "MbdCalib.h"
0003 
0004 #include <phool/phool.h>
0005 
0006 #include <TF1.h>
0007 #include <TFile.h>
0008 #include <TGraphErrors.h>
0009 #include <TH2.h>
0010 #include <TMath.h>
0011 #include <TPad.h>
0012 #include <TSpectrum.h>
0013 #include <TSpline.h>
0014 #include <TTree.h>
0015 
0016 #include <algorithm>
0017 #include <fstream>
0018 #include <iostream>
0019 #include <iomanip>
0020 #include <limits>
0021 
0022 MbdSig::MbdSig(const int chnum, const int nsamp)
0023   : _ch{chnum}
0024   , _nsamples{nsamp}
0025 {
0026   // std::cout << "In MbdSig::MbdSig(" << _ch << "," << _nsamples << ")" << std::endl;
0027 }
0028 
0029 void MbdSig::Init()
0030 {
0031   TString name;
0032 
0033   name = "hrawpulse";
0034   name += _ch;
0035   hRawPulse = new TH1F(name, name, _nsamples, -0.5, _nsamples - 0.5);
0036   name = "hsubpulse";
0037   name += _ch;
0038   hSubPulse = new TH1F(name, name, _nsamples, -0.5, _nsamples - 0.5);
0039 
0040   // gRawPulse = new TGraphErrors(_nsamples);
0041   gRawPulse = new TGraphErrors();
0042   name = "grawpulse";
0043   name += _ch;
0044   gRawPulse->SetName(name);
0045   // gSubPulse = new TGraphErrors(_nsamples);
0046   gSubPulse = new TGraphErrors();
0047   name = "gsubpulse";
0048   name += _ch;
0049   gSubPulse->SetName(name);
0050 
0051   hpulse = hRawPulse;  // hpulse,gpulse point to raw by default
0052   gpulse = gRawPulse;  // we switch to sub for default if ped is applied
0053 
0054   //ped0stats = std::make_unique<MbdRunningStats>(100);  // use the last 100 events for running pedestal
0055   ped0stats = new MbdRunningStats(100);  // use the last 100 events for running pedestal
0056   name = "hPed0_";
0057   name += _ch;
0058   hPed0 = new TH1F(name, name, 3000, -0.5, 2999.5);
0059   // hPed0 = new TH1F(name,name,10000,1,0); // automatically determine the range
0060   name = "hPedEvt_";
0061   name += _ch;
0062   hPedEvt = new TH1F(name, name, 3000, -0.5, 2999.5);
0063 
0064   SetTemplateSize(900, 1000, -10., 20.);
0065   // SetTemplateSize(300,300,0.,15.);
0066 
0067   // Set pedestal function
0068   ped_fcn = new TF1("ped_fcn","[0]",0,2);
0069   ped_fcn->SetLineColor(3);
0070 
0071   // Set tail function
0072   ped_tail = new TF1("ped_tail","[0]+[1]*exp(-[2]*x)",0,2);
0073   ped_tail->SetLineColor(2);
0074 
0075   // uncomment this to write out waveforms from events that have pileup from prev. crossing
0076   /*
0077   name = "mbdsig"; name += _ch; name += ".txt";
0078   _pileupfile = new ofstream(name);
0079   */
0080 }
0081 
0082 void MbdSig::SetTemplateSize(const Int_t nptsx, const Int_t nptsy, const Double_t begt, const Double_t endt)
0083 {
0084   template_npointsx = nptsx;
0085   template_npointsy = nptsy;
0086   template_begintime = begt;
0087   template_endtime = endt;
0088 
0089   template_y.resize(template_npointsx);
0090   template_yrms.resize(template_npointsx);
0091 
0092   Double_t xbinwid = (template_endtime - template_begintime) / (template_npointsx - 1);
0093   Double_t ybinwid = (1.1 + 0.1) / template_npointsy;  // yscale... should we vary this?
0094   
0095   
0096   delete h2Template;
0097   
0098   delete h2Residuals;
0099   
0100   TString name = "h2Template";
0101   name += _ch;
0102   h2Template = new TH2F(name, name, template_npointsx, template_begintime - (xbinwid / 2.), template_endtime + (xbinwid / 2),
0103                         template_npointsy, -0.1 + (ybinwid / 2.0), 1.1 + (ybinwid / 2.0));
0104 
0105   name = "h2Residuals";
0106   name += _ch;
0107   h2Residuals = new TH2F(name, name, template_npointsx, template_begintime - (xbinwid / 2.), template_endtime + (xbinwid / 2),
0108                          80, -20, 20);
0109 
0110   /*
0111   int nbins[] = { template_npointsx, nbinsy };
0112   Double_t lowrange[] = { template_begintime-xbinwid/2.0, -0.1+ybinwid/2.0 };
0113   Double_t highrange[] = { template_endtime+xbinwid/2.0, 1.1+ybinwid/2.0 };
0114   h2Template = new THnSparseF(name,name,2,nbins,lowrange,highrange);
0115   */
0116   // h2Template->cd( gDirectory );
0117 }
0118 
0119 void MbdSig::SetMinMaxFitTime(const Double_t mintime, const Double_t maxtime)
0120 {
0121   fit_min_time = mintime;
0122   fit_max_time = maxtime;
0123 }
0124 
0125 MbdSig::~MbdSig()
0126 {
0127   if ( _pileupfile )
0128   {
0129     _pileupfile->close();
0130   }
0131   delete hRawPulse;
0132   delete hSubPulse;
0133   delete gRawPulse;
0134   delete gSubPulse;
0135   delete hPed0;
0136   delete ped0stats;
0137   // h2Template->Write();
0138   delete h2Template;
0139   delete h2Residuals;
0140   delete hAmpl;
0141   delete hTime;
0142   delete template_fcn;
0143   delete ped_fcn;
0144   delete ped_tail;
0145 }
0146 
0147 void MbdSig::SetEventPed0PreSamp(const Int_t presample, const Int_t nsamps, const int max_samp)
0148 {
0149   ped_presamp = presample;
0150   ped_presamp_nsamps = nsamps;
0151   if ( ped_presamp_nsamps < 0 && max_samp > 0 )
0152   {
0153     ped_presamp_nsamps = max_samp - presample + 1;
0154   }
0155   ped_presamp_maxsamp = max_samp; // max of waveform
0156 }
0157 
0158 void MbdSig::SetCalib(MbdCalib *m)
0159 {
0160   _mbdcal = m;
0161   _pileup_p0 = m->get_pileup(_ch,0);
0162   _pileup_p1 = m->get_pileup(_ch,1);
0163   _pileup_p2 = m->get_pileup(_ch,2);
0164 }
0165 
0166 // This sets y, and x to sample number (starts at 0)
0167 void MbdSig::SetY(const Float_t* y, const int invert)
0168 {
0169   if (hRawPulse == nullptr)
0170   {
0171     Init();
0172   }
0173 
0174   hpulse->Reset();
0175   f_ampl = -9999.;
0176   f_time = -9999.;
0177 
0178   for (int isamp = 0; isamp < _nsamples; isamp++)
0179   {
0180     hRawPulse->SetBinContent(isamp + 1, y[isamp]);
0181     gRawPulse->SetPoint(isamp, Double_t(isamp), y[isamp]);
0182   }
0183 
0184   // Apply pedestal
0185   if (use_ped0 != 0 || minped0samp >= 0 || minped0x != maxped0x || ped_presamp != 0)
0186   {
0187     // std::cout << "sub" << std::endl;
0188  
0189     int ispileup = 0; // whether pileup event or not (from prev crossing)
0190                       // note pileup correction not implemented for CalcEventPed0(),
0191                       // but that is not used
0192 
0193     if (minped0samp >= 0)
0194     {
0195       CalcEventPed0(minped0samp, maxped0samp);
0196     }
0197     else if (minped0x != maxped0x)
0198     {
0199       CalcEventPed0(minped0x, maxped0x);
0200     }
0201     else if (ped_presamp != 0)
0202     {
0203       ispileup = CalcEventPed0_PreSamp(ped_presamp, ped_presamp_nsamps);
0204     }
0205 
0206     for (int isamp = 0; isamp < _nsamples; isamp++)
0207     {
0208       hSubPulse->SetBinContent(isamp + 1, invert * (y[isamp] - ped0));
0209       hSubPulse->SetBinError(isamp + 1, ped0rms);
0210       gSubPulse->SetPoint(isamp, (Double_t) isamp, invert * (y[isamp] - ped0));
0211       gSubPulse->SetPointError(isamp, 0., ped0rms);
0212     }
0213 
0214     if ( ispileup==1 && !std::isnan(_pileup_p0) )
0215     {
0216       Remove_Pileup();
0217     }
0218   }
0219 
0220   _evt_counter++;
0221 }
0222 
0223 void MbdSig::SetXY(const Float_t* x, const Float_t* y, const int invert)
0224 {
0225   //_verbose = 100;
0226   if (hRawPulse == nullptr)
0227   {
0228     Init();
0229   }
0230 
0231   hRawPulse->Reset();
0232   hSubPulse->Reset();
0233   _status = 0;
0234 
0235   f_ampl = -9999.;
0236   f_time = -9999.;
0237 
0238   // std::cout << "_nsamples " << _nsamples << std::endl;
0239   // std::cout << "use_ped0 " << use_ped0 << "\t" << ped0 << std::endl;
0240 
0241   for (int isamp = 0; isamp < _nsamples; isamp++)
0242   {
0243     // std::cout << "aaa\t" << isamp << "\t" << x[isamp] << "\t" << y[isamp] << std::endl;
0244     hRawPulse->SetBinContent(isamp + 1, y[isamp]);
0245     gRawPulse->SetPoint(isamp, x[isamp], y[isamp]);
0246     gRawPulse->SetPointError(isamp, 0, 4.0);
0247   }
0248   if ( _verbose && _ch==9 )
0249   {
0250     gRawPulse->Draw("ap");
0251     gRawPulse->GetHistogram()->SetTitle(gRawPulse->GetName());
0252     gPad->SetGridy(1);
0253     PadUpdate();
0254   }
0255 
0256   if (use_ped0 != 0 || minped0samp >= 0 || minped0x != maxped0x || ped_presamp != 0)
0257   {
0258     int ispileup = 0; // whether pileup event or not (from prev crossing)
0259                       // note pileup correction not implemented for CalcEventPed0(),
0260                       // but that is not used
0261 
0262     if (minped0samp >= 0)
0263     {
0264       CalcEventPed0(minped0samp, maxped0samp);
0265     }
0266     else if (minped0x != maxped0x)
0267     {
0268       CalcEventPed0(minped0x, maxped0x);
0269     }
0270     else if (ped_presamp != 0)
0271     {
0272       ispileup = CalcEventPed0_PreSamp(ped_presamp, ped_presamp_nsamps);
0273     }
0274 
0275     for (int isamp = 0; isamp < _nsamples; isamp++)
0276     {
0277       if ( _verbose && isamp==(_nsamples-1) )
0278       {
0279         std::cout << "bbb ch " << _ch << "\t" << isamp << "\t" << x[isamp] << "\t" << invert*(y[isamp]-ped0) << std::endl;
0280       }
0281       hSubPulse->SetBinContent(isamp + 1, invert * (y[isamp] - ped0));
0282       hSubPulse->SetBinError(isamp + 1, ped0rms);
0283       gSubPulse->SetPoint(isamp, x[isamp], invert * (y[isamp] - ped0));
0284       gSubPulse->SetPointError(isamp, 0., ped0rms);
0285     }
0286 
0287     if ( ispileup==1 )
0288     {
0289       Remove_Pileup();
0290     }
0291 
0292     if ( _verbose && _ch==9 )
0293     {
0294       std::cout << "SetXY: ch " << _ch << std::endl;
0295       gSubPulse->Print("ALL");
0296     }
0297   }
0298 
0299   _evt_counter++;
0300   _verbose = 0;
0301 }
0302 
0303 void MbdSig::Remove_Pileup()
0304 {
0305   //_verbose = 100;
0306   _verbose = 0;
0307 
0308   if ( (_ch/8)%2 == 0 )   // time ch
0309   {
0310     float offset = _pileup_p0*gSubPulse->GetPointY(0);
0311 
0312     for (int isamp = 0; isamp < _nsamples; isamp++)
0313     {
0314       double x = gSubPulse->GetPointX(isamp);
0315       double y = gSubPulse->GetPointY(isamp);
0316 
0317       hSubPulse->SetBinContent( isamp + 1, y - offset );
0318       gSubPulse->SetPoint( isamp, x, y - offset );
0319     }
0320   }
0321   else
0322   {
0323     if ( fit_pileup == nullptr )
0324     {
0325       TString name = "fit_pileup"; name += _ch;
0326       fit_pileup = new TF1(name,"gaus",-0.1,4.1);
0327       fit_pileup->SetLineColor(2);
0328     }
0329 
0330     fit_pileup->SetRange(-0.1,4.1);
0331     fit_pileup->SetParameters( _pileup_p0*gSubPulse->GetPointY(0), _pileup_p1, _pileup_p2 );
0332     
0333     if ( _verbose )
0334     {
0335       gSubPulse->Fit( fit_pileup, "R" );
0336     }
0337     else
0338     {
0339       gSubPulse->Fit( fit_pileup, "RNQ" );
0340     }
0341 
0342     for (int isamp = 0; isamp < _nsamples; isamp++)
0343     {
0344 
0345       double bkg = fit_pileup->Eval(isamp);
0346 
0347       double x = gSubPulse->GetPointX(isamp);
0348       double y = gSubPulse->GetPointY(isamp);
0349 
0350       float newval = static_cast<float>( y - bkg );
0351 
0352       hSubPulse->SetBinContent( isamp + 1, newval );
0353       gSubPulse->SetPoint( isamp, x, newval );
0354     }
0355   }
0356 
0357   if ( _verbose )
0358   {
0359     gSubPulse->Draw("ap");
0360     PadUpdate();
0361   }
0362 
0363   _verbose = 0;
0364 }
0365 
0366 Double_t MbdSig::GetSplineAmpl()
0367 {
0368   if (gSubPulse == nullptr)
0369   {
0370     std::cout << "gsub bad " << (uint64_t) gSubPulse << std::endl;
0371     return 0.;
0372   }
0373 
0374   TSpline3 s3("s3", gSubPulse);
0375 
0376   // First find maximum, to rescale
0377   f_ampl = -999999.;
0378   double step_size = 0.01;
0379   // std::cout << "step size " << step_size << std::endl;
0380   for (double ix = 0; ix < _nsamples; ix += step_size)
0381   {
0382     Double_t val = s3.Eval(ix);
0383     f_ampl = std::max(val, f_ampl);
0384   }
0385 
0386   return f_ampl;
0387 }
0388 
0389 void MbdSig::WritePedHist()
0390 {
0391   hPed0->Write();
0392 }
0393 
0394 void MbdSig::FillPed0(const Int_t sampmin, const Int_t sampmax)
0395 {
0396   Double_t x;
0397   Double_t y;
0398   for (int isamp = sampmin; isamp <= sampmax; isamp++)
0399   {
0400     gRawPulse->GetPoint(isamp, x, y);
0401     // gRawPulse->Print("all");
0402     hPed0->Fill(y);
0403 
0404     /*
0405     // chiu taken out
0406     ped0stats->Push( y );
0407     ped0 = ped0stats->Mean();
0408     ped0rms = ped0stats->RMS();
0409     */
0410 
0411     // std::cout << "ped0 " << _ch << " " << n << "\t" << ped0 << std::endl;
0412     // std::cout << "ped0 " << _ch << "\t" << ped0 << std::endl;
0413   }
0414 }
0415 
0416 void MbdSig::FillPed0(const Double_t begin, const Double_t end)
0417 {
0418   Double_t x;
0419   Double_t y;
0420   Int_t n = gRawPulse->GetN();
0421   for (int isamp = 0; isamp < n; isamp++)
0422   {
0423     gRawPulse->GetPoint(isamp, x, y);
0424     if (x >= begin && x <= end)
0425     {
0426       hPed0->Fill(y);
0427 
0428       /*
0429          ped0stats->Push( y );
0430          ped0 = ped0stats->Mean();
0431          ped0rms = ped0stats->RMS();
0432          */
0433 
0434       // std::cout << "ped0 " << _ch << " " << n << "\t" << x << "\t" << y << std::endl;
0435     }
0436 
0437     // quit if we are past the ped region
0438     if (x > end)
0439     {
0440       break;
0441     }
0442   }
0443 }
0444 
0445 void MbdSig::SetPed0(const Double_t mean, const Double_t rms)
0446 {
0447   ped0 = mean;
0448   ped0rms = rms;
0449   use_ped0 = 1;
0450   hpulse = hSubPulse;
0451   gpulse = gSubPulse;
0452   // if ( _ch==8 ) std::cout << "_ch " << _ch << " Ped = " << ped0 << std::endl;
0453 }
0454 
0455 // Get Event by Event Ped0 if requested
0456 void MbdSig::CalcEventPed0(const Int_t minpedsamp, const Int_t maxpedsamp)
0457 {
0458   // if (_ch==8) std::cout << "In MbdSig::CalcEventPed0(int,int)" << std::endl;
0459   hPedEvt->Reset();
0460 
0461   Double_t x;
0462   Double_t y;
0463   for (int isamp = minpedsamp; isamp <= maxpedsamp; isamp++)
0464   {
0465     gRawPulse->GetPoint(isamp, x, y);
0466 
0467     hPed0->Fill(y);
0468     hPedEvt->Fill(y);
0469     // ped0stats->Push( y );
0470     // if ( _ch==8 ) std::cout << "ped0stats " << isamp << "\t" << y << std::endl;
0471   }
0472 
0473   // use straight mean for pedestal
0474   // Could consider using fit to hPed0 to remove outliers
0475   float mean = hPedEvt->GetMean();
0476   float rms = hPedEvt->GetRMS();
0477 
0478   SetPed0(mean, rms);
0479   // if (_ch==8) std::cout << "ped0stats mean, rms " << mean << "\t" << rms << std::endl;
0480 }
0481 
0482 // Get Event by Event Ped0 if requested
0483 void MbdSig::CalcEventPed0(const Double_t minpedx, const Double_t maxpedx)
0484 {
0485   hPedEvt->Reset();
0486 
0487   Double_t x;
0488   Double_t y;
0489   Int_t n = gRawPulse->GetN();
0490 
0491   for (int isamp = 0; isamp < n; isamp++)
0492   {
0493     gRawPulse->GetPoint(isamp, x, y);
0494 
0495     if (x >= minpedx && x <= maxpedx)
0496     {
0497       hPed0->Fill(y);
0498       hPedEvt->Fill(y);
0499       // ped0stats->Push( y );
0500     }
0501   }
0502 
0503   // use straight mean for pedestal
0504   // Could consider using fit to hPed0 to remove outliers
0505   SetPed0(hPedEvt->GetMean(), hPedEvt->GetRMS());
0506 }
0507 
0508 // Get Event by Event Ped0, num samples before peak
0509 // presample is number of samples before peak, nsamps is how many samples
0510 // If difficult to get pedestal, use running mean from previous 100 events
0511 // If a prev event pileup is detected, return 1, otherwise, return 0
0512 int MbdSig::CalcEventPed0_PreSamp(const int presample, const int nsamps)
0513 {
0514   //std::cout << PHWHERE << std::endl;  //chiu
0515   //_verbose = 100;
0516   //ped0stats->Clear();
0517 
0518   int status = 0; // assume no pileup
0519 
0520   // Int_t n = gRawPulse->GetN();
0521   // Int_t max = gRawPulse->GetHistogram()->GetMaximumBin();
0522   Long64_t max = ped_presamp_maxsamp;
0523 
0524   // actual max from event
0525   Long64_t actual_max = TMath::LocMax(gRawPulse->GetN(), gRawPulse->GetY());
0526 
0527   if ( ped_presamp_maxsamp == -1 ) // if there is no maxsamp set, use the max found in this event
0528   {
0529     max = actual_max;
0530   }
0531 
0532   Int_t minsamp = max - presample - nsamps + 1;
0533   Int_t maxsamp = max - presample;
0534   //std::cout << "CalcEventPed0_PreSamp: ch ctr max " << _ch << "\t" << _evt_counter << "\t" << max << std::endl;
0535 
0536   if (minsamp < 0)
0537   {
0538     static int counter = 0;
0539     if ( counter<1 )
0540     {
0541       std::cout << PHWHERE << " ped minsamp " << minsamp << "\t" << max << "\t" << presample << "\t" << nsamps << std::endl;
0542       counter++;
0543     }
0544     minsamp = 0;
0545   }
0546   if (maxsamp < 0)
0547   {
0548     static int counter = 0;
0549     if ( counter<1 )
0550     {
0551       std::cout << PHWHERE << " ped maxsamp " << maxsamp << "\t" << max << "\t" << presample << std::endl;
0552       counter++;
0553     }
0554     maxsamp = minsamp;
0555   }
0556 
0557   if ( _verbose )
0558   {
0559     std::cout << "CalcEventPed0_Presamp(), ch " << _ch << std::endl;
0560   }
0561 
0562   double mean = ped0stats->Mean();
0563   //double rms = ped0stats->RMS();
0564   double rms = _mbdcal->get_pedrms(_ch);
0565   if ( std::isnan(rms) )
0566   {
0567     // ped calib doesn't exist, set to average
0568     rms = 5.0;
0569   }
0570 
0571   ped_fcn->SetRange(minsamp-0.1,maxsamp+0.1);
0572   ped_fcn->SetParameter(0,1500.);
0573 
0574   if ( gRawPulse->GetN()==0 )//chiu
0575   {
0576     std::cout << PHWHERE << " gRawPulse 0" << std::endl;
0577   }
0578 
0579   if ( _verbose )
0580   {
0581     gRawPulse->Fit( ped_fcn, "RQ" );
0582 
0583     double chi2ndf = ped_fcn->GetChisquare()/ped_fcn->GetNDF();
0584     if ( chi2ndf > 4.0 )
0585     {
0586       gRawPulse->Draw("ap");
0587       ped_fcn->Draw("same");
0588       PadUpdate();
0589     }
0590   }
0591   else
0592   {
0593     //std::cout << PHWHERE << std::endl;
0594     gRawPulse->Fit( ped_fcn, "RNQ" );
0595 
0596     double chi2ndf = ped_fcn->GetChisquare()/ped_fcn->GetNDF();
0597     if ( _pileupfile != nullptr && chi2ndf > 4.0 )
0598     {
0599       *_pileupfile << "ped " << _ch << " mean " << mean << "\t";
0600       for ( int i=0; i<gRawPulse->GetN(); i++)
0601       {
0602         *_pileupfile << std::setw(6) << gRawPulse->GetPointY(i);
0603       }
0604       *_pileupfile << std::endl;
0605     }
0606   }
0607 
0608   double chi2 = ped_fcn->GetChisquare();
0609   double ndf = ped_fcn->GetNDF();
0610 
0611   if ( chi2/ndf < 4.0 )
0612   {
0613     mean = ped_fcn->GetParameter(0);
0614 
0615     Double_t x;
0616     Double_t y;
0617 
0618     for (int isamp = minsamp; isamp <= maxsamp; isamp++)
0619     {
0620       gRawPulse->GetPoint(isamp, x, y);
0621 
0622       // exclude outliers
0623       if ( fabs(y-mean) < 4.0*rms )
0624       {
0625         ped0stats->Push( y );
0626       }
0627 
0628       if ( _verbose )
0629       {
0630         std::cout << "CalcEventPed0_PreSamp: ped0stats " << _ch << "\t" << _evt_counter << "\t" 
0631           << "isamp " << isamp << "\t" << x << "\t" << y << std::endl;
0632       }
0633     }
0634   }
0635   else
0636   {
0637     // ped fit was bad, we have signal contamination in ped region
0638     // or other thing going on
0639     status = 1;
0640 
0641     if ( ped0stats->Size() < ped0stats->MaxNum() && !std::isnan(_mbdcal->get_ped(_ch)) ) // use pre-calib for early events
0642     {
0643       mean = _mbdcal->get_ped(_ch);
0644     }
0645     else  // use running mean, once enough stats are accumulated, or if no pre-calib pedestals
0646     {
0647       mean = ped0stats->Mean();
0648 
0649       if ( _verbose )
0650       {
0651         gRawPulse->Draw("ap");
0652         PadUpdate();
0653 
0654         if ( _evt_counter<100 )
0655         {
0656           double ped0statsrms = ped0stats->RMS();
0657           std::cout << "aaa ch " << _ch << "\t" << _evt_counter << "\t" << mean << "\t" << ped0statsrms << std::endl;
0658         }
0659         std::string junk;
0660         std::cin >> junk;
0661       }
0662     }
0663 
0664     // use straight mean for pedestal
0665     // Could consider using fit to hPed0 to remove outliers
0666     //rms = ped0stats->RMS();
0667     //Double_t mean = hPed0->GetMean();
0668     //Double_t rms = hPed0->GetRMS();
0669   }
0670 
0671   SetPed0(mean, rms);
0672 
0673   if (_verbose > 0 && _evt_counter < 10)
0674   {
0675     std::cout << "CalcEventPed0_PreSamp: ped0stats " << _ch << "\t" << _evt_counter << "\t" << mean << "\t" << rms << std::endl;
0676     std::cout << "CalcEventPed0_PreSamp: ped0stats " << ped0stats->RMS() << std::endl;
0677   }
0678 
0679   _verbose = 0;
0680 
0681   return status;
0682 }
0683 
0684 Double_t MbdSig::LeadingEdge(const Double_t threshold)
0685 {
0686   // Find first point above threshold
0687   // We also make sure the next point is above threshold
0688   // to get rid of a high fluctuation
0689   int n = gSubPulse->GetN();
0690   Double_t* x = gSubPulse->GetX();
0691   Double_t* y = gSubPulse->GetY();
0692 
0693   int sample = -1;
0694   for (int isamp = 0; isamp < n; isamp++)
0695   {
0696     if (y[isamp] > threshold)
0697     {
0698       if (isamp == (n - 1) || y[isamp + 1] > threshold)
0699       {
0700         sample = isamp;
0701         break;
0702       }
0703     }
0704   }
0705   if (sample < 1)
0706   {
0707     return -9999.;  // no signal above threshold
0708   }
0709 
0710   // Linear Interpolation of start time
0711   Double_t dx = x[sample] - x[sample - 1];
0712   Double_t dy = y[sample] - y[sample - 1];
0713   Double_t dt1 = y[sample] - threshold;
0714 
0715   Double_t t0 = x[sample] - (dt1 * (dx / dy));
0716 
0717   return t0;
0718 }
0719 
0720 Double_t MbdSig::dCFD(const Double_t fraction_threshold)
0721 {
0722   // Find first point above threshold
0723   // We also make sure the next point is above threshold
0724   // to get rid of a high fluctuation
0725   int n = gSubPulse->GetN();
0726   Double_t* x = gSubPulse->GetX();
0727   Double_t* y = gSubPulse->GetY();
0728 
0729   // Get max amplitude
0730   Double_t ymax = TMath::MaxElement(n, y);
0731   if (f_ampl == -9999.)
0732   {
0733     f_ampl = ymax;
0734   }
0735 
0736   Double_t threshold = fraction_threshold * ymax;  // get fraction of amplitude
0737   // std::cout << "threshold = " << threshold << "\tymax = " << ymax <<std::endl;
0738 
0739   int sample = -1;
0740   for (int isamp = 0; isamp < n; isamp++)
0741   {
0742     if (y[isamp] > threshold)
0743     {
0744       if (isamp == (n - 1) || y[isamp + 1] > threshold)
0745       {
0746         sample = isamp;
0747         break;
0748       }
0749     }
0750   }
0751   if (sample < 1)
0752   {
0753     return -9999.;  // no signal above threshold
0754   }
0755 
0756   // Linear Interpolation of start time
0757   Double_t dx = x[sample] - x[sample - 1];
0758   Double_t dy = y[sample] - y[sample - 1];
0759   Double_t dt1 = y[sample] - threshold;
0760 
0761   Double_t t0 = x[sample] - (dt1 * (dx / dy));
0762 
0763   return t0;
0764 }
0765 
0766 Double_t MbdSig::MBDTDC(const Int_t max_samp)
0767 {
0768   // Get the amplitude of a fixed sample (max_samp) to get time
0769   // Used in MBD Time Channels
0770   Double_t* y = gSubPulse->GetY();
0771 
0772   if (y == nullptr)
0773   {
0774     std::cout << "ERROR y == 0" << std::endl;
0775     return std::numeric_limits<Double_t>::quiet_NaN();
0776   }
0777 
0778   f_time = y[max_samp];
0779 
0780   if ( y[2]>100. )
0781   {
0782     f_time = 0.;
0783   }
0784 
0785   /*
0786   if ( _ch==128 )
0787   {
0788     std::cout << "msig\t" << _ch << "\t" << f_time << "\t" << f_ampl << std::endl;
0789     gSubPulse->Print("ALL");
0790   }
0791   */
0792 
0793   return f_time;
0794 }
0795 
0796 Double_t MbdSig::Integral(const Double_t xmin, const Double_t xmax)
0797 {
0798   Int_t n = gSubPulse->GetN();
0799   Double_t* x = gSubPulse->GetX();
0800   Double_t* y = gSubPulse->GetY();
0801 
0802   f_integral = 0.;
0803   for (int ix = 0; ix < n; ix++)
0804   {
0805     if (x[ix] >= xmin && x[ix] <= xmax)
0806     {
0807       // Get dx
0808       Double_t dx = (x[ix + 1] - x[ix - 1]) / 2.0;
0809       f_integral += (y[ix] * dx);
0810     }
0811   }
0812 
0813   return f_integral;
0814 }
0815 
0816 void MbdSig::LocMax(Double_t& x_at_max, Double_t& ymax, Double_t xminrange, Double_t xmaxrange)
0817 {
0818   _verbose = 0;
0819   if ( _verbose && _ch==250 )
0820   {
0821     gSubPulse->Draw("ap");
0822     gPad->Modified();
0823     gPad->Update();
0824   }
0825 
0826   // Find index of maximum peak
0827   Int_t n = gSubPulse->GetN();
0828   Double_t* x = gSubPulse->GetX();
0829   Double_t* y = gSubPulse->GetY();
0830 
0831   // if flipped or equal, we search the whole range
0832   if (xmaxrange <= xminrange)
0833   {
0834     xminrange = -DBL_MAX;
0835     xmaxrange = DBL_MAX;
0836   }
0837 
0838   x_at_max = -DBL_MAX;
0839   ymax = -DBL_MAX;
0840 
0841   for (int i = 0; i < n; i++)
0842   {
0843     // Skip if out of range
0844     if (x[i] < xminrange)
0845     {
0846       continue;
0847     }
0848     if (x[i] > xmaxrange)
0849     {
0850       break;
0851     }
0852 
0853     if (y[i] > ymax)
0854     {
0855       ymax = y[i];
0856       x_at_max = x[i];
0857     }
0858   }
0859 
0860   if ( _verbose )
0861   {
0862     std::string junk;
0863     std::cin >> junk;
0864     std::cout << "MbdSig::LocMax, " << x_at_max << "\t" << ymax << std::endl;
0865     _verbose = 0;
0866   }
0867 }
0868 
0869 void MbdSig::LocMin(Double_t& x_at_min, Double_t& ymin, Double_t xminrange, Double_t xmaxrange)
0870 {
0871   // Find index of minimum peak (for neg signals)
0872   Int_t n = gSubPulse->GetN();
0873   Double_t* x = gSubPulse->GetX();
0874   Double_t* y = gSubPulse->GetY();
0875 
0876   // if flipped or equal, we search the whole range
0877   if (xmaxrange <= xminrange)
0878   {
0879     xminrange = -DBL_MAX;
0880     xmaxrange = DBL_MAX;
0881   }
0882 
0883   ymin = DBL_MAX;
0884 
0885   for (int i = 0; i < n; i++)
0886   {
0887     // Skip if out of range
0888     if (x[i] < xminrange)
0889     {
0890       continue;
0891     }
0892     if (x[i] > xmaxrange)
0893     {
0894       break;
0895     }
0896 
0897     if (y[i] < ymin)
0898     {
0899       ymin = y[i];
0900       x_at_min = x[i];
0901     }
0902   }
0903 
0904   // old way of getting locmax
0905   // int locmax = TMath::LocMin(n,y);
0906 }
0907 
0908 void MbdSig::Print()
0909 {
0910   Double_t x;
0911   Double_t y;
0912   std::cout << "CH " << _ch << std::endl;
0913   for (int isamp = 0; isamp < _nsamples; isamp++)
0914   {
0915     gpulse->GetPoint(isamp, x, y);
0916     std::cout << isamp << "\t" << x << "\t" << y << std::endl;
0917   }
0918 }
0919 
0920 void MbdSig::DrawWaveform()
0921 {
0922   int orig_verbose = _verbose;
0923   if ( _verbose <= 5 )
0924   {
0925     _verbose = 6;
0926   }
0927 
0928   gSubPulse->Draw("ap");
0929   gSubPulse->GetHistogram()->SetTitle(gSubPulse->GetName());
0930   gPad->SetGridy(1);
0931   if ( template_fcn!=nullptr ) template_fcn->Draw("same");  // should check if fit was made
0932   PadUpdate();
0933 
0934   _verbose = orig_verbose;
0935 }
0936 
0937 void MbdSig::PadUpdate() const
0938 {
0939   // Make sure TCanvas is created externally!
0940   std::cout << PHWHERE << " PadUpdate\t_verbose = " << _verbose << std::endl;
0941   if ( _verbose>5 )
0942   {
0943     gPad->Modified();
0944     gPad->Update();
0945     std::cout << _ch << " ? ";
0946     if ( _verbose>10 )
0947     {
0948       std::string junk;
0949       std::cin >> junk;
0950 
0951       if (junk[0] == 'w' || junk[0] == 's')
0952       {
0953         TString name = "ch";
0954         name += _ch;
0955         name += ".png";
0956         gPad->SaveAs(name);
0957       }
0958     }
0959   }
0960 }
0961 
0962 Double_t MbdSig::TemplateFcn(const Double_t* x, const Double_t* par)
0963 {
0964   // par[0] is the amplitude (relative to the spline amplitude)
0965   // par[1] is the start time (in sample number)
0966   // x[0] units are in sample number
0967   Double_t xx = x[0] - par[1];
0968   Double_t f = 0.;
0969 
0970   int verbose = 0;
0971   if ( verbose )
0972   {
0973     std::cout << PHWHERE << " " << _ch << " x par0 par1 " << x[0] << "\t" << par[0] << "\t" << par[1] << std::endl;
0974   }
0975 
0976 
0977   // When fit is out of limits of good part of spline, ignore fit
0978   if (xx < template_begintime || xx > template_endtime || std::isnan(xx) )
0979   {
0980     TF1::RejectPoint();
0981 
0982     if ( std::isnan(xx) )
0983     {
0984       if (verbose > 0)
0985       {
0986     std::cout << PHWHERE << " " << _evt_counter << ", " << _ch
0987           << " ERROR x par0 par1 " << x[0] << "\t" << par[0] << "\t" << par[1] << std::endl;
0988     if ( x[0] == 0. )
0989     {
0990       gSubPulse->Print("ALL");
0991     }
0992       }
0993       return 0.;
0994     }
0995 
0996     if (xx < template_begintime)
0997     {
0998       // Double_t x0,y0;
0999       Double_t y0 = template_y[0];
1000       return par[0] * y0;
1001     }
1002     if (xx > template_endtime)
1003     {
1004       // Double_t x0,y0;
1005       Double_t y0 = template_y[template_npointsx - 1];
1006       return par[0] * y0;
1007     }
1008   }
1009 
1010   // Linear Interpolation of template
1011   Double_t x0 = 0.;
1012   Double_t y0 = 0.;
1013   Double_t x1 = 0.;
1014   Double_t y1 = 0.;
1015 
1016   // find the index in the vector which is closest to xx
1017   Double_t step = (template_endtime - template_begintime) / (template_npointsx - 1);
1018   Double_t index = (xx - template_begintime) / step;
1019   //std::cout << "xxx " << index << "\t" << xx << "\t" << template_begintime << "\t" << template_endtime << "\t" << step << std::endl;
1020 
1021   int ilow = TMath::FloorNint(index);
1022   int ihigh = TMath::CeilNint(index);
1023   if (ilow < 0 || ihigh >= template_npointsx)
1024   {
1025     if (verbose > 0)
1026     {
1027       std::cout << "ERROR, ilow ihigh " << ilow << "\t" << ihigh << std::endl;
1028       std::cout << " " << xx << " " << x[0] << " " << par[1] << std::endl;
1029     }
1030 
1031     if (ilow < 0)
1032     {
1033       ilow = 0;
1034     }
1035     else if (ihigh >= template_npointsx)
1036     {
1037       ihigh = template_npointsx - 1;
1038     }
1039   }
1040 
1041   if (ilow == ihigh)
1042   {
1043     f = par[0] * template_y[ilow];
1044   }
1045   else
1046   {
1047     x0 = template_begintime + ilow * step;
1048     y0 = template_y[ilow];
1049     x1 = template_begintime + ihigh * step;
1050     y1 = template_y[ihigh];
1051     f = par[0] * (y0 + ((y1 - y0) / (x1 - x0)) * (xx - x0));  // linear interpolation
1052   }
1053 
1054   // reject points with very bad rms in shape
1055   if (template_yrms[ilow] >= 1.0 || template_yrms[ihigh] >= 1.0)
1056   {
1057     TF1::RejectPoint();
1058     // return f;
1059   }
1060 
1061   // Reject points where ADC saturates
1062   int samp_point = static_cast<int>(x[0]);
1063   Double_t temp_x;
1064   Double_t temp_y;
1065   gRawPulse->GetPoint(samp_point, temp_x, temp_y);
1066   if (temp_y > 16370)
1067   {
1068     //if ( _ch==185 ) std::cout << "ADCSATURATED " << _ch << "\t" << samp_point << "\t" << temp_x << "\t" << temp_y << std::endl;
1069     TF1::RejectPoint();
1070   }
1071 
1072   //verbose = 0;
1073   return f;
1074 }
1075 
1076 // sampmax>0 means fit to the peak near sampmax
1077 int MbdSig::FitTemplate( const Int_t sampmax )
1078 {
1079   //std::cout << PHWHERE << std::endl;  //chiu
1080   //_verbose = 100; // uncomment to see fits
1081   //_verbose = 12;        // don't see pedestal fits
1082 
1083   // Check if channel is empty
1084   if (gSubPulse->GetN() == 0)
1085   {
1086     f_ampl = 0.;
1087     f_time = std::numeric_limits<Float_t>::quiet_NaN();
1088     std::cout << "ERROR, gSubPulse empty" << std::endl;
1089     return 1;
1090   }
1091 
1092   // Determine if channel is saturated
1093   Double_t *rawsamps = gRawPulse->GetY();
1094   Int_t nrawsamps = gRawPulse->GetN();
1095   int nsaturated = 0;
1096   for (int ipt=0; ipt<nrawsamps; ipt++)
1097   {
1098     if ( rawsamps[ipt] > 16370. )
1099     {
1100       nsaturated++;
1101     }
1102   }
1103   /*
1104   if ( nsaturated>2 && _ch==185 )
1105   {
1106     _verbose = 12;
1107   }
1108   */
1109 
1110 
1111   if (_verbose > 0)
1112   {
1113     std::cout << "Fitting ch " << _ch << std::endl;
1114   }
1115 
1116   // Get x and y of maximum
1117   Double_t x_at_max{-1.};
1118   Double_t ymax{0.};
1119   if ( sampmax>=0 )
1120   {
1121     gSubPulse->GetPoint(sampmax, x_at_max, ymax);
1122     if ( nsaturated<=3 )
1123     {
1124       x_at_max -= 2.0;
1125     }
1126     else
1127     {
1128       x_at_max -= 1.5;
1129       ymax = 16370.+nsaturated*2000.;
1130     }
1131   }
1132   else
1133   {
1134     ymax = TMath::MaxElement( gSubPulse->GetN(), gSubPulse->GetY() );
1135     x_at_max = TMath::LocMax( gSubPulse->GetN(), gSubPulse->GetY() );
1136   }
1137 
1138   // Threshold cut
1139   if ( ymax < 20. )
1140   {
1141     f_ampl = 0.;
1142     f_time = std::numeric_limits<Float_t>::quiet_NaN();
1143     if ( _verbose>20 )
1144     {
1145       // for checking pedestal
1146       std::cout << "skipping, ymax < 20" << std::endl;
1147       gSubPulse->Draw("ap");
1148       gSubPulse->GetHistogram()->SetTitle(gSubPulse->GetName());
1149       gPad->SetGridy(1);
1150       PadUpdate();
1151     }
1152 
1153     _verbose = 0;
1154     return 1;
1155   }
1156 
1157   template_fcn->SetParameters(ymax, x_at_max);
1158   // template_fcn->SetParLimits(1, fit_min_time, fit_max_time);
1159   // template_fcn->SetParLimits(1, 3, 15);
1160   // template_fcn->SetRange(template_min_xrange,template_max_xrange);
1161   if ( nsaturated<=3 )
1162   {
1163     template_fcn->SetRange(0, _nsamples);
1164   }
1165   else
1166   {
1167     template_fcn->SetRange(0, sampmax + nsaturated - 0.5);
1168   }
1169 
1170   if ( gSubPulse->GetN()==0 )//chiu
1171   {
1172     std::cout << PHWHERE << " gSubPulse 0" << std::endl;
1173   }
1174 
1175   if (_verbose == 0)
1176   {
1177     //std::cout << PHWHERE << std::endl;
1178     gSubPulse->Fit(template_fcn, "RNQ");
1179   }
1180   else
1181   {
1182     std::cout << "doing fit1 " << x_at_max << "\t" << ymax << std::endl;
1183     gSubPulse->Fit(template_fcn, "R");
1184     gSubPulse->Draw("ap");
1185     gSubPulse->GetHistogram()->SetTitle(gSubPulse->GetName());
1186     gPad->SetGridy(1);
1187     PadUpdate();
1188     //std::cout << "doing fit2 " << _verbose << std::endl;
1189     //std::cout << "doing fit3 " << _verbose << std::endl;
1190     //gSubPulse->Print("ALL");
1191   }
1192 
1193   // Get fit parameters
1194   f_ampl = template_fcn->GetParameter(0);
1195   f_time = template_fcn->GetParameter(1);
1196   if ( f_time<0. || f_time>_nsamples )
1197   {
1198     f_time = _nsamples*0.5;  // bad fit last time
1199   }
1200 
1201   // refit with new range to exclude after-pulses
1202   template_fcn->SetParameters( f_ampl, f_time );
1203   if ( nsaturated<=3 )
1204   {
1205     template_fcn->SetRange( 0., f_time+4.0 );
1206   }
1207   else
1208   {
1209     template_fcn->SetRange( 0., f_time+nsaturated+0.8 );
1210   }
1211 
1212   if (_verbose == 0)
1213   {
1214     //std::cout << PHWHERE << std::endl;
1215     int fit_status = gSubPulse->Fit(template_fcn, "RNQ");
1216     if ( fit_status<0 )
1217     {
1218       std::cout << PHWHERE << "\t" << fit_status << std::endl;
1219       gSubPulse->Print("ALL");
1220       gSubPulse->Draw("ap");
1221       gSubPulse->Fit(template_fcn, "R");
1222       std::cout << "ampl time before refit " << f_ampl << "\t" << f_time << std::endl;
1223       f_ampl = template_fcn->GetParameter(0);
1224       f_time = template_fcn->GetParameter(1);
1225       std::cout << "ampl time after  refit " << f_ampl << "\t" << f_time << std::endl;
1226       PadUpdate();
1227       std::string junk;
1228       std::cin >> junk;
1229     }
1230   }
1231   else
1232   {
1233     gSubPulse->Fit(template_fcn, "R");
1234     //gSubPulse->Print("ALL");
1235     std::cout << "ampl time before refit " << f_ampl << "\t" << f_time << std::endl;
1236     f_ampl = template_fcn->GetParameter(0);
1237     f_time = template_fcn->GetParameter(1);
1238     std::cout << "ampl time after  refit " << f_ampl << "\t" << f_time << std::endl;
1239   }
1240 
1241   f_ampl = template_fcn->GetParameter(0);
1242   f_time = template_fcn->GetParameter(1);
1243 
1244   //if ( f_time<0 || f_time>30 )
1245   //if ( (_ch==185||_ch==155||_ch==249) && (fabs(f_ampl) > 44000.) )
1246   //double chi2 = template_fcn->GetChisquare();
1247   //double ndf = template_fcn->GetNDF();
1248   //if ( (_ch==185||_ch==155||_ch==249) && (fabs(chi2/ndf) > 100.) && nsaturated > 3)
1249   if (_verbose > 0 && fabs(f_ampl) > 0.)
1250   {
1251     _verbose = 12;
1252     std::cout << "FitTemplate " << _ch << "\t" << f_ampl << "\t" << f_time << std::endl;
1253     std::cout << "            " << template_fcn->GetChisquare()/template_fcn->GetNDF() << std::endl;
1254     gSubPulse->Draw("ap");
1255     gSubPulse->GetHistogram()->SetTitle(gSubPulse->GetName());
1256     gPad->SetGridy(1);
1257     template_fcn->SetLineColor(4);
1258     template_fcn->Draw("same");
1259     PadUpdate();
1260   }
1261 
1262   _verbose = 0;
1263   return 1;
1264 }
1265 
1266 int MbdSig::SetTemplate(const std::vector<float>& shape, const std::vector<float>& sherr)
1267 {
1268   template_y = shape;
1269   template_yrms = sherr;
1270 
1271   if (_verbose)
1272   {
1273     std::cout << "SHAPE " << _ch << "\t" << template_y.size() << std::endl;
1274     for (size_t i = 0; i < template_y.size(); i++)
1275     {
1276       if (i % 10 == 0)
1277       {
1278         std::cout << i << ":\t" << std::endl;
1279       }
1280       std::cout << " " << template_y[i];
1281     }
1282     std::cout << std::endl;
1283   }
1284 
1285   if (template_fcn == nullptr)
1286   {
1287     TString name = "template_fcn";
1288     name += _ch;
1289     // template_fcn = new TF1(name,this,&MbdSig::TemplateFcn,template_min_xrange,template_max_xrange,2,"MbdSig","TemplateFcn");
1290     // template_fcn = new TF1(name,this,&MbdSig::TemplateFcn,-10,20,2,"MbdSig","TemplateFcn");
1291     template_fcn = new TF1(name, this, &MbdSig::TemplateFcn, 0, _nsamples, 2, "MbdSig", "TemplateFcn");
1292     template_fcn->SetParameters(1, 10);
1293     template_fcn->SetParName(0, "ampl");
1294     template_fcn->SetParName(1, "time");
1295     SetTemplateSize(900, 1000, -10., 20.);
1296 
1297     if (_verbose)
1298     {
1299       std::cout << "SHAPE " << _ch << std::endl;
1300       template_fcn->Draw("acp");
1301       gPad->Modified();
1302       gPad->Update();
1303     }
1304   }
1305 
1306   return 1;
1307 }