Back to home page

sPhenix code displayed by LXR

 
 

    


Warning, file /coresoftware/offline/packages/mbd/MbdCalib.cc was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 #include "MbdCalib.h"
0002 
0003 #include "MbdDefs.h"
0004 #include "MbdGeomV1.h"
0005 
0006 // Database Includes
0007 #ifndef ONLINE
0008 #include <ffamodules/CDBInterface.h>
0009 #include <cdbobjects/CDBTTree.h>
0010 #endif
0011 
0012 #include <phool/phool.h>
0013 
0014 #include <cmath>
0015 #include <fstream>
0016 #include <iostream>
0017 #include <regex>
0018 #include <string>
0019 
0020 #include <TString.h>
0021 
0022 
0023 MbdCalib::MbdCalib()
0024 {
0025   Reset();
0026   _mbdgeom = std::make_unique<MbdGeomV1>();
0027 
0028 #ifndef ONLINE
0029   // clang-tidy suggests to put this into the ctor but then ONLINE
0030   // does not compile anymore
0031   // NOLINTNEXTLINE(cppcoreguidelines-prefer-member-initializer)
0032   _rc = recoConsts::instance();
0033   if ( _rc->FlagExist("MBD_TEMPLATEFIT") )
0034   {
0035     do_templatefit = _rc->get_IntFlag("MBD_TEMPLATEFIT");
0036   }
0037   else
0038   {
0039     do_templatefit = 1;
0040   }
0041 #else
0042   do_templatefit = 0;
0043 #endif
0044 
0045 }
0046 
0047 int MbdCalib::Download_All()
0048 {
0049   if ( Verbosity()>0 )
0050   {
0051     std::cout << PHWHERE << " In MbdCalib::Download_All()" << std::endl;
0052   }
0053   _status = 0;
0054 
0055   std::string bbc_caldir;
0056 
0057 #ifdef ONLINE
0058   bbc_caldir = getenv("BBCCALIB");
0059   if (bbc_caldir.size()==0)
0060   {
0061     std::cout << "BBCCALIB environment variable not set" << std::endl;
0062     exit(1);
0063   }
0064 #else
0065   if (Verbosity() > 0)
0066   {
0067     std::cout << "MBD CDB " << _rc->get_StringFlag("CDB_GLOBALTAG") << "\t" << _rc->get_uint64Flag("TIMESTAMP") << std::endl;
0068   }
0069   _cdb = CDBInterface::instance();
0070 
0071   if (_rc->FlagExist("MBD_CALDIR"))
0072   {
0073     bbc_caldir = _rc->get_StringFlag("MBD_CALDIR");
0074     std::cout << "Reading MBD Calibrations from " << bbc_caldir << std::endl;
0075   }
0076 
0077   // if rc flag MBD_CALDIR does not exist, we create it and set it to an empty string
0078   if (!_rc->FlagExist("MBD_CALDIR"))
0079   {
0080     std::string sampmax_url = _cdb->getUrl("MBD_SAMPMAX");
0081     if (Verbosity() > 0)
0082     {
0083       std::cout << "sampmax_url " << sampmax_url << std::endl;
0084     }
0085     Download_SampMax(sampmax_url);
0086 
0087     if ( !_rawdstflag )
0088     {
0089       std::string ped_url = _cdb->getUrl("MBD_PED");
0090       if (Verbosity() > 0)
0091       {
0092         std::cout << "ped_url " << ped_url << std::endl;
0093       }
0094       Download_Ped(ped_url);
0095 
0096     
0097       std::string pileup_url = _cdb->getUrl("MBD_PILEUP");
0098       if (Verbosity() > 0)
0099       {
0100         std::cout << "pileup_url " << pileup_url << std::endl;
0101       }
0102       Download_Pileup(pileup_url);
0103 
0104       if (do_templatefit)
0105       {
0106         std::string shape_url = _cdb->getUrl("MBD_SHAPES");
0107         if (Verbosity() > 0)
0108         {
0109           std::cout << "shape_url " << shape_url << std::endl;
0110         }
0111         Download_Shapes(shape_url);
0112       }
0113     }
0114 
0115     std::string qfit_url = _cdb->getUrl("MBD_QFIT");
0116     if (Verbosity() > 0)
0117     {
0118       std::cout << "qfit_url " << qfit_url << std::endl;
0119     }
0120     Download_Gains(qfit_url);
0121 
0122     std::string tt_t0_url = _cdb->getUrl("MBD_TT_T0");
0123     if ( Verbosity() > 0 )
0124     {
0125       std::cout << "tt_t0_url " << tt_t0_url << std::endl;
0126     }
0127     Download_TTT0(tt_t0_url);
0128 
0129     std::string tq_t0_url = _cdb->getUrl("MBD_TQ_T0");
0130     if (Verbosity() > 0)
0131     {
0132       std::cout << "tq_t0_url " << tq_t0_url << std::endl;
0133     }
0134     Download_TQT0(tq_t0_url);
0135 
0136     if ( !_fitsonly )
0137     {
0138       std::string t0corr_url = _cdb->getUrl("MBD_T0CORR");
0139       if ( Verbosity() > 0 )
0140       {
0141         std::cout << "t0corr_url " << t0corr_url << std::endl;
0142       }
0143       Download_T0Corr(t0corr_url);
0144 
0145       std::string timecorr_url = _cdb->getUrl("MBD_TIMECORR");
0146       if ( Verbosity() > 0 )
0147       {
0148         std::cout << "timecorr_url " << timecorr_url << std::endl;
0149       }
0150       Download_TimeCorr(timecorr_url);
0151 
0152       std::string slew_url = _cdb->getUrl("MBD_SLEWCORR");
0153       if ( Verbosity() > 0 )
0154       {
0155         std::cout << "slew_url " << slew_url << std::endl;
0156       }
0157       Download_SlewCorr(slew_url);
0158     }
0159 
0160     Verbosity(0);
0161   }
0162 #endif
0163 
0164   // download local calibs (text versions)
0165   if ( !bbc_caldir.empty() )
0166   {
0167     std::string sampmax_file = bbc_caldir + "/mbd_sampmax.calib";
0168     Download_SampMax(sampmax_file);
0169 
0170     if ( !_rawdstflag )
0171     {
0172       std::string ped_file = bbc_caldir + "/mbd_ped.calib";
0173       Download_Ped(ped_file);
0174 
0175       std::string pileup_file = bbc_caldir + "/mbd_pileup.calib";
0176       Download_Pileup(pileup_file);
0177 
0178       if (do_templatefit)
0179       {
0180         std::string shape_file = bbc_caldir + "/mbd_shape.calib";
0181         Download_Shapes(shape_file);
0182       }
0183     }
0184 
0185     std::string qfit_file = bbc_caldir + "/mbd_qfit.calib";
0186     Download_Gains(qfit_file);
0187 
0188     std::string tq_t0_file = bbc_caldir + "/mbd_tq_t0.calib";
0189     Download_TQT0(tq_t0_file);
0190 
0191     std::string tt_t0_file = bbc_caldir + "/mbd_tt_t0.calib";
0192     Download_TTT0(tt_t0_file);
0193 
0194     if ( !_fitsonly )
0195     {
0196       std::string t0corr_file = bbc_caldir + "/mbd_t0corr.calib";
0197       Download_T0Corr(t0corr_file);
0198 
0199       std::string tt_tcorr_file = bbc_caldir + "/mbd_timecorr.calib";
0200       Download_TimeCorr(tt_tcorr_file);
0201 
0202       std::string slew_file = bbc_caldir + "/mbd_slewcorr.calib";
0203       Download_SlewCorr(slew_file);
0204     }
0205   }
0206 
0207   return _status;
0208 }
0209 
0210 int MbdCalib::Download_Gains(const std::string& dbase_location)
0211 {
0212   // Reset All Values
0213   _qfit_integ.fill(std::numeric_limits<float>::quiet_NaN());
0214   _qfit_mpv.fill(std::numeric_limits<float>::quiet_NaN());
0215   _qfit_sigma.fill(std::numeric_limits<float>::quiet_NaN());
0216   _qfit_integerr.fill(std::numeric_limits<float>::quiet_NaN());
0217   _qfit_mpverr.fill(std::numeric_limits<float>::quiet_NaN());
0218   _qfit_sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
0219   _qfit_chi2ndf.fill(std::numeric_limits<float>::quiet_NaN());
0220 
0221   if (Verbosity() > 0)
0222   {
0223     std::cout << "Opening " << dbase_location << std::endl;
0224   }
0225   TString dbase_file = dbase_location;
0226 
0227 #ifndef ONLINE
0228   if (dbase_file.EndsWith(".root"))  // read from database
0229   {
0230     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0231     cdbttree->LoadCalibrations();
0232 
0233     for (int ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
0234     {
0235       _qfit_integ[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_integ");
0236       _qfit_mpv[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_mpv");
0237       _qfit_sigma[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_sigma");
0238       _qfit_integerr[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_integerr");
0239       _qfit_mpverr[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_mpverr");
0240       _qfit_sigmaerr[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_sigmaerr");
0241       _qfit_chi2ndf[ipmt] = cdbttree->GetFloatValue(ipmt, "qfit_chi2ndf");
0242       if (Verbosity() > 0)
0243       {
0244         if (ipmt < 5)
0245         {
0246           std::cout << ipmt << "\t" << _qfit_mpv[ipmt] << std::endl;
0247         }
0248       }
0249     }
0250     delete cdbttree;
0251   }
0252 #endif
0253 
0254   if (dbase_file.EndsWith(".calib"))  // read from text file
0255   {
0256     std::ifstream infile(dbase_location);
0257     if (!infile.is_open())
0258     {
0259       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0260       _status = -3;
0261       return _status;
0262     }
0263 
0264     int pmt = -1;
0265     while (infile >> pmt)
0266     {
0267       infile >> _qfit_integ[pmt] >> _qfit_mpv[pmt] >> _qfit_sigma[pmt] >> _qfit_integerr[pmt] >> _qfit_mpverr[pmt] >> _qfit_sigmaerr[pmt] >> _qfit_chi2ndf[pmt];
0268       if (Verbosity() > 0)
0269       {
0270         if (pmt < 5 || pmt >= MbdDefs::MBD_N_PMT - 5)
0271         {
0272           std::cout << pmt << "\t" << _qfit_integ[pmt] << "\t" << _qfit_mpv[pmt] << "\t" << _qfit_sigma[pmt]
0273                     << "\t" << _qfit_integerr[pmt] << "\t" << _qfit_mpverr[pmt] << "\t" << _qfit_sigmaerr[pmt]
0274                     << "\t" << _qfit_chi2ndf[pmt] << std::endl;
0275         }
0276       }
0277     }
0278   }
0279   
0280   if ( std::isnan(_qfit_mpv[0]) )
0281   {
0282     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
0283     _status = -1;
0284     return _status;
0285   }
0286 
0287   return 1;
0288 }
0289 
0290 int MbdCalib::Download_TQT0(const std::string& dbase_location)
0291 {
0292   // Reset All Values
0293   _tqfit_t0mean.fill(std::numeric_limits<float>::quiet_NaN());
0294   _tqfit_t0meanerr.fill(std::numeric_limits<float>::quiet_NaN());
0295   _tqfit_t0sigma.fill(std::numeric_limits<float>::quiet_NaN());
0296   _tqfit_t0sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
0297 
0298   if (Verbosity() > 0)
0299   {
0300     std::cout << "Opening " << dbase_location << std::endl;
0301   }
0302   TString dbase_file = dbase_location;
0303 
0304 #ifndef ONLINE
0305   if (dbase_file.EndsWith(".root"))  // read from database
0306   {
0307     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0308     cdbttree->LoadCalibrations();
0309 
0310     for (int ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
0311     {
0312       _tqfit_t0mean[ipmt] = cdbttree->GetFloatValue(ipmt, "tqfit_t0mean");
0313       _tqfit_t0meanerr[ipmt] = cdbttree->GetFloatValue(ipmt, "tqfit_t0meanerr");
0314       _tqfit_t0sigma[ipmt] = cdbttree->GetFloatValue(ipmt, "tqfit_t0sigma");
0315       _tqfit_t0sigmaerr[ipmt] = cdbttree->GetFloatValue(ipmt, "tqfit_t0sigmaerr");
0316       if (Verbosity() > 0)
0317       {
0318         if (ipmt < 5 || ipmt >= MbdDefs::MBD_N_PMT - 5)
0319         {
0320           std::cout << ipmt << "\t" << _tqfit_t0mean[ipmt] << std::endl;
0321         }
0322       }
0323     }
0324     delete cdbttree;
0325   }
0326 #endif
0327 
0328   if (dbase_file.EndsWith(".calib"))  // read from text file
0329   {
0330     std::ifstream infile(dbase_location);
0331     if (!infile.is_open())
0332     {
0333       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0334       _status = -3;
0335       return _status;
0336     }
0337 
0338     int pmt = -1;
0339     while (infile >> pmt)
0340     {
0341       infile >> _tqfit_t0mean[pmt] >> _tqfit_t0meanerr[pmt] >> _tqfit_t0sigma[pmt] >> _tqfit_t0sigmaerr[pmt];
0342 
0343       if (Verbosity() > 0)
0344       {
0345         if (pmt < 5 || pmt >= MbdDefs::MBD_N_PMT - 5)
0346         {
0347           std::cout << pmt << "\t" << _tqfit_t0mean[pmt] << "\t" << _tqfit_t0meanerr[pmt]
0348                     << "\t" << _tqfit_t0sigma[pmt] << "\t" << _tqfit_t0sigmaerr[pmt] << std::endl;
0349         }
0350       }
0351     }
0352     infile.close();
0353   }
0354 
0355   if ( std::isnan(_tqfit_t0mean[0]) )
0356   {
0357     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
0358     _status = -1;
0359     return _status;
0360   }
0361 
0362   return 1;
0363 }
0364 
0365 int MbdCalib::Download_TTT0(const std::string& dbase_location)
0366 {
0367   // Reset All Values
0368   _ttfit_t0mean.fill(std::numeric_limits<float>::quiet_NaN());
0369   _ttfit_t0meanerr.fill(std::numeric_limits<float>::quiet_NaN());
0370   _ttfit_t0sigma.fill(std::numeric_limits<float>::quiet_NaN());
0371   _ttfit_t0sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
0372 
0373   if (Verbosity() > 0)
0374   {
0375     std::cout << "Opening " << dbase_location << std::endl;
0376   }
0377   TString dbase_file = dbase_location;
0378 
0379 #ifndef ONLINE
0380   if (dbase_file.EndsWith(".root"))  // read from database
0381   {
0382     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0383     cdbttree->LoadCalibrations();
0384 
0385     for (int ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
0386     {
0387       _ttfit_t0mean[ipmt] = cdbttree->GetFloatValue(ipmt, "ttfit_t0mean");
0388       _ttfit_t0meanerr[ipmt] = cdbttree->GetFloatValue(ipmt, "ttfit_t0meanerr");
0389       _ttfit_t0sigma[ipmt] = cdbttree->GetFloatValue(ipmt, "ttfit_t0sigma");
0390       _ttfit_t0sigmaerr[ipmt] = cdbttree->GetFloatValue(ipmt, "ttfit_t0sigmaerr");
0391 
0392       if (Verbosity() > 0)
0393       {
0394         if (ipmt < 5 || ipmt >= MbdDefs::MBD_N_PMT - 5)
0395         {
0396           std::cout << ipmt << "\t" << _ttfit_t0mean[ipmt] << std::endl;
0397         }
0398       }
0399     }
0400     delete cdbttree;
0401   }
0402 #endif
0403 
0404   if (dbase_file.EndsWith(".calib"))  // read from text file
0405   {
0406     std::ifstream infile(dbase_location);
0407     if (!infile.is_open())
0408     {
0409       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0410       _status = -3;
0411       return _status;
0412     }
0413 
0414     int pmt = -1;
0415     while (infile >> pmt)
0416     {
0417       infile >> _ttfit_t0mean[pmt] >> _ttfit_t0meanerr[pmt] >> _ttfit_t0sigma[pmt] >> _ttfit_t0sigmaerr[pmt];
0418 
0419       if (Verbosity() > 0)
0420       {
0421         if (pmt < 5 || pmt >= MbdDefs::MBD_N_PMT - 5)
0422         {
0423           std::cout << pmt << "\t" << _ttfit_t0mean[pmt] << "\t" << _ttfit_t0meanerr[pmt]
0424                     << "\t" << _ttfit_t0sigma[pmt] << "\t" << _ttfit_t0sigmaerr[pmt] << std::endl;
0425         }
0426       }
0427     }
0428     infile.close();
0429   }
0430 
0431   if ( std::isnan(_ttfit_t0mean[0]) )
0432   {
0433     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
0434     _status = -1;
0435     return _status;
0436   }
0437 
0438   return 1;
0439 }
0440 
0441 
0442 int MbdCalib::Download_T0Corr(const std::string& dbase_location)
0443 {
0444   // Reset All Values
0445   _t0corrmean = 0.;
0446   _t0corrmeanerr = 0.;
0447   _t0corr_fitmean.fill(std::numeric_limits<float>::quiet_NaN());
0448   _t0corr_fitmeanerr.fill(std::numeric_limits<float>::quiet_NaN());
0449   _t0corr_fitsigma.fill(std::numeric_limits<float>::quiet_NaN());
0450   _t0corr_fitsigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
0451   _t0corr_hmean.fill(std::numeric_limits<float>::quiet_NaN());
0452   _t0corr_hmeanerr.fill(std::numeric_limits<float>::quiet_NaN());
0453   _t0corr_hstddev.fill(std::numeric_limits<float>::quiet_NaN());
0454   _t0corr_hstddeverr.fill(std::numeric_limits<float>::quiet_NaN());
0455 
0456   if ( dbase_location.empty() )
0457   {
0458     return 0;
0459   }
0460 
0461   if (Verbosity() > 0)
0462   {
0463     std::cout << "Opening " << dbase_location << std::endl;
0464   }
0465   TString dbase_file = dbase_location;
0466 
0467 #ifndef ONLINE
0468   if (dbase_file.EndsWith(".root"))  // read from database
0469   {
0470     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0471     if ( cdbttree == nullptr )
0472     {
0473       std::cerr << "T0Corr not found, skipping" << std::endl;
0474       _status = -1;
0475       return _status;
0476     }
0477     cdbttree->LoadCalibrations();
0478 
0479     _t0corrmean = cdbttree->GetSingleFloatValue("_t0corrmean");
0480     _t0corrmeanerr = cdbttree->GetSingleFloatValue("_t0corrmeanerr");
0481     _t0corr_fitmean[0] = cdbttree->GetSingleFloatValue("_t0corr_fitmean0");
0482     _t0corr_fitmeanerr[0] = cdbttree->GetSingleFloatValue("_t0corr_fitmeanerr0");
0483     _t0corr_fitsigma[0] = cdbttree->GetSingleFloatValue("_t0corr_fitsigma0");
0484     _t0corr_fitsigmaerr[0] = cdbttree->GetSingleFloatValue("_t0corr_fitsigmaerr0");
0485     _t0corr_fitmean[1] = cdbttree->GetSingleFloatValue("_t0corr_fitmean1");
0486     _t0corr_fitmeanerr[1] = cdbttree->GetSingleFloatValue("_t0corr_fitmeanerr1");
0487     _t0corr_fitsigma[1] = cdbttree->GetSingleFloatValue("_t0corr_fitsigma1");
0488     _t0corr_fitsigmaerr[1] = cdbttree->GetSingleFloatValue("_t0corr_fitsigmaerr1");
0489     _t0corr_hmean[0] = cdbttree->GetSingleFloatValue("_t0corr_hmean0");
0490     _t0corr_hmeanerr[0] = cdbttree->GetSingleFloatValue("_t0corr_hmeanerr0");
0491     _t0corr_hstddev[0] = cdbttree->GetSingleFloatValue("_t0corr_hstddev0");
0492     _t0corr_hstddeverr[0] = cdbttree->GetSingleFloatValue("_t0corr_hstddeverr0");
0493     _t0corr_hmean[1] = cdbttree->GetSingleFloatValue("_t0corr_hmean1");
0494     _t0corr_hmeanerr[1] = cdbttree->GetSingleFloatValue("_t0corr_hmeanerr1");
0495     _t0corr_hstddev[1] = cdbttree->GetSingleFloatValue("_t0corr_hstddev1");
0496     _t0corr_hstddeverr[1] = cdbttree->GetSingleFloatValue("_t0corr_hstddeverr1");
0497 
0498     if (Verbosity() > 0)
0499     {
0500       std::cout << "T0Corr\t" << _t0corrmean << std::endl;
0501     }
0502     delete cdbttree;
0503   }
0504 #endif
0505 
0506   if (dbase_file.EndsWith(".calib"))  // read from text file
0507   {
0508     std::ifstream infile(dbase_location);
0509     if (!infile.is_open())
0510     {
0511       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0512       _status = -3;
0513       return _status;
0514     }
0515 
0516     infile >> _t0corrmean >> _t0corrmeanerr
0517       >> _t0corr_fitmean[0] >> _t0corr_fitmeanerr[0] >> _t0corr_fitsigma[0] >> _t0corr_fitsigmaerr[0]
0518       >> _t0corr_fitmean[1] >> _t0corr_fitmeanerr[1] >> _t0corr_fitsigma[1] >> _t0corr_fitsigmaerr[1]
0519       >> _t0corr_hmean[0] >> _t0corr_hmeanerr[0] >> _t0corr_hstddev[0] >> _t0corr_hstddeverr[0]
0520       >> _t0corr_hmean[1] >> _t0corr_hmeanerr[1] >> _t0corr_hstddev[1] >> _t0corr_hstddeverr[1];
0521 
0522     if (Verbosity() > 0)
0523     {
0524       std::cout << "T0Corr\t" << _t0corrmean << "\t" << _t0corrmeanerr << std::endl;
0525     }
0526     infile.close();
0527   }
0528 
0529   if ( std::isnan(_t0corrmean) )
0530   {
0531     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
0532     _status = -1;
0533     return _status;
0534   }
0535 
0536   return 1;
0537 }
0538 
0539 int MbdCalib::Download_Ped(const std::string& dbase_location)
0540 {
0541   // Reset All Values
0542   _pedmean.fill(std::numeric_limits<float>::quiet_NaN());
0543   _pedmeanerr.fill(std::numeric_limits<float>::quiet_NaN());
0544   _pedsigma.fill(std::numeric_limits<float>::quiet_NaN());
0545   _pedsigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
0546 
0547   if (Verbosity() > 0)
0548   {
0549     std::cout << "Opening " << dbase_location << std::endl;
0550   }
0551   TString dbase_file = dbase_location;
0552 
0553 #ifndef ONLINE
0554   if (dbase_file.EndsWith(".root"))  // read from database
0555   {
0556     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0557     cdbttree->LoadCalibrations();
0558 
0559     for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
0560     {
0561       _pedmean[ifeech] = cdbttree->GetFloatValue(ifeech, "pedmean");
0562       _pedmeanerr[ifeech] = cdbttree->GetFloatValue(ifeech, "pedmeanerr");
0563       _pedsigma[ifeech] = cdbttree->GetFloatValue(ifeech, "pedsigma");
0564       _pedsigmaerr[ifeech] = cdbttree->GetFloatValue(ifeech, "pedsigmaerr");
0565 
0566       if (Verbosity() > 0)
0567       {
0568         if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
0569         {
0570           std::cout << ifeech << "\t" << _pedmean[ifeech] << std::endl;
0571         }
0572       }
0573     }
0574     delete cdbttree;
0575   }
0576 #endif
0577 
0578   if (dbase_file.EndsWith(".calib"))  // read from text file
0579   {
0580     std::ifstream infile(dbase_location);
0581     if (!infile.is_open())
0582     {
0583       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0584       _status = -3;
0585       return _status;
0586     }
0587 
0588     int feech = -1;
0589     while (infile >> feech)
0590     {
0591       infile >> _pedmean[feech] >> _pedmeanerr[feech] >> _pedsigma[feech] >> _pedsigmaerr[feech];
0592 
0593       if (Verbosity() > 0)
0594       {
0595         if (feech < 5 || feech >= MbdDefs::MBD_N_FEECH - 5)
0596         {
0597           std::cout << feech << "\t" << _pedmean[feech] << "\t" << _pedmeanerr[feech]
0598                     << "\t" << _pedsigma[feech] << "\t" << _pedsigmaerr[feech] << std::endl;
0599         }
0600       }
0601     }
0602     infile.close();
0603   }
0604 
0605   if ( std::isnan(_pedmean[0]) )
0606   {
0607     std::cout << PHWHERE << ", WARNING, ped calib missing, " << dbase_location << std::endl;
0608     _status = -1;
0609     return _status;
0610   }
0611 
0612   return 1;
0613 }
0614 
0615 int MbdCalib::Download_SampMax(const std::string& dbase_location)
0616 {
0617   // Reset All Values
0618   _sampmax.fill(-1);
0619 
0620   TString dbase_file = dbase_location;
0621 
0622 #ifndef ONLINE
0623   if (dbase_file.EndsWith(".root"))  // read from database
0624   {
0625     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0626     cdbttree->LoadCalibrations();
0627 
0628     for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
0629     {
0630       _sampmax[ifeech] = cdbttree->GetIntValue(ifeech, "sampmax");
0631       if (Verbosity() > 0)
0632       {
0633         if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
0634         {
0635           std::cout << ifeech << "\t" << _sampmax[ifeech] << std::endl;
0636         }
0637       }
0638     }
0639     delete cdbttree;
0640   }
0641 #endif
0642 
0643   if (dbase_file.EndsWith(".calib"))  // read from text file
0644   {
0645     std::ifstream infile(dbase_location);
0646     if (!infile.is_open())
0647     {
0648       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0649       _status = -3;
0650       return _status;
0651     }
0652 
0653     int feech = -1;
0654     while (infile >> feech)
0655     {
0656       infile >> _sampmax[feech];
0657       if (Verbosity() > 0)
0658       {
0659         if (feech < 5 || feech >= MbdDefs::MBD_N_FEECH - 5)
0660         {
0661           std::cout << "sampmax\t" << feech << "\t" << _sampmax[feech] << std::endl;
0662         }
0663       }
0664     }
0665     infile.close();
0666   }
0667   
0668 
0669   if ( _sampmax[0] == -1 )
0670   {
0671     std::cout << PHWHERE << ", WARNING, sampmax calib missing, " << dbase_location << std::endl;
0672     _status = -1;
0673     return _status;  // file not found
0674   }
0675 
0676   return 1;
0677 }
0678 
0679 int MbdCalib::Download_Shapes(const std::string& dbase_location)
0680 {
0681   // Verbosity(100);
0682   if (Verbosity())
0683   {
0684     std::cout << "In MbdCalib::Download_Shapes" << std::endl;
0685   }
0686   // Reset All Values
0687   for (auto& shape : _shape_y)
0688   {
0689     shape.clear();
0690   }
0691   for (auto& sherr : _sherr_yerr)
0692   {
0693     sherr.clear();
0694   }
0695 
0696   TString dbase_file = dbase_location;
0697 
0698 #ifndef ONLINE
0699   if (dbase_file.EndsWith(".root"))  // read from database
0700   {
0701     if (Verbosity())
0702     {
0703       std::cout << "Reading from CDB " << dbase_location << std::endl;
0704     }
0705     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0706     cdbttree->LoadCalibrations();
0707 
0708     for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
0709     {
0710       if (_mbdgeom->get_type(ifeech) == 0)
0711       {
0712         continue;  // skip t-channels
0713       }
0714 
0715       _shape_npts[ifeech] = cdbttree->GetIntValue(ifeech, "shape_npts");
0716       _shape_minrange[ifeech] = cdbttree->GetFloatValue(ifeech, "shape_min");
0717       _shape_maxrange[ifeech] = cdbttree->GetFloatValue(ifeech, "shape_max");
0718 
0719       _sherr_npts[ifeech] = cdbttree->GetIntValue(ifeech, "sherr_npts");
0720       _sherr_minrange[ifeech] = cdbttree->GetFloatValue(ifeech, "sherr_min");
0721       _sherr_maxrange[ifeech] = cdbttree->GetFloatValue(ifeech, "sherr_max");
0722 
0723       for (int ipt = 0; ipt < _shape_npts[ifeech]; ipt++)
0724       {
0725         int chtemp = (1000 * ipt) + ifeech;
0726 
0727         float val = cdbttree->GetFloatValue(chtemp, "shape_val");
0728         _shape_y[ifeech].push_back(val);
0729 
0730         val = cdbttree->GetFloatValue(chtemp, "sherr_val");
0731         _sherr_yerr[ifeech].push_back(val);
0732       }
0733 
0734       if (Verbosity() > 0)
0735       {
0736         if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
0737         {
0738           std::cout << ifeech << "\t" << _shape_y[ifeech][0] << std::endl;
0739         }
0740       }
0741     }
0742     delete cdbttree;
0743   }
0744 #endif
0745 
0746   if (dbase_file.EndsWith(".calib"))  // read from text file
0747   {
0748     if (Verbosity())
0749     {
0750       std::cout << "Reading from " << dbase_location << std::endl;
0751     }
0752     std::ifstream infile(dbase_location);
0753     if (!infile.is_open())
0754     {
0755       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0756       _status = -3;
0757       return _status;
0758     }
0759 
0760     int temp_feech = -1;
0761     int temp_npoints = -1;
0762     float temp_begintime = -1;
0763     float temp_endtime = -1;
0764     while (infile >> temp_feech >> temp_npoints >> temp_begintime >> temp_endtime)
0765     {
0766       if (Verbosity())
0767       {
0768         std::cout << "shape " << temp_feech << "\t" << temp_npoints << "\t" << temp_begintime << "\t" << temp_endtime << std::endl;
0769       }
0770       if (temp_feech < 0 || temp_feech > 255)
0771       {
0772         std::cout << "ERROR, invalid FEECH " << temp_feech << " in MBD waveforms calibration" << std::endl;
0773         _status = -2;
0774         return _status;
0775       }
0776 
0777       _shape_npts[temp_feech] = temp_npoints;
0778       _shape_minrange[temp_feech] = temp_begintime;
0779       _shape_maxrange[temp_feech] = temp_endtime;
0780 
0781       float temp_val{0.};
0782       for (int isamp = 0; isamp < temp_npoints; isamp++)
0783       {
0784         infile >> temp_val;
0785         _shape_y[temp_feech].push_back(temp_val);
0786         if (Verbosity() && (temp_feech == 8 || temp_feech == 255))
0787         {
0788           std::cout << _shape_y[temp_feech][isamp] << " ";
0789           if (isamp % 10 == 9)
0790           {
0791             std::cout << std::endl;
0792           }
0793         }
0794       }
0795       if (Verbosity())
0796       {
0797         std::cout << std::endl;
0798       }
0799     }
0800 
0801     infile.close();
0802 
0803     // Now read in the sherr file
0804     std::string sherr_dbase_location = std::regex_replace(dbase_location, std::regex("bbc_shape.calib"), "bbc_sherr.calib");
0805     if (Verbosity())
0806     {
0807       std::cout << "Reading from " << sherr_dbase_location << std::endl;
0808     }
0809     infile.open(sherr_dbase_location);
0810     if (!infile.is_open())
0811     {
0812       std::cout << PHWHERE << "unable to open " << sherr_dbase_location << std::endl;
0813       _status = -3;
0814       return _status;
0815     }
0816 
0817     temp_feech = -1;
0818     temp_npoints = -1;
0819     temp_begintime = -1;
0820     temp_endtime = -1;
0821     while (infile >> temp_feech >> temp_npoints >> temp_begintime >> temp_endtime)
0822     {
0823       if (Verbosity())
0824       {
0825         std::cout << "sheer " << temp_feech << "\t" << temp_npoints << "\t" << temp_begintime << "\t" << temp_endtime << std::endl;
0826       }
0827       if (temp_feech < 0 || temp_feech > 255)
0828       {
0829         std::cout << "ERROR, invalid FEECH " << temp_feech << " in MBD waveforms calibration" << std::endl;
0830         _status = -2;
0831         return _status;
0832       }
0833 
0834       _sherr_npts[temp_feech] = temp_npoints;
0835       _sherr_minrange[temp_feech] = temp_begintime;
0836       _sherr_maxrange[temp_feech] = temp_endtime;
0837 
0838       float temp_val{0.};
0839       for (int isamp = 0; isamp < temp_npoints; isamp++)
0840       {
0841         infile >> temp_val;
0842         _sherr_yerr[temp_feech].push_back(temp_val);
0843         if (Verbosity() && (temp_feech == 8 || temp_feech == 255))
0844         {
0845           std::cout << _sherr_yerr[temp_feech][isamp] << " ";
0846           if (isamp % 10 == 9)
0847           {
0848             std::cout << std::endl;
0849           }
0850         }
0851       }
0852       if (Verbosity())
0853       {
0854         std::cout << std::endl;
0855       }
0856     }
0857 
0858     infile.close();
0859   }
0860 
0861   if ( _shape_y[8].empty() )
0862   {
0863     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
0864     _status = -1;
0865     return _status;  // file not found
0866   }
0867 
0868   // Verbosity(0);
0869   return 1;
0870 }
0871 
0872 
0873 int MbdCalib::Download_TimeCorr(const std::string& dbase_location)
0874 {
0875   //Verbosity(100);
0876   if ( Verbosity() )
0877   {
0878     std::cout << "In MbdCalib::Download_TimeCorr" << std::endl;
0879   }
0880   // Reset All Values
0881   for(auto& tcorr : _tcorr_y) {
0882     tcorr.clear();
0883   }
0884 
0885   TString dbase_file = dbase_location;
0886 
0887 #ifndef ONLINE
0888   if (dbase_file.EndsWith(".root"))  // read from CDB database file
0889   {
0890     if ( Verbosity() )
0891     {
0892       std::cout << "Reading from CDB " << dbase_location << std::endl;
0893     }
0894     CDBTTree* cdbttree = new CDBTTree(dbase_location);
0895     cdbttree->LoadCalibrations();
0896 
0897     for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
0898     {
0899       if ( _mbdgeom->get_type(ifeech) == 1 )
0900       {
0901         continue;  // skip q-channels
0902       }
0903 
0904       _tcorr_npts[ifeech] = cdbttree->GetIntValue(ifeech, "tcorr_npts");
0905       _tcorr_minrange[ifeech] = cdbttree->GetFloatValue(ifeech, "tcorr_min");
0906       _tcorr_maxrange[ifeech] = cdbttree->GetFloatValue(ifeech, "tcorr_max");
0907 
0908       for (int ipt=0; ipt<_tcorr_npts[ifeech]; ipt++)
0909       {
0910         int chtemp = (1000*ipt) + ifeech; // in cdbtree, entry has id = 1000*datapoint + ifeech
0911 
0912         float val = cdbttree->GetFloatValue(chtemp, "tcorr_val");
0913         _tcorr_y[ifeech].push_back( val );
0914       }
0915 
0916       if (Verbosity() > 0)
0917       {
0918         if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
0919         {
0920           std::cout << ifeech << "\t" << _tcorr_y[ifeech][0] << std::endl;
0921         }
0922       }
0923     }
0924     delete cdbttree;
0925   }
0926 #endif
0927 
0928   if (dbase_file.EndsWith(".calib"))  // read from text file
0929   {
0930     if ( Verbosity() )
0931     {
0932       std::cout << "Reading from " << dbase_location << std::endl;
0933     }
0934     std::ifstream infile(dbase_location);
0935     if (!infile.is_open())
0936     {
0937       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
0938       _status = -3;
0939       return _status;
0940     }
0941 
0942     int temp_feech = -1;
0943     int temp_npoints = -1;
0944     float temp_begintime = -1;
0945     float temp_endtime = -1;
0946     while ( infile >> temp_feech >> temp_npoints >> temp_begintime >> temp_endtime )
0947     {
0948       if ( Verbosity() )
0949       {
0950         std::cout << "tcorr " << temp_feech << "\t" <<  temp_npoints << "\t" <<  temp_begintime << "\t" <<  temp_endtime << std::endl;
0951       }
0952       if ( temp_feech<0 || temp_feech>255 )
0953       {
0954         std::cout << "ERROR, invalid FEECH " << temp_feech << " in MBD timecorr calibration" << std::endl;
0955         _status = -2;
0956         return _status;
0957       }
0958 
0959       _tcorr_npts[temp_feech] = temp_npoints;
0960       _tcorr_minrange[temp_feech] = temp_begintime;
0961       _tcorr_maxrange[temp_feech] = temp_endtime;
0962 
0963       float temp_val{0.};
0964       for (int isamp=0; isamp<temp_npoints; isamp++)
0965       {
0966         infile >> temp_val;
0967         _tcorr_y[temp_feech].push_back( temp_val );
0968         if ( Verbosity() && (temp_feech==0 || temp_feech==64) )
0969         {
0970           std::cout << _tcorr_y[temp_feech][isamp] << " ";
0971           if ( isamp%10==9 )
0972           {
0973             std::cout << std::endl;
0974           }
0975         }
0976       }
0977       if ( Verbosity() )
0978       {
0979         std::cout << std::endl;
0980       }
0981     }
0982 
0983     infile.close();
0984   }
0985 
0986   if ( _tcorr_y[0].empty() )
0987   {
0988     std::cout << PHWHERE << ", ERROR, MBD tcorr not loaded, " << dbase_location << std::endl;
0989     _status = -1;
0990     return _status;  // file not loaded
0991   }
0992 
0993   // Now we interpolate the timecorr
0994   for (size_t ifeech=0; ifeech<MbdDefs::MBD_N_FEECH; ifeech++) 
0995   {
0996     if ( _mbdgeom->get_type(ifeech) == 1 )
0997     {
0998       continue;  // skip q-channels
0999     }
1000 
1001     int step = static_cast<int>( (_tcorr_maxrange[ifeech] - _tcorr_minrange[ifeech]) / (_tcorr_npts[ifeech]-1) );
1002     //std::cout << ifeech << " step = " << step << std::endl;
1003 
1004     for (int itdc=0; itdc<=_tcorr_maxrange[ifeech]; itdc++)
1005     {
1006       int calib_index = itdc/step;
1007       int interp = itdc%step;
1008 
1009       // simple linear interpolation for now
1010       double slope = (_tcorr_y[ifeech][calib_index+1] - _tcorr_y[ifeech][calib_index])/step;
1011       float tcorr_interp = _tcorr_y[ifeech][calib_index] + (interp*slope);
1012  
1013       _tcorr_y_interp[ifeech].push_back( tcorr_interp );
1014 
1015       if ( ifeech==0 && itdc<2*step && Verbosity() )
1016       {
1017         std::cout << _tcorr_y_interp[ifeech][itdc] << " ";
1018         if ( itdc%step==(step-1) )
1019         {
1020           std::cout << std::endl;
1021         }
1022       }
1023     }
1024 
1025   }
1026 
1027   //Verbosity(0);
1028   return 1;
1029 }
1030 
1031 int MbdCalib::Download_SlewCorr(const std::string& dbase_location)
1032 {
1033   //Verbosity(100);
1034   if ( Verbosity() )
1035   {
1036     std::cout << "In MbdCalib::Download_SlewCorr" << std::endl;
1037   }
1038   // Reset All Values
1039   for(auto& scorr : _scorr_y) {
1040     scorr.clear();
1041   }
1042   std::fill(_scorr_npts.begin(), _scorr_npts.end(), 0);
1043   
1044   TString dbase_file = dbase_location;
1045 
1046 #ifndef ONLINE
1047   if (dbase_file.EndsWith(".root"))  // read from CDB database file
1048   {
1049     if ( Verbosity() )
1050     {
1051       std::cout << "Reading from CDB " << dbase_location << std::endl;
1052     }
1053     CDBTTree* cdbttree = new CDBTTree(dbase_location);
1054     cdbttree->LoadCalibrations();
1055 
1056     for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1057     {
1058       if ( _mbdgeom->get_type(ifeech) == 1 )
1059       {
1060         continue;  // skip q-channels
1061       }
1062 
1063       _scorr_npts[ifeech] = cdbttree->GetIntValue(ifeech, "scorr_npts");
1064       _scorr_minrange[ifeech] = cdbttree->GetFloatValue(ifeech, "scorr_min");
1065       _scorr_maxrange[ifeech] = cdbttree->GetFloatValue(ifeech, "scorr_max");
1066 
1067       for (int ipt=0; ipt<_scorr_npts[ifeech]; ipt++)
1068       {
1069         int chtemp = (1000*ipt) + ifeech; // in cdbtree, entry has id = 1000*datapoint + ifeech
1070 
1071         float val = cdbttree->GetFloatValue(chtemp, "scorr_val");
1072         _scorr_y[ifeech].push_back( val );
1073       }
1074 
1075       if (Verbosity() > 0)
1076       {
1077         if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
1078         {
1079           std::cout << ifeech << "\t" << _scorr_y[ifeech][0] << std::endl;
1080         }
1081       }
1082     }
1083     delete cdbttree;
1084   }
1085 #endif
1086 
1087   if (dbase_file.EndsWith(".calib"))  // read from text file
1088   {
1089     if ( Verbosity() )
1090     {
1091       std::cout << "Reading from " << dbase_location << std::endl;
1092     }
1093 
1094     std::ifstream infile(dbase_location);
1095     if (!infile.is_open())
1096     {
1097       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
1098       _status = -3;
1099       return _status;
1100     }
1101 
1102     int temp_feech = -1;
1103     int temp_npoints = 0;
1104     float temp_beginadc = -1;
1105     float temp_endadc = -1;
1106     while ( infile >> temp_feech >> temp_npoints >> temp_beginadc >> temp_endadc )
1107     {
1108       if ( Verbosity() )
1109       {
1110         std::cout << "scorr " << temp_feech << "\t" <<  temp_npoints << "\t" <<  temp_beginadc << "\t" <<  temp_endadc << std::endl;
1111       }
1112 
1113       if ( temp_feech<0 || temp_feech>255 )
1114       {
1115         std::cout << "ERROR, invalid FEECH " << temp_feech << " in MBD slewcorr calibration" << std::endl;
1116         _status = -2;
1117         return _status;
1118       }
1119 
1120       _scorr_npts[temp_feech] = temp_npoints;
1121       _scorr_minrange[temp_feech] = temp_beginadc;
1122       _scorr_maxrange[temp_feech] = temp_endadc;
1123 
1124       float temp_val{0.};
1125       for (int isamp=0; isamp<temp_npoints; isamp++)
1126       {
1127         infile >> temp_val;
1128         _scorr_y[temp_feech].push_back( temp_val );
1129         if ( Verbosity() && (temp_feech==0 || temp_feech==64) )
1130         {
1131           std::cout << _scorr_y[temp_feech][isamp] << " ";
1132           if ( isamp%10==9 )
1133           {
1134             std::cout << std::endl;
1135           }
1136         }
1137       }
1138       if ( Verbosity() )
1139       {
1140         std::cout << std::endl;
1141       }
1142     }
1143 
1144     infile.close();
1145   }
1146 
1147   if ( _scorr_y[0].empty() )
1148   {
1149     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
1150     _status = -1;
1151     return _status;  // file not found
1152   }
1153 
1154   // Now we interpolate the slewcorr
1155   for (size_t ifeech=0; ifeech<MbdDefs::MBD_N_FEECH; ifeech++) 
1156   {
1157     if ( _mbdgeom->get_type(ifeech) == 1 )
1158     {
1159       continue;  // skip q-channels
1160     }
1161     // skip bad t-channels
1162     if ( _scorr_npts[ifeech] == 0 )
1163     {
1164       //std::cout << "skipping " << ifeech << std::endl;
1165       continue;
1166     }
1167 
1168     int step = static_cast<int>( (_scorr_maxrange[ifeech] - _scorr_minrange[ifeech]) / (_scorr_npts[ifeech]-1) );
1169     //std::cout << ifeech << " step = " << step << std::endl;
1170 
1171     for (int iadc=0; iadc<=_scorr_maxrange[ifeech]; iadc++)
1172     {
1173       int calib_index = iadc/step;
1174       int interp = iadc%step;
1175 
1176       // simple linear interpolation for now
1177       double slope = (_scorr_y[ifeech][calib_index+1] - _scorr_y[ifeech][calib_index])/step;
1178       float scorr_interp = _scorr_y[ifeech][calib_index] + (interp*slope);
1179  
1180       _scorr_y_interp[ifeech].push_back( scorr_interp );
1181 
1182 
1183       if ( ifeech==4 && iadc<12 && Verbosity() )
1184       {
1185         if ( iadc==0 )
1186         {
1187           std::cout << "slewcorr " << ifeech << "\t" << _scorr_npts[ifeech] << "\t"
1188             << _scorr_minrange[ifeech] << "\t" << _scorr_maxrange[ifeech] << std::endl;
1189         }
1190         std::cout << _scorr_y_interp[ifeech][iadc] << " ";
1191         if ( iadc%step==(step-1) )
1192         {
1193           std::cout << std::endl;
1194         }
1195       }
1196     }
1197 
1198   }
1199 
1200   //Verbosity(0);
1201   return 1;
1202 }
1203 
1204 int MbdCalib::Download_Pileup(const std::string& dbase_location)
1205 {
1206   // Reset All Values
1207   Reset_Pileup();
1208 
1209   if (Verbosity() > 0)
1210   {
1211     std::cout << "Opening " << dbase_location << std::endl;
1212   }
1213   TString dbase_file = dbase_location;
1214 
1215 #ifndef ONLINE
1216   if (dbase_file.EndsWith(".root"))  // read from database
1217   {
1218     CDBTTree* cdbttree = new CDBTTree(dbase_location);
1219     if ( cdbttree == nullptr )
1220     {
1221       std::cerr << "MBD pileup calib not found, skipping" << std::endl;
1222       _status = -1;
1223       return _status;
1224     }
1225     cdbttree->LoadCalibrations();
1226 
1227     for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1228     {
1229       _pileup_p0[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_p0");
1230       _pileup_p0err[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_p0err");
1231       _pileup_p1[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_p1");
1232       _pileup_p1err[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_p1err");
1233       _pileup_p2[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_p2");
1234       _pileup_p2err[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_p2err");
1235       _pileup_chi2ndf[ifeech] = cdbttree->GetFloatValue(ifeech, "pileup_chi2ndf");
1236       if (Verbosity() > 0)
1237       {
1238         if (ifeech < 2 || ifeech >= (MbdDefs::MBD_N_FEECH-2) )
1239         {
1240           std::cout << ifeech << "\t" << _pileup_p0[ifeech] << std::endl;
1241         }
1242       }
1243     }
1244     delete cdbttree;
1245   }
1246 #endif
1247 
1248   if (dbase_file.EndsWith(".calib"))  // read from text file
1249   {
1250     std::ifstream infile(dbase_location);
1251     if (!infile.is_open())
1252     {
1253       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
1254       _status = -3;
1255       return _status;
1256     }
1257 
1258     int feech = -1;
1259     while (infile >> feech)
1260     {
1261       infile >> _pileup_p0[feech] >> _pileup_p0err[feech]
1262              >> _pileup_p1[feech] >> _pileup_p1err[feech]
1263              >> _pileup_p2[feech] >> _pileup_p2err[feech]
1264              >> _pileup_chi2ndf[feech];
1265 
1266       if (Verbosity() > 0)
1267       {
1268         if (feech < 2 || feech >= MbdDefs::MBD_N_PMT - 2)
1269         {
1270           std::cout << feech << "\t" << _pileup_p0[feech] << "\t" << _pileup_p0err[feech]
1271                              << "\t" << _pileup_p1[feech] << "\t" << _pileup_p1err[feech]
1272                              << "\t" << _pileup_p2[feech] << "\t" << _pileup_p2err[feech]
1273                              << "\t" << _pileup_chi2ndf[feech] << std::endl;
1274         }
1275       }
1276     }
1277   }
1278   
1279   if ( std::isnan(_pileup_p0[0]) )
1280   {
1281     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
1282     _status = -1;
1283     return _status;
1284   }
1285 
1286   return 1;
1287 }
1288 
1289 int MbdCalib::Download_Thresholds(const std::string& dbase_location)
1290 {
1291   // Reset All Values
1292   Reset_Thresholds();
1293 
1294   if (Verbosity() > 0)
1295   {
1296     std::cout << "Opening " << dbase_location << std::endl;
1297   }
1298   TString dbase_file = dbase_location;
1299 
1300 #ifndef ONLINE
1301   if (dbase_file.EndsWith(".root"))  // read from database
1302   {
1303     CDBTTree* cdbttree = new CDBTTree(dbase_location);
1304     cdbttree->LoadCalibrations();
1305 
1306     for (int ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1307     {
1308       _thresh_mean[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_mean");
1309       _thresh_meanerr[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_meanerr");
1310       _thresh_width[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_width");
1311       _thresh_widtherr[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_widtherr");
1312       _thresh_eff[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_eff");
1313       _thresh_efferr[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_efferr");
1314       _thresh_chi2ndf[ipmt] = cdbttree->GetFloatValue(ipmt, "thresh_chi2ndf");
1315       if (Verbosity() > 0)
1316       {
1317         if (ipmt < 5)
1318         {
1319           std::cout << ipmt << "\t" << _thresh_mean[ipmt] << std::endl;
1320         }
1321       }
1322     }
1323     delete cdbttree;
1324   }
1325 #endif
1326 
1327   if (dbase_file.EndsWith(".calib"))  // read from text file
1328   {
1329     std::ifstream infile(dbase_location);
1330     if (!infile.is_open())
1331     {
1332       std::cout << PHWHERE << "unable to open " << dbase_location << std::endl;
1333       _status = -3;
1334       return _status;
1335     }
1336 
1337     int pmt = -1;
1338     while (infile >> pmt)
1339     {
1340       infile >> _thresh_mean[pmt] >> _thresh_meanerr[pmt] >> _thresh_width[pmt] >> _thresh_widtherr[pmt] >> _thresh_eff[pmt] >> _thresh_efferr[pmt] >> _thresh_chi2ndf[pmt];
1341       if (Verbosity() > 0)
1342       {
1343         if (pmt < 5 || pmt >= MbdDefs::MBD_N_PMT - 5)
1344         {
1345           std::cout << pmt << "\t" << _thresh_mean[pmt] << "\t" << _thresh_meanerr[pmt] << "\t" << _thresh_width[pmt]
1346                     << "\t" << _thresh_widtherr[pmt] << "\t" << _thresh_eff[pmt] << "\t" << _thresh_efferr[pmt]
1347                     << "\t" << _thresh_chi2ndf[pmt] << std::endl;
1348         }
1349       }
1350     }
1351   }
1352   
1353   if ( std::isnan(_thresh_mean[0]) )
1354   {
1355     std::cout << PHWHERE << ", ERROR, unknown file type, " << dbase_location << std::endl;
1356     _status = -1;
1357     return _status;
1358   }
1359 
1360   return 1;
1361 }
1362 
1363 #ifndef ONLINE
1364 int MbdCalib::Write_CDB_SampMax(const std::string& dbfile)
1365 {
1366   CDBTTree* cdbttree{ nullptr };
1367 
1368   std::cout << "Creating " << dbfile << std::endl;
1369   cdbttree = new CDBTTree( dbfile );
1370   cdbttree->SetSingleIntValue("version", 1);
1371   cdbttree->CommitSingle();
1372 
1373   std::cout << "SAMPMAX" << std::endl;
1374   for (size_t ifeech = 0; ifeech < _sampmax.size(); ifeech++)
1375   {
1376     // store in a CDBTree
1377     cdbttree->SetIntValue(ifeech, "sampmax", _sampmax[ifeech]);
1378 
1379     if (ifeech < 12 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
1380     {
1381       std::cout << ifeech << "\t" << cdbttree->GetIntValue(ifeech, "sampmax") << std::endl;
1382     }
1383   }
1384 
1385   cdbttree->Commit();
1386   // cdbttree->Print();
1387 
1388   // for now we create the tree after reading it
1389   cdbttree->WriteCDBTTree();
1390   delete cdbttree;
1391 
1392   return 1;
1393 }
1394 #endif
1395 
1396 int MbdCalib::Write_SampMax(const std::string& dbfile)
1397 {
1398   std::ofstream cal_file;
1399   cal_file.open(dbfile);
1400   for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1401   {
1402     cal_file << ifeech << "\t" << _sampmax[ifeech] << std::endl;
1403   }
1404   cal_file.close();
1405 
1406   return 1;
1407 }
1408 
1409 #ifndef ONLINE
1410 int MbdCalib::Write_CDB_TTT0(const std::string& dbfile)
1411 {
1412   CDBTTree* cdbttree{ nullptr };
1413 
1414   std::cout << "Creating " << dbfile << std::endl;
1415   cdbttree = new CDBTTree( dbfile );
1416   cdbttree->SetSingleIntValue("version", 1);
1417   cdbttree->CommitSingle();
1418 
1419   std::cout << "TTT0" << std::endl;
1420   for (size_t ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1421   {
1422     // store in a CDBTree
1423     cdbttree->SetFloatValue(ipmt, "ttfit_t0mean", _ttfit_t0mean[ipmt]);
1424     cdbttree->SetFloatValue(ipmt, "ttfit_t0meanerr", _ttfit_t0meanerr[ipmt]);
1425     cdbttree->SetFloatValue(ipmt, "ttfit_t0sigma", _ttfit_t0sigma[ipmt]);
1426     cdbttree->SetFloatValue(ipmt, "ttfit_t0sigmaerr", _ttfit_t0sigmaerr[ipmt]);
1427 
1428     if (ipmt < 5 || ipmt >= MbdDefs::MBD_N_PMT - 5)
1429     {
1430       std::cout << ipmt << "\t" << cdbttree->GetFloatValue(ipmt, "ttfit_t0mean") << std::endl;
1431     }
1432   }
1433 
1434   cdbttree->Commit();
1435   // cdbttree->Print();
1436 
1437   // for now we create the tree after reading it
1438   cdbttree->WriteCDBTTree();
1439   delete cdbttree;
1440 
1441   return 1;
1442 }
1443 #endif
1444 
1445 int MbdCalib::Write_TTT0(const std::string& dbfile)
1446 {
1447   std::ofstream cal_t0_file;
1448   cal_t0_file.open(dbfile);
1449   for (int ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1450   {
1451     cal_t0_file << ipmt << "\t" << _ttfit_t0mean[ipmt] << "\t" << _ttfit_t0meanerr[ipmt]
1452       << "\t" << _ttfit_t0sigma[ipmt] << "\t" << _ttfit_t0sigmaerr[ipmt] << std::endl;
1453   }
1454   cal_t0_file.close();
1455 
1456   return 1;
1457 }
1458 
1459 #ifndef ONLINE
1460 int MbdCalib::Write_CDB_TQT0(const std::string& dbfile)
1461 {
1462   CDBTTree* cdbttree{ nullptr };
1463 
1464   std::cout << "Creating " << dbfile << std::endl;
1465   cdbttree = new CDBTTree( dbfile );
1466   cdbttree->SetSingleIntValue("version", 1);
1467   cdbttree->CommitSingle();
1468 
1469   std::cout << "TQT0" << std::endl;
1470   for (size_t ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1471   {
1472     // store in a CDBTree
1473     cdbttree->SetFloatValue(ipmt, "tqfit_t0mean", _tqfit_t0mean[ipmt]);
1474     cdbttree->SetFloatValue(ipmt, "tqfit_t0meanerr", _tqfit_t0meanerr[ipmt]);
1475     cdbttree->SetFloatValue(ipmt, "tqfit_t0sigma", _tqfit_t0sigma[ipmt]);
1476     cdbttree->SetFloatValue(ipmt, "tqfit_t0sigmaerr", _tqfit_t0sigmaerr[ipmt]);
1477 
1478     if (ipmt < 5 || ipmt >= MbdDefs::MBD_N_PMT - 5)
1479     {
1480       std::cout << ipmt << "\t" << cdbttree->GetFloatValue(ipmt, "tqfit_t0mean") << std::endl;
1481     }
1482   }
1483 
1484   cdbttree->Commit();
1485   // cdbttree->Print();
1486 
1487   // for now we create the tree after reading it
1488   cdbttree->WriteCDBTTree();
1489   delete cdbttree;
1490 
1491   return 1;
1492 }
1493 #endif
1494 
1495 int MbdCalib::Write_TQT0(const std::string& dbfile)
1496 {
1497   std::ofstream cal_t0_file;
1498   cal_t0_file.open(dbfile);
1499   for (int ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1500   {
1501     cal_t0_file << ipmt << "\t" << _tqfit_t0mean[ipmt] << "\t" << _tqfit_t0meanerr[ipmt]
1502       << "\t" << _tqfit_t0sigma[ipmt] << "\t" << _tqfit_t0sigmaerr[ipmt] << std::endl;
1503   }
1504   cal_t0_file.close();
1505 
1506   return 1;
1507 }
1508 
1509 #ifndef ONLINE
1510 int MbdCalib::Write_CDB_T0Corr(const std::string& dbfile)
1511 {
1512   CDBTTree* cdbttree{ nullptr };
1513 
1514   std::cout << "Creating " << dbfile << std::endl;
1515   cdbttree = new CDBTTree( dbfile );
1516   cdbttree->SetSingleIntValue("version", 1);
1517 
1518   std::cout << "T0Corr" << std::endl;
1519 
1520   // store in a CDBTree
1521   cdbttree->SetSingleFloatValue("_t0corrmean", _t0corrmean);
1522   cdbttree->SetSingleFloatValue("_t0corrmeanerr", _t0corrmeanerr);
1523   cdbttree->SetSingleFloatValue("_t0corr_fitmean0", _t0corr_fitmean[0]);
1524   cdbttree->SetSingleFloatValue("_t0corr_fitmeanerr0", _t0corr_fitmeanerr[0]);
1525   cdbttree->SetSingleFloatValue("_t0corr_fitsigma0", _t0corr_fitsigma[0]);
1526   cdbttree->SetSingleFloatValue("_t0corr_fitsigmaerr0", _t0corr_fitsigmaerr[0]);
1527   cdbttree->SetSingleFloatValue("_t0corr_fitmean1", _t0corr_fitmean[1]);
1528   cdbttree->SetSingleFloatValue("_t0corr_fitmeanerr1", _t0corr_fitmeanerr[1]);
1529   cdbttree->SetSingleFloatValue("_t0corr_fitsigma1", _t0corr_fitsigma[1]);
1530   cdbttree->SetSingleFloatValue("_t0corr_fitsigmaerr1", _t0corr_fitsigmaerr[1]);
1531   cdbttree->SetSingleFloatValue("_t0corr_hmean0", _t0corr_hmean[0]);
1532   cdbttree->SetSingleFloatValue("_t0corr_hmeanerr0", _t0corr_hmeanerr[0]);
1533   cdbttree->SetSingleFloatValue("_t0corr_hstddev0", _t0corr_hstddev[0]);
1534   cdbttree->SetSingleFloatValue("_t0corr_hstddeverr0", _t0corr_hstddeverr[0]);
1535   cdbttree->SetSingleFloatValue("_t0corr_hmean1", _t0corr_hmean[1]);
1536   cdbttree->SetSingleFloatValue("_t0corr_hmeanerr1", _t0corr_hmeanerr[1]);
1537   cdbttree->SetSingleFloatValue("_t0corr_hstddev1", _t0corr_hstddev[1]);
1538   cdbttree->SetSingleFloatValue("_t0corr_hstddeverr1", _t0corr_hstddeverr[1]);
1539 
1540   std::cout << "T0Corr\t" << cdbttree->GetSingleFloatValue("_t0corrmean") << std::endl;
1541 
1542   cdbttree->CommitSingle();
1543   //cdbttree->Commit();
1544   //cdbttree->Print();
1545 
1546   // for now we create the tree after reading it
1547   cdbttree->WriteCDBTTree();
1548   delete cdbttree;
1549 
1550   return 1;
1551 }
1552 #endif
1553 
1554 int MbdCalib::Write_T0Corr(const std::string& dbfile)
1555 {
1556   std::ofstream cal_t0corr_file;
1557   cal_t0corr_file.open(dbfile);
1558   cal_t0corr_file << _t0corrmean << "\t" << _t0corrmeanerr << std::endl;
1559   cal_t0corr_file << _t0corr_fitmean[0] << "\t" << _t0corr_fitmeanerr[0] << "\t" << _t0corr_fitsigma[0] << "\t" << _t0corr_fitsigmaerr[0] << std::endl;
1560   cal_t0corr_file << _t0corr_fitmean[1] << "\t" << _t0corr_fitmeanerr[1] << "\t" << _t0corr_fitsigma[1] << "\t" << _t0corr_fitsigmaerr[1] << std::endl;
1561   cal_t0corr_file << _t0corr_hmean[0] << "\t" << _t0corr_hmeanerr[0] << "\t" << _t0corr_hstddev[0] << "\t" << _t0corr_hstddeverr[0] << std::endl;
1562   cal_t0corr_file << _t0corr_hmean[1] << "\t" << _t0corr_hmeanerr[1] << "\t" << _t0corr_hstddev[1] << "\t" << _t0corr_hstddeverr[1] << std::endl;
1563   cal_t0corr_file.close();
1564 
1565   return 1;
1566 }
1567 
1568 #ifndef ONLINE
1569 int MbdCalib::Write_CDB_Ped(const std::string& dbfile)
1570 {
1571   CDBTTree* cdbttree{ nullptr };
1572 
1573   std::cout << "Creating " << dbfile << std::endl;
1574   cdbttree = new CDBTTree( dbfile );
1575   cdbttree->SetSingleIntValue("version", 1);
1576 
1577   std::cout << "Ped" << std::endl;
1578   for (size_t ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1579   {
1580     // store in a CDBTree
1581     cdbttree->SetFloatValue(ifeech, "pedmean", _pedmean[ifeech]);
1582     cdbttree->SetFloatValue(ifeech, "pedmeanerr", _pedmeanerr[ifeech]);
1583     cdbttree->SetFloatValue(ifeech, "pedsigma", _pedsigma[ifeech]);
1584     cdbttree->SetFloatValue(ifeech, "pedsigmaerr", _pedsigmaerr[ifeech]);
1585 
1586     if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
1587     {
1588       std::cout << ifeech << "\t" << cdbttree->GetFloatValue(ifeech, "pedmean") << std::endl;
1589     }
1590   }
1591   cdbttree->CommitSingle();
1592 
1593   cdbttree->Commit();
1594   // cdbttree->Print();
1595 
1596   // for now we create the tree after reading it
1597   cdbttree->WriteCDBTTree();
1598   delete cdbttree;
1599 
1600   return 1;
1601 }
1602 #endif
1603 
1604 int MbdCalib::Write_Ped(const std::string& dbfile)
1605 {
1606   std::ofstream cal_ped_file;
1607   cal_ped_file.open(dbfile);
1608   for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1609   {
1610     cal_ped_file << ifeech << "\t" << _pedmean[ifeech] << "\t" << _pedmeanerr[ifeech]
1611       << "\t" << _pedsigma[ifeech] << "\t" << _pedsigmaerr[ifeech] << std::endl;
1612   }
1613   cal_ped_file.close();
1614 
1615   return 1;
1616 }
1617 
1618 #ifndef ONLINE
1619 int MbdCalib::Write_CDB_Shapes(const std::string& dbfile)
1620 {
1621   // store in a CDBTree
1622   CDBTTree* cdbttree{nullptr};
1623 
1624   std::cout << "Creating " << dbfile << std::endl;
1625   cdbttree = new CDBTTree(dbfile);
1626   cdbttree->SetSingleIntValue("version", 1);
1627   cdbttree->CommitSingle();
1628 
1629   std::cout << "SHAPES" << std::endl;
1630   for (unsigned int ifeech = 0; ifeech < _sampmax.size(); ifeech++)
1631   {
1632     if (_mbdgeom->get_type(ifeech) == 0)
1633     {
1634       continue;  // skip t-channels
1635     }
1636 
1637     cdbttree->SetIntValue(ifeech, "shape_npts", _shape_npts[ifeech]);
1638     cdbttree->SetFloatValue(ifeech, "shape_min", _shape_minrange[ifeech]);
1639     cdbttree->SetFloatValue(ifeech, "shape_max", _shape_maxrange[ifeech]);
1640 
1641     for (int ipt = 0; ipt < _shape_npts[ifeech]; ipt++)
1642     {
1643       int temp_ch = (ipt * 1000) + ifeech;
1644       cdbttree->SetFloatValue(temp_ch, "shape_val", _shape_y[ifeech][ipt]);
1645     }
1646 
1647     cdbttree->SetIntValue(ifeech, "sherr_npts", _sherr_npts[ifeech]);
1648     cdbttree->SetFloatValue(ifeech, "sherr_min", _sherr_minrange[ifeech]);
1649     cdbttree->SetFloatValue(ifeech, "sherr_max", _sherr_maxrange[ifeech]);
1650 
1651     for (int ipt = 0; ipt < _shape_npts[ifeech]; ipt++)
1652     {
1653       int temp_ch = (ipt * 1000) + ifeech;
1654       cdbttree->SetFloatValue(temp_ch, "sherr_val", _sherr_yerr[ifeech][ipt]);
1655     }
1656   }
1657 
1658   cdbttree->Commit();
1659   // cdbttree->Print();
1660 
1661   for (unsigned int ifeech = 0; ifeech < _sampmax.size(); ifeech++)
1662   {
1663     if (_mbdgeom->get_type(ifeech) == 0)
1664     {
1665       continue;  // skip t-channels
1666     }
1667 
1668     if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_FEECH - 5)
1669     {
1670       std::cout << ifeech << "\t" << cdbttree->GetIntValue(ifeech, "shape_npts") << std::endl;
1671       for (int ipt = 0; ipt < 10; ipt++)
1672       {
1673         int temp_ch = (ipt * 1000) + (int)ifeech;
1674         std::cout << cdbttree->GetFloatValue(temp_ch, "shape_val") << "  ";
1675       }
1676       std::cout << std::endl;
1677     }
1678   }
1679 
1680   // for now we create the tree after reading it
1681   cdbttree->WriteCDBTTree();
1682   delete cdbttree;
1683 
1684   return 1;
1685 }
1686 #endif
1687 
1688 #ifndef ONLINE
1689 int MbdCalib::Write_CDB_TimeCorr(const std::string& dbfile)
1690 {
1691   // store in a CDBTree
1692   CDBTTree *cdbttree {nullptr};
1693 
1694   std::cout << "Creating " << dbfile << std::endl;
1695   cdbttree = new CDBTTree( dbfile );
1696   cdbttree->SetSingleIntValue("version",1);
1697   cdbttree->CommitSingle();
1698 
1699   std::cout << "TIMECORR" << std::endl;
1700   for (int ifeech=0; ifeech<MbdDefs::MBD_N_FEECH; ifeech++) 
1701   {
1702     if ( _mbdgeom->get_type(ifeech) == 1 )
1703     {
1704       continue;  // skip q-channels
1705     }
1706 
1707     cdbttree->SetIntValue(ifeech,"tcorr_npts",_tcorr_npts[ifeech]);
1708     cdbttree->SetFloatValue(ifeech,"tcorr_min",_tcorr_minrange[ifeech]);
1709     cdbttree->SetFloatValue(ifeech,"tcorr_max",_tcorr_maxrange[ifeech]);
1710 
1711     for (int ipt=0; ipt<_tcorr_npts[ifeech]; ipt++)
1712     {
1713       int temp_ch = (ipt*1000) + ifeech;
1714       cdbttree->SetFloatValue(temp_ch,"tcorr_val",_tcorr_y[ifeech][ipt]);
1715     }
1716   }
1717 
1718   cdbttree->Commit();
1719   //cdbttree->Print();
1720 
1721   for (size_t ifeech=0; ifeech<MbdDefs::MBD_N_FEECH; ifeech++) 
1722   {
1723     if ( _mbdgeom->get_type(ifeech) == 1 )
1724     {
1725       continue;  // skip q-channels
1726     }
1727 
1728     if ( ifeech<5 || ifeech>=MbdDefs::MBD_N_FEECH-5-8 )
1729     {
1730       std::cout << ifeech << "\t" <<  cdbttree->GetIntValue(ifeech,"tcorr_npts") << std::endl;
1731       for (int ipt=0; ipt<10; ipt++)
1732       {
1733         int temp_ch = (ipt*1000) + (int)ifeech;
1734         std::cout << cdbttree->GetFloatValue(temp_ch,"tcorr_val") << "  ";
1735       }
1736       std::cout << std::endl;
1737     }
1738   }
1739 
1740   // for now we create the tree after reading it
1741   cdbttree->WriteCDBTTree();
1742   delete cdbttree;
1743 
1744   return 1;
1745 }
1746 #endif
1747 
1748 #ifndef ONLINE
1749 int MbdCalib::Write_CDB_SlewCorr(const std::string& dbfile)
1750 {
1751   // store in a CDBTree
1752   CDBTTree *cdbttree {nullptr};
1753 
1754   std::cout << "Creating " << dbfile << std::endl;
1755   cdbttree = new CDBTTree( dbfile );
1756   cdbttree->SetSingleIntValue("version",1);
1757   cdbttree->CommitSingle();
1758 
1759   std::cout << "SLEWCORR" << std::endl;
1760   //for (size_t ifeech=0; ifeech<_sampmax.size(); ifeech++) 
1761   for (size_t ifeech=0; ifeech<MbdDefs::MBD_N_FEECH; ifeech++) 
1762   {
1763     if ( _mbdgeom->get_type(ifeech) == 1 )
1764     {
1765       continue;  // skip q-channels
1766     }
1767 
1768     cdbttree->SetIntValue(ifeech,"scorr_npts",_scorr_npts[ifeech]);
1769     cdbttree->SetFloatValue(ifeech,"scorr_min",_scorr_minrange[ifeech]);
1770     cdbttree->SetFloatValue(ifeech,"scorr_max",_scorr_maxrange[ifeech]);
1771 
1772     for (int ipt=0; ipt<_scorr_npts[ifeech]; ipt++)
1773     {
1774       int temp_ch = (ipt*1000) + (int)ifeech;
1775       cdbttree->SetFloatValue(temp_ch,"scorr_val",_scorr_y[ifeech][ipt]);
1776     }
1777   }
1778 
1779   cdbttree->Commit();
1780   //cdbttree->Print();
1781 
1782   for (size_t ifeech=0; ifeech<MbdDefs::MBD_N_FEECH; ifeech++) 
1783   {
1784     if ( _mbdgeom->get_type(ifeech) == 1 )
1785     {
1786       continue;  // skip q-channels
1787     }
1788 
1789     if ( ifeech<5 || ifeech>=MbdDefs::MBD_N_FEECH-8-5 )
1790     {
1791       std::cout << ifeech << "\t" <<  cdbttree->GetIntValue(ifeech,"scorr_npts") << std::endl;
1792       for (int ipt=0; ipt<10; ipt++)
1793       {
1794         int temp_ch = (ipt*1000) + (int)ifeech;
1795         std::cout << cdbttree->GetFloatValue(temp_ch,"scorr_val") << "  ";
1796       }
1797       std::cout << std::endl;
1798     }
1799   }
1800 
1801   // for now we create the tree after reading it
1802   cdbttree->WriteCDBTTree();
1803   delete cdbttree;
1804 
1805   return 1;
1806 }
1807 #endif
1808 
1809 #ifndef ONLINE
1810 int MbdCalib::Write_CDB_Gains(const std::string& dbfile)
1811 {
1812   CDBTTree* cdbttree{ nullptr };
1813 
1814   std::cout << "Creating " << dbfile << std::endl;
1815   cdbttree = new CDBTTree( dbfile );
1816   cdbttree->SetSingleIntValue("version", 1);
1817   cdbttree->CommitSingle();
1818 
1819   std::cout << "MBD_QFIT" << std::endl;
1820   for (size_t ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1821   {
1822     // store in a CDBTree
1823     cdbttree->SetFloatValue(ipmt, "qfit_integ", _qfit_integ[ipmt]);
1824     cdbttree->SetFloatValue(ipmt, "qfit_mpv", _qfit_mpv[ipmt]);
1825     cdbttree->SetFloatValue(ipmt, "qfit_sigma", _qfit_sigma[ipmt]);
1826     cdbttree->SetFloatValue(ipmt, "qfit_integerr", _qfit_integerr[ipmt]);
1827     cdbttree->SetFloatValue(ipmt, "qfit_mpverr", _qfit_mpverr[ipmt]);
1828     cdbttree->SetFloatValue(ipmt, "qfit_sigmaerr", _qfit_sigmaerr[ipmt]);
1829     cdbttree->SetFloatValue(ipmt, "qfit_chi2ndf", _qfit_chi2ndf[ipmt]);
1830 
1831     if (ipmt < 5 || ipmt >= MbdDefs::MBD_N_PMT - 5)
1832     {
1833       std::cout << ipmt << "\t" << cdbttree->GetFloatValue(ipmt, "qfit_mpv") << std::endl;
1834     }
1835   }
1836 
1837   cdbttree->Commit();
1838   // cdbttree->Print();
1839 
1840   // for now we create the tree after reading it
1841   cdbttree->WriteCDBTTree();
1842   delete cdbttree;
1843 
1844   return 1;
1845 }
1846 #endif
1847 
1848 int MbdCalib::Write_Gains(const std::string& dbfile)
1849 {
1850   std::ofstream cal_gains_file;
1851   cal_gains_file.open(dbfile);
1852   for (int ipmtch = 0; ipmtch < MbdDefs::MBD_N_PMT; ipmtch++)
1853   {
1854     cal_gains_file << ipmtch << "\t" << _qfit_integ[ipmtch] << "\t" << _qfit_mpv[ipmtch]
1855       << "\t" << _qfit_sigma[ipmtch] << "\t" << _qfit_integerr[ipmtch]
1856       << "\t" << _qfit_mpverr[ipmtch] << "\t" << _qfit_sigmaerr[ipmtch]
1857       << "\t" << _qfit_chi2ndf[ipmtch] << std::endl;
1858   }
1859   cal_gains_file.close();
1860 
1861   return 1;
1862 }
1863 
1864 #ifndef ONLINE
1865 int MbdCalib::Write_CDB_Pileup(const std::string& dbfile)
1866 {
1867   CDBTTree* cdbttree{ nullptr };
1868 
1869   std::cout << "Creating " << dbfile << std::endl;
1870   cdbttree = new CDBTTree( dbfile );
1871   cdbttree->SetSingleIntValue("version", 1);
1872   cdbttree->CommitSingle();
1873 
1874   std::cout << "MBD_PILEUP" << std::endl;
1875   for (size_t ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1876   {
1877     // store in a CDBTree
1878     cdbttree->SetFloatValue(ifeech, "pileup_p0", _pileup_p0[ifeech]);
1879     cdbttree->SetFloatValue(ifeech, "pileup_p0err", _pileup_p0err[ifeech]);
1880     cdbttree->SetFloatValue(ifeech, "pileup_p1", _pileup_p1[ifeech]);
1881     cdbttree->SetFloatValue(ifeech, "pileup_p1err", _pileup_p1err[ifeech]);
1882     cdbttree->SetFloatValue(ifeech, "pileup_p2", _pileup_p2[ifeech]);
1883     cdbttree->SetFloatValue(ifeech, "pileup_p2err", _pileup_p2err[ifeech]);
1884     cdbttree->SetFloatValue(ifeech, "pileup_chi2ndf", _pileup_chi2ndf[ifeech]);
1885 
1886     if (ifeech < 5 || ifeech >= MbdDefs::MBD_N_PMT - 5)
1887     {
1888       std::cout << ifeech << "\t" << cdbttree->GetFloatValue(ifeech, "pileup_p0") << std::endl;
1889     }
1890   }
1891 
1892   cdbttree->Commit();
1893   // cdbttree->Print();
1894 
1895   // for now we create the tree after reading it
1896   cdbttree->WriteCDBTTree();
1897   delete cdbttree;
1898 
1899   return 1;
1900 }
1901 #endif
1902 
1903 int MbdCalib::Write_Pileup(const std::string& dbfile)
1904 {
1905   std::ofstream cal_pileup_file;
1906   cal_pileup_file.open(dbfile);
1907   for (int ifeech = 0; ifeech < MbdDefs::MBD_N_FEECH; ifeech++)
1908   {
1909     cal_pileup_file << ifeech << "\t" << _pileup_p0[ifeech] << "\t" << _pileup_p0err[ifeech]
1910       << "\t" << _pileup_p1[ifeech] << "\t" << _pileup_p1err[ifeech]
1911       << "\t" << _pileup_p2[ifeech] << "\t" << _pileup_p2err[ifeech]
1912       << "\t" << _pileup_chi2ndf[ifeech] << std::endl;
1913   }
1914   cal_pileup_file.close();
1915 
1916   return 1;
1917 }
1918 
1919 #ifndef ONLINE
1920 int MbdCalib::Write_CDB_Thresholds(const std::string& dbfile)
1921 {
1922   CDBTTree* cdbttree{ nullptr };
1923 
1924   std::cout << "Creating " << dbfile << std::endl;
1925   cdbttree = new CDBTTree( dbfile );
1926   cdbttree->SetSingleIntValue("version", 1);
1927   cdbttree->CommitSingle();
1928 
1929   std::cout << "MBD_THRESHOLDS" << std::endl;
1930   for (size_t ipmt = 0; ipmt < MbdDefs::MBD_N_PMT; ipmt++)
1931   {
1932     // store in a CDBTree
1933     cdbttree->SetFloatValue(ipmt, "thresh_mean", _thresh_mean[ipmt]);
1934     cdbttree->SetFloatValue(ipmt, "thresh_meanerr", _thresh_meanerr[ipmt]);
1935     cdbttree->SetFloatValue(ipmt, "thresh_width", _thresh_width[ipmt]);
1936     cdbttree->SetFloatValue(ipmt, "thresh_widtherr", _thresh_widtherr[ipmt]);
1937     cdbttree->SetFloatValue(ipmt, "thresh_eff", _thresh_eff[ipmt]);
1938     cdbttree->SetFloatValue(ipmt, "thresh_efferr", _thresh_efferr[ipmt]);
1939     cdbttree->SetFloatValue(ipmt, "thresh_chi2ndf", _thresh_chi2ndf[ipmt]);
1940 
1941     if (ipmt < 5 || ipmt >= MbdDefs::MBD_N_PMT - 5)
1942     {
1943       std::cout << ipmt << "\t" << cdbttree->GetFloatValue(ipmt, "thresh_mpv") << std::endl;
1944     }
1945   }
1946 
1947   cdbttree->Commit();
1948   // cdbttree->Print();
1949 
1950   // for now we create the tree after reading it
1951   cdbttree->WriteCDBTTree();
1952   delete cdbttree;
1953 
1954   return 1;
1955 }
1956 #endif
1957 
1958 int MbdCalib::Write_Thresholds(const std::string& dbfile)
1959 {
1960   std::ofstream cal_thresh_file;
1961   cal_thresh_file.open(dbfile);
1962   for (int ipmtch = 0; ipmtch < MbdDefs::MBD_N_PMT; ipmtch++)
1963   {
1964     cal_thresh_file << ipmtch << "\t" << _thresh_mean[ipmtch] << "\t" << _thresh_meanerr[ipmtch]
1965       << "\t" << _thresh_width[ipmtch] << "\t" << _thresh_widtherr[ipmtch]
1966       << "\t" << _thresh_eff[ipmtch] << "\t" << _thresh_efferr[ipmtch]
1967       << "\t" << _thresh_chi2ndf[ipmtch] << std::endl;
1968   }
1969   cal_thresh_file.close();
1970 
1971   return 1;
1972 }
1973 
1974 #ifndef ONLINE
1975 int MbdCalib::Write_CDB_All()
1976 {
1977   return 1;
1978 }
1979 #endif
1980 
1981 // dz is what we need to move the MBD z by
1982 // dt is what we change the MBD t0 by
1983 void MbdCalib::Update_TQT0(const float dz, const float dt)
1984 {
1985   const float z_dt = dz/MbdDefs::C;
1986 
1987   for (int ipmt=0; ipmt<MbdDefs::MBD_N_PMT; ipmt++)
1988   {
1989     // update zvtx
1990     if ( ipmt<64 )  // south
1991     {
1992       _tqfit_t0mean[ipmt] -= z_dt;
1993     }
1994     else
1995     {
1996       _tqfit_t0mean[ipmt] += z_dt;
1997     }
1998 
1999     // update t0
2000     _tqfit_t0mean[ipmt] -= dt;
2001   }
2002 }
2003 
2004 void MbdCalib::Update_TTT0(const float dz, const float dt)
2005 {
2006   // dz is what we need to move the MBD z by
2007   const float z_dt = dz/MbdDefs::C;
2008 
2009   for (int ipmt=0; ipmt<MbdDefs::MBD_N_PMT; ipmt++)
2010   {
2011     if ( ipmt<64 )  // south
2012     {
2013       _ttfit_t0mean[ipmt] -= z_dt;
2014     }
2015     else
2016     {
2017       _ttfit_t0mean[ipmt] += z_dt;
2018     }
2019 
2020     // update t0
2021     _ttfit_t0mean[ipmt] -= dt;
2022   }
2023 }
2024 
2025 void MbdCalib::Reset_TTT0()
2026 {
2027   _ttfit_t0mean.fill( 0. );
2028   _ttfit_t0meanerr.fill( 0. );
2029   _ttfit_t0sigma.fill(std::numeric_limits<float>::quiet_NaN());
2030   _ttfit_t0sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
2031 }
2032 
2033 void MbdCalib::Reset_TQT0()
2034 {
2035   _tqfit_t0mean.fill( 0. );
2036   _tqfit_t0meanerr.fill( 0. );
2037   _tqfit_t0sigma.fill(std::numeric_limits<float>::quiet_NaN());
2038   _tqfit_t0sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
2039 }
2040 
2041 void MbdCalib::Reset_T0Corr()
2042 {
2043   _t0corrmean = 0.;
2044   _t0corrmeanerr = 0.;
2045   _t0corr_fitmean.fill(std::numeric_limits<float>::quiet_NaN());
2046   _t0corr_fitmeanerr.fill(std::numeric_limits<float>::quiet_NaN());
2047   _t0corr_fitsigma.fill(std::numeric_limits<float>::quiet_NaN());
2048   _t0corr_fitsigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
2049   _t0corr_hmean.fill(std::numeric_limits<float>::quiet_NaN());
2050   _t0corr_hmeanerr.fill(std::numeric_limits<float>::quiet_NaN());
2051   _t0corr_hstddev.fill(std::numeric_limits<float>::quiet_NaN());
2052   _t0corr_hstddeverr.fill(std::numeric_limits<float>::quiet_NaN());
2053 }
2054 
2055 void MbdCalib::Reset_Ped()
2056 {
2057   _tqfit_t0mean.fill( 0. );
2058   _tqfit_t0meanerr.fill( 0. );
2059   _tqfit_t0sigma.fill(std::numeric_limits<float>::quiet_NaN());
2060   _tqfit_t0sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
2061 }
2062 
2063 void MbdCalib::Reset_Gains()
2064 {
2065   // Set all initial values
2066   _qfit_integ.fill(std::numeric_limits<float>::quiet_NaN());
2067   _qfit_mpv.fill( 1.0 );
2068   _qfit_sigma.fill(std::numeric_limits<float>::quiet_NaN());
2069   _qfit_integerr.fill(std::numeric_limits<float>::quiet_NaN());
2070   _qfit_mpverr.fill(std::numeric_limits<float>::quiet_NaN());
2071   _qfit_sigmaerr.fill(std::numeric_limits<float>::quiet_NaN());
2072   _qfit_chi2ndf.fill(std::numeric_limits<float>::quiet_NaN());
2073 }
2074 
2075 void MbdCalib::Reset_Pileup()
2076 {
2077   // Set all initial values
2078   _pileup_p0.fill(std::numeric_limits<float>::quiet_NaN());
2079   _pileup_p1.fill(std::numeric_limits<float>::quiet_NaN());
2080   _pileup_p2.fill(std::numeric_limits<float>::quiet_NaN());
2081   _pileup_p0err.fill(std::numeric_limits<float>::quiet_NaN());
2082   _pileup_p1err.fill(std::numeric_limits<float>::quiet_NaN());
2083   _pileup_p2err.fill(std::numeric_limits<float>::quiet_NaN());
2084   _qfit_chi2ndf.fill(std::numeric_limits<float>::quiet_NaN());
2085 }
2086 
2087 void MbdCalib::Reset_Thresholds()
2088 {
2089   // Set all initial values
2090   _thresh_mean.fill(std::numeric_limits<float>::quiet_NaN());
2091   _thresh_meanerr.fill(std::numeric_limits<float>::quiet_NaN());
2092   _thresh_width.fill(std::numeric_limits<float>::quiet_NaN());
2093   _thresh_widtherr.fill(std::numeric_limits<float>::quiet_NaN());
2094   _thresh_eff.fill(std::numeric_limits<float>::quiet_NaN());
2095   _thresh_efferr.fill(std::numeric_limits<float>::quiet_NaN());
2096   _thresh_chi2ndf.fill(std::numeric_limits<float>::quiet_NaN());
2097 }
2098 
2099 void MbdCalib::Reset()
2100 {
2101   Reset_TTT0();
2102   Reset_TQT0();
2103   Reset_Ped();
2104   Reset_Gains();
2105   Reset_T0Corr();
2106   Reset_Pileup();
2107   Reset_Thresholds();
2108 
2109   _sampmax.fill(-1);
2110 }
2111 
2112 void MbdCalib::set_ped(const int ifeech, const float m, const float merr, const float s, const float serr)
2113 {
2114   _pedmean[ifeech] = m;
2115   _pedmeanerr[ifeech] = merr;
2116   _pedsigma[ifeech] = s;
2117   _pedsigmaerr[ifeech] = serr;
2118 }
2119 
2120 
2121 float MbdCalib::get_threshold(const int pmtch, const int rel_or_abs)
2122 {
2123   if ( rel_or_abs==0 )
2124   {
2125     return _thresh_mean[pmtch]/_qfit_mpv[pmtch];
2126   }
2127  
2128   return _thresh_mean[pmtch];
2129 }