Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:16:42

0001 #ifndef MBD_MBDCALIB_H
0002 #define MBD_MBDCALIB_H
0003 
0004 #include "MbdDefs.h"
0005 #include "MbdGeom.h"
0006 
0007 #ifndef ONLINE
0008 #include <fun4all/Fun4AllBase.h>
0009 #include <phool/recoConsts.h>
0010 #endif
0011 
0012 
0013 #include <array>
0014 #include <vector>
0015 #include <string>
0016 #include <memory>
0017 
0018 class TTree;
0019 class CDBInterface;
0020 
0021 class MbdCalib 
0022 {
0023  public:
0024   MbdCalib();
0025 
0026   // MbdCalib(MbdCalib &other) = delete;
0027   // void operator=(const MbdCalib &) = delete;
0028 
0029   virtual ~MbdCalib() {}
0030 
0031   float get_qgain(const int ipmt) const { return _qfit_mpv[ipmt]; }
0032   float get_tt0(const int ipmt) const { return _ttfit_t0mean[ipmt]; }
0033   float get_tq0(const int ipmt) const { return _tqfit_t0mean[ipmt]; }
0034   float get_t0corr() const { return _t0corrmean; }
0035   float get_ped(const int ifeech) const { return _pedmean[ifeech]; }
0036   float get_pedrms(const int ifeech) const { return _pedsigma[ifeech]; }
0037   int   get_sampmax(const int ifeech) const { return _sampmax[ifeech]; }
0038   float get_tcorr(const int ifeech, const int tdc) const {
0039     if (tdc<0)
0040     {
0041       //std::cout << "bad tdc " << ifeech << " " << tdc << std::endl;
0042       return _tcorr_y_interp[ifeech][0];
0043     }
0044     if (tdc>=_tcorr_maxrange[ifeech])
0045     {
0046       //std::cout << "bad tdc " << ifeech << " " << tdc << " " << _tcorr_maxrange[ifeech] << std::endl;
0047       return _tcorr_y_interp[ifeech][_tcorr_maxrange[ifeech]-1];
0048     }
0049     //std::cout << "aaa " << _tcorr_y_interp[ifeech].size() << std::endl;
0050     return _tcorr_y_interp[ifeech][tdc];
0051   }
0052   float get_scorr(const int ifeech, const int adc) const {
0053     if ( _scorr_y_interp[ifeech].size() == 0 ) return 0.; // return 0 if calib doesn't exist
0054     if (adc<0)
0055     {
0056       //std::cout << "bad adc " << ifeech << " " << adc << std::endl;
0057       return _scorr_y_interp[ifeech][0];
0058     }
0059     if (adc>=_scorr_maxrange[ifeech])
0060     {
0061       //std::cout << "high adc " << ifeech << " " << adc << std::endl;
0062       return _scorr_y_interp[ifeech][_scorr_maxrange[ifeech]-1];
0063     }
0064     return _scorr_y_interp[ifeech][adc];
0065   }
0066 
0067   std::vector<float> get_shape(const int ifeech) const { return _shape_y[ifeech]; }
0068   std::vector<float> get_sherr(const int ifeech) const { return _sherr_yerr[ifeech]; }
0069 
0070   float get_pileup(const int ifeech, const int ipar) const {
0071 
0072     if (ipar==0)
0073     {
0074       return _pileup_p0[ifeech];
0075     }
0076     else if (ipar==1)
0077     {
0078       return _pileup_p1[ifeech];
0079     }
0080     else if (ipar==2)
0081     {
0082       return _pileup_p2[ifeech];
0083     }
0084 
0085     return std::numeric_limits<float>::quiet_NaN();
0086   }
0087 
0088   float get_threshold(const int pmtch, const int rel_or_abs = 0);
0089 
0090   void set_sampmax(const int ifeech, const int val) { _sampmax[ifeech] = val; }
0091   void set_ped(const int ifeech, const float m, const float merr, const float s, const float serr);
0092   void set_tt0(const int ipmt, const float t0) { _ttfit_t0mean[ipmt] = t0; }
0093   void set_tq0(const int ipmt, const float t0) { _tqfit_t0mean[ipmt] = t0; }
0094 
0095   int Download_Gains(const std::string& dbase_location);
0096   int Download_TQT0(const std::string& dbase_location);
0097   int Download_TTT0(const std::string& dbase_location);
0098   int Download_T0Corr(const std::string& dbase_location);
0099   int Download_Ped(const std::string& dbase_location);
0100   int Download_SampMax(const std::string& dbase_location);
0101   int Download_Shapes(const std::string& dbase_location);
0102   int Download_TimeCorr(const std::string& dbase_location);
0103   int Download_SlewCorr(const std::string& dbase_location);
0104   int Download_Pileup(const std::string& dbase_location);
0105   int Download_Thresholds(const std::string& dbase_location);
0106   int Download_All();
0107 
0108 #ifndef ONLINE
0109   int Write_CDB_SampMax(const std::string& dbfile);
0110   int Write_CDB_TTT0(const std::string& dbfile);
0111   int Write_CDB_TQT0(const std::string& dbfile);
0112   int Write_CDB_T0Corr(const std::string& dbfile);
0113   int Write_CDB_Ped(const std::string& dbfile);
0114   int Write_CDB_Shapes(const std::string& dbfile);
0115   int Write_CDB_TimeCorr(const std::string& dbfile);
0116   int Write_CDB_SlewCorr(const std::string& dbfile);
0117   int Write_CDB_Gains(const std::string& dbfile);
0118   int Write_CDB_Pileup(const std::string& dbfile);
0119   int Write_CDB_Thresholds(const std::string& dbfile);
0120   static int Write_CDB_All();
0121 #endif
0122 
0123   int Write_SampMax(const std::string& dbfile);
0124   int Write_TQT0(const std::string& dbfile);
0125   int Write_TTT0(const std::string& dbfile);
0126   int Write_T0Corr(const std::string& dbfile);
0127   int Write_Ped(const std::string& dbfile);
0128   int Write_Gains(const std::string& dbfile);
0129   int Write_Pileup(const std::string& dbfile);
0130   int Write_Thresholds(const std::string& dbfile);
0131 
0132   void Reset_TQT0();
0133   void Reset_TTT0();
0134   void Reset_T0Corr();
0135   void Reset_Ped();
0136   void Reset_Gains();
0137   void Reset_Pileup();
0138   void Reset_Thresholds();
0139 
0140   void Update_TQT0(const float dz, const float dt = 0.); // update with new z-vertex, t0
0141   void Update_TTT0(const float dz, const float dt = 0.);
0142 
0143   // void Dump_to_file(const std::string& what = "ALL");
0144 
0145   void Reset();
0146   // void Print(Option_t* option) const;
0147 
0148   int  Verbosity()            { return _verbose; }
0149   void Verbosity(const int v) { _verbose = v; }
0150 
0151  private:
0152 #ifndef ONLINE
0153   CDBInterface* _cdb{nullptr};
0154   recoConsts* _rc{nullptr};
0155 #endif
0156 
0157   std::unique_ptr<MbdGeom> _mbdgeom{nullptr};
0158 
0159   int _status{0};
0160   int _verbose{0};
0161   // int          _run_number {0};
0162   // uint64_t     _timestamp {0};
0163   std::string _dbfilename;
0164 
0165   // Assumes Landau fit
0166   std::array<float, MbdDefs::MBD_N_PMT> _qfit_integ{};
0167   std::array<float, MbdDefs::MBD_N_PMT> _qfit_mpv{};
0168   std::array<float, MbdDefs::MBD_N_PMT> _qfit_sigma{};
0169   std::array<float, MbdDefs::MBD_N_PMT> _qfit_integerr{};
0170   std::array<float, MbdDefs::MBD_N_PMT> _qfit_mpverr{};
0171   std::array<float, MbdDefs::MBD_N_PMT> _qfit_sigmaerr{};
0172   std::array<float, MbdDefs::MBD_N_PMT> _qfit_chi2ndf{};
0173 
0174   // T0 offsets, time channels
0175   std::array<float, MbdDefs::MBD_N_PMT> _ttfit_t0mean{};
0176   std::array<float, MbdDefs::MBD_N_PMT> _ttfit_t0meanerr{};
0177   std::array<float, MbdDefs::MBD_N_PMT> _ttfit_t0sigma{};
0178   std::array<float, MbdDefs::MBD_N_PMT> _ttfit_t0sigmaerr{};
0179 
0180   // T0 offsets, charge channels
0181   std::array<float, MbdDefs::MBD_N_PMT> _tqfit_t0mean{};
0182   std::array<float, MbdDefs::MBD_N_PMT> _tqfit_t0meanerr{};
0183   std::array<float, MbdDefs::MBD_N_PMT> _tqfit_t0sigma{};
0184   std::array<float, MbdDefs::MBD_N_PMT> _tqfit_t0sigmaerr{};
0185 
0186   // Overall T0 offset
0187   Float_t _t0corrmean{ 0. };
0188   Float_t _t0corrmeanerr{ 0. };
0189   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_fitmean{};
0190   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_fitmeanerr{};
0191   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_fitsigma{};
0192   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_fitsigmaerr{};
0193   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_hmean{};
0194   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_hmeanerr{};
0195   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_hstddev{};
0196   std::array<float, MbdDefs::MBD_N_ARMS> _t0corr_hstddeverr{};
0197 
0198   // Pedestals
0199   std::array<float, MbdDefs::MBD_N_FEECH> _pedmean{};
0200   std::array<float, MbdDefs::MBD_N_FEECH> _pedmeanerr{};
0201   std::array<float, MbdDefs::MBD_N_FEECH> _pedsigma{};
0202   std::array<float, MbdDefs::MBD_N_FEECH> _pedsigmaerr{};
0203 
0204   // SampMax (Peak of waveform)
0205   std::array<int, MbdDefs::MBD_N_FEECH> _sampmax{};
0206 
0207   // Pileup waveform correction
0208   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_p0{};
0209   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_p0err{};
0210   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_p1{};
0211   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_p1err{};
0212   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_p2{};
0213   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_p2err{};
0214   std::array<float, MbdDefs::MBD_N_FEECH> _pileup_chi2ndf{};
0215 
0216   // Waveform Template
0217   int do_templatefit{0};
0218   std::array<int, MbdDefs::MBD_N_FEECH>   _shape_npts{};      // num points in template
0219   std::array<float, MbdDefs::MBD_N_FEECH> _shape_minrange{};  // in template units (samples)
0220   std::array<float, MbdDefs::MBD_N_FEECH> _shape_maxrange{};  // in template units (samples)
0221   std::array<std::vector<float>, MbdDefs::MBD_N_FEECH> _shape_y{};
0222 
0223   std::array<int, MbdDefs::MBD_N_FEECH>   _sherr_npts{};      // num points in template
0224   std::array<float, MbdDefs::MBD_N_FEECH> _sherr_minrange{};  // in template units (samples)
0225   std::array<float, MbdDefs::MBD_N_FEECH> _sherr_maxrange{};  // in template units (samples)
0226   std::array<std::vector<float>, MbdDefs::MBD_N_FEECH> _sherr_yerr{};
0227 
0228   // Fine Timing Corrections
0229   std::array<int, MbdDefs::MBD_N_FEECH>   _tcorr_npts{};      // num points in timing corr
0230   std::array<float, MbdDefs::MBD_N_FEECH> _tcorr_minrange{};  // in units (delta-TDC)
0231   std::array<float, MbdDefs::MBD_N_FEECH> _tcorr_maxrange{};  // in units (detta-TDC)
0232   std::array<std::vector<float>, MbdDefs::MBD_N_FEECH> _tcorr_y{};
0233   std::array<std::vector<float>, MbdDefs::MBD_N_FEECH> _tcorr_y_interp{}; // interpolated tcorr
0234 
0235   // Slew Correction
0236   std::array<int, MbdDefs::MBD_N_FEECH>   _scorr_npts{};      // num points in slew corr
0237   std::array<float, MbdDefs::MBD_N_FEECH> _scorr_minrange{};  // in units (delta-TDC)
0238   std::array<float, MbdDefs::MBD_N_FEECH> _scorr_maxrange{};  // in units (detta-TDC)
0239   std::array<std::vector<float>, MbdDefs::MBD_N_FEECH> _scorr_y{};
0240   std::array<std::vector<float>, MbdDefs::MBD_N_FEECH> _scorr_y_interp{}; // interpolated scorr
0241 
0242   // Thresolds
0243   std::array<float, MbdDefs::MBD_N_PMT> _thresh_mean{};
0244   std::array<float, MbdDefs::MBD_N_PMT> _thresh_meanerr{};
0245   std::array<float, MbdDefs::MBD_N_PMT> _thresh_width{};
0246   std::array<float, MbdDefs::MBD_N_PMT> _thresh_widtherr{};
0247   std::array<float, MbdDefs::MBD_N_PMT> _thresh_eff{};
0248   std::array<float, MbdDefs::MBD_N_PMT> _thresh_efferr{};
0249   std::array<float, MbdDefs::MBD_N_PMT> _thresh_chi2ndf{};
0250 };
0251 
0252 #endif  // MBD_MBDCALIB_H