Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:51

0001 #define MELaserPrim_cxx
0002 #include <cassert>
0003 #include "CalibCalorimetry/EcalLaserAnalyzer/interface/MELaserPrim.h"
0004 #include "CalibCalorimetry/EcalLaserAnalyzer/interface/MEGeom.h"
0005 #include <cassert>
0006 #include <cstdlib>
0007 
0008 const TString MELaserPrim::apdpn_arrayName[MELaserPrim::iSizeArray_apdpn] = {
0009     "APD", "APDoPN", "APDoPNA", "APDoPNB", "APDoAPD", "APDoAPDA", "APDoAPDB", "Time"};
0010 const TString MELaserPrim::apdpn_varName[MELaserPrim::iSize_apdpn] = {"Mean", "RMS", "M3", "Nevt", "Min", "Max"};
0011 const TString MELaserPrim::apdpn_varUnit[MELaserPrim::iSizeArray_apdpn][MELaserPrim::iSize_apdpn] =
0012 
0013     {{" (ADC Counts)", " (ADC Counts)", " (ADC Counts)", "", " (ADC Counts)", " (ADC Counts)"},
0014      {"", "", "", "", "", ""},
0015      {"", "", "", "", "", ""},
0016      {"", "", "", "", "", ""},
0017      {"", "", "", "", "", ""},
0018      {"", "", "", "", "", ""},
0019      {" (25 ns)", " (25 ns)", " (25 ns)", "", " (25 ns)", " (25 ns)"}};
0020 const TString MELaserPrim::apdpn_extraVarName[MELaserPrim::iSizeExtra_apdpn] = {"ShapeCor"};
0021 const TString MELaserPrim::apdpn_extraVarUnit[MELaserPrim::iSizeExtra_apdpn] = {""};
0022 const TString MELaserPrim::ab_varName[MELaserPrim::iSize_ab] = {"alpha", "beta", "width", "chi2"};
0023 const TString MELaserPrim::mtq_varName[MELaserPrim::iSize_mtq] = {
0024     "peak", "sigma", "fit", "ampl", "trise", "fwhm", "fw20", "fw80", "sliding"};
0025 const TString MELaserPrim::mtq_varUnit[MELaserPrim::iSize_mtq] = {"(nanoseconds)",
0026                                                                   "(nanoseconds)",
0027                                                                   "(nanoseconds)",
0028                                                                   "(ADC counts)",
0029                                                                   "(nanoseconds)",
0030                                                                   "(nanoseconds)",
0031                                                                   "(nanoseconds)",
0032                                                                   "(nanoseconds)",
0033                                                                   "(ADC counts)"};
0034 
0035 const TString MELaserPrim::separator = "__";
0036 
0037 //ClassImp( MELaserPrim )
0038 
0039 MELaserPrim::MELaserPrim(ME::Header header, ME::Settings settings, const char* inpath, const char* outfile)
0040     : init_ok(false), _isBarrel(true), _inpath(inpath), _outfile(outfile) {
0041   apdpn_file = nullptr;
0042   ab_file = nullptr;
0043   mtq_file = nullptr;
0044   tpapd_file = nullptr;
0045   apdpn_tree = nullptr;
0046   ab_tree = nullptr;
0047   pn_tree = nullptr;
0048   mtq_tree = nullptr;
0049   tpapd_tree = nullptr;
0050   tppn_tree = nullptr;
0051   ixmin = 0;
0052   ixmax = 0;
0053   iymin = 0;
0054   iymax = 0;
0055 
0056   _dcc = header.dcc;
0057   _side = header.side;
0058   _run = header.run;
0059   _lb = header.lb;
0060   _events = header.events;
0061   _ts = header.ts_beg;
0062   _ts_beg = header.ts_beg;
0063   _ts_end = header.ts_end;
0064 
0065   _type = settings.type;
0066   _color = settings.wavelength;
0067   _power = settings.power;
0068   _filter = settings.filter;
0069   _delay = settings.delay;
0070   _mgpagain = settings.mgpagain;
0071   _memgain = settings.memgain;
0072 
0073   if (_type == ME::iLaser) {
0074     _primStr = lmfLaserName(ME::iLmfLaserPrim, _type, _color) + separator;
0075     _pnPrimStr = lmfLaserName(ME::iLmfLaserPnPrim, _type, _color) + separator;
0076     _pulseStr = lmfLaserName(ME::iLmfLaserPulse, _type, _color) + separator;
0077   } else if (_type == ME::iTestPulse) {
0078     _tpPrimStr = lmfLaserName(ME::iLmfTestPulsePrim, _type) + separator;
0079     _tpPnPrimStr = lmfLaserName(ME::iLmfTestPulsePnPrim, _type) + separator;
0080   }
0081 
0082   _lmr = ME::lmr(_dcc, _side);
0083   ME::regionAndSector(_lmr, _reg, _sm, _dcc, _side);
0084   _isBarrel = (_reg == ME::iEBM || _reg == ME::iEBP);
0085   _sectorStr = ME::smName(_lmr);
0086   _regionStr = _sectorStr;
0087   _regionStr += "_";
0088   _regionStr += _side;
0089 
0090   init();
0091   bookHistograms();
0092   //fillHistograms();
0093   //writeHistograms();
0094 }
0095 
0096 TString MELaserPrim::channelViewName(int iname) {
0097   switch (iname) {
0098     case iECAL:
0099       return "ECAL";
0100     case iECAL_LMR:
0101       return "ECAL_LMR";
0102     case iEB_crystal_number:
0103       return "EB_crystal_number";
0104     case iEB_LM_LMM:
0105       return "EB_LM_LMM";
0106     case iEB_LM_PN:
0107       return "EB_LM_PN";
0108     case iEE_crystal_number:
0109       return "EE_crystal_number";
0110     case iEE_LM_LMM:
0111       return "EE_LM_LMM";
0112     case iEE_LM_PN:
0113       return "EE_LM_PN";
0114     default:
0115       abort();
0116   }
0117   return "";
0118 }
0119 
0120 int MELaserPrim::logicId(int channelView, int id1, int id2) {
0121   assert(channelView >= iECAL && channelView < iSize_cv);
0122   return 1000000 * channelView + 10000 * id1 + id2;
0123 }
0124 
0125 bool MELaserPrim::getViewIds(int logic_id, int& channelView, int& id1, int& id2) {
0126   bool out = true;
0127   int channelView_ = logic_id / 1000000;
0128   if (channelView != 0 && channelView_ != channelView)
0129     out = false;
0130   channelView = channelView_;
0131   id1 = (logic_id % 1000000) / 10000;
0132   id2 = logic_id % 10000;
0133   return out;
0134 }
0135 
0136 void MELaserPrim::init() {
0137   bool verbose_ = false;
0138 
0139   if (_inpath == "0") {
0140     if (verbose_)
0141       std::cout << "no input file" << std::endl;
0142     init_ok = true;
0143     return;  // GHM
0144   }
0145 
0146   TString cur(_inpath);
0147   if (!cur.EndsWith("/"))
0148     cur += "/";
0149 
0150   if (_type == ME::iLaser) {
0151     TString _APDPN_fname = cur;
0152     _APDPN_fname += "APDPN_LASER.root";
0153     TString _AB_fname = cur;
0154     _AB_fname += "AB.root";
0155     TString _MTQ_fname = cur;
0156     _MTQ_fname += "MATACQ.root";
0157 
0158     bool apdpn_ok, ab_ok, pn_ok, mtq_ok;
0159     apdpn_ok = false;
0160     ab_ok = false;
0161     pn_ok = false;
0162     mtq_ok = false;
0163 
0164     FILE* test;
0165     test = fopen(_APDPN_fname, "r");
0166     if (test) {
0167       apdpn_ok = true;
0168       pn_ok = true;
0169       fclose(test);
0170     }
0171     test = fopen(_AB_fname, "r");
0172     if (test) {
0173       ab_ok = true;
0174       fclose(test);
0175     }
0176     test = fopen(_MTQ_fname, "r");
0177     if (test) {
0178       mtq_ok = true;
0179       fclose(test);
0180     }
0181 
0182     if (apdpn_ok)
0183       apdpn_file = TFile::Open(_APDPN_fname);
0184     if (ab_ok)
0185       ab_file = TFile::Open(_AB_fname);
0186     if (mtq_ok)
0187       mtq_file = TFile::Open(_MTQ_fname);
0188 
0189     if (verbose_) {
0190       std::cout << _APDPN_fname << " ok=" << apdpn_ok << std::endl;
0191       std::cout << _AB_fname << " ok=" << ab_ok << std::endl;
0192       std::cout << _MTQ_fname << " ok=" << mtq_ok << std::endl;
0193     }
0194     if (!apdpn_ok || !pn_ok)
0195       return;  // FIXME !
0196 
0197     TString apdpn_tree_name;
0198     TString ab_tree_name;
0199     TString pn_tree_name;
0200     TString mtq_tree_name;
0201 
0202     apdpn_tree_name = "APDCol";
0203     ab_tree_name = "ABCol";
0204     pn_tree_name = "PNCol";
0205     mtq_tree_name = "MatacqCol";
0206 
0207     apdpn_tree_name += _color;
0208     ab_tree_name += _color;
0209     pn_tree_name += _color;
0210     mtq_tree_name += _color;
0211 
0212     if (mtq_ok) {
0213       TTree* ckeckMtq = (TTree*)mtq_file->Get(mtq_tree_name);
0214       if (ckeckMtq == nullptr)
0215         mtq_ok = false;
0216     }
0217 
0218     if (_color != ME::iIRed && _color != ME::iBlue) {
0219       std::cout << "MELaserPrim::init() -- Fatal Error -- Wrong Laser Color : " << _color << " ---- Abort "
0220                 << std::endl;
0221       return;
0222     }
0223 
0224     apdpn_tree = (TTree*)apdpn_file->Get(apdpn_tree_name);
0225     assert(apdpn_tree != nullptr);
0226     apdpn_tree->SetMakeClass(1);
0227     apdpn_tree->SetBranchAddress("dccID", &apdpn_dccID, &b_apdpn_dccID);
0228     apdpn_tree->SetBranchAddress("towerID", &apdpn_towerID, &b_apdpn_towerID);
0229     apdpn_tree->SetBranchAddress("channelID", &apdpn_channelID, &b_apdpn_channelID);
0230     apdpn_tree->SetBranchAddress("moduleID", &apdpn_moduleID, &b_apdpn_moduleID);
0231     apdpn_tree->SetBranchAddress("side", &apdpn_side, &b_apdpn_side);
0232     apdpn_tree->SetBranchAddress("ieta", &apdpn_ieta, &b_apdpn_ieta);
0233     apdpn_tree->SetBranchAddress("iphi", &apdpn_iphi, &b_apdpn_iphi);
0234     apdpn_tree->SetBranchAddress("flag", &apdpn_flag, &b_apdpn_flag);
0235     if (apdpn_tree->GetBranchStatus("ShapeCor"))
0236       apdpn_tree->SetBranchAddress("ShapeCor", &apdpn_ShapeCor, &b_apdpn_ShapeCor);
0237     else
0238       apdpn_ShapeCor = 0.0;
0239     for (int jj = 0; jj < iSizeArray_apdpn; jj++) {
0240       TString name_ = apdpn_arrayName[jj];
0241       apdpn_tree->SetBranchAddress(name_, apdpn_apdpn[jj], &b_apdpn_apdpn[jj]);
0242     }
0243 
0244     if (ab_ok) {
0245       ab_tree = (TTree*)ab_file->Get(ab_tree_name);
0246       assert(ab_tree != nullptr);
0247       ab_tree->SetMakeClass(1);
0248       ab_tree->SetBranchAddress("dccID", &ab_dccID, &b_ab_dccID);
0249       ab_tree->SetBranchAddress("towerID", &ab_towerID, &b_ab_towerID);
0250       ab_tree->SetBranchAddress("channelID", &ab_channelID, &b_ab_channelID);
0251       ab_tree->SetBranchAddress("ieta", &ab_ieta, &b_ab_ieta);
0252       ab_tree->SetBranchAddress("iphi", &ab_iphi, &b_ab_iphi);
0253       ab_tree->SetBranchAddress("flag", &ab_flag, &b_ab_flag);
0254       for (int ii = 0; ii < iSize_ab; ii++) {
0255         ab_tree->SetBranchAddress(ab_varName[ii], &ab_ab[ii], &b_ab_ab[ii]);
0256       }
0257     }
0258 
0259     pn_tree = (TTree*)apdpn_file->Get(pn_tree_name);
0260     assert(pn_tree != nullptr);
0261     pn_tree->SetMakeClass(1);
0262     pn_tree->SetBranchAddress("side", &pn_side, &b_pn_side);
0263     pn_tree->SetBranchAddress("pnID", &pn_pnID, &b_pn_pnID);
0264     pn_tree->SetBranchAddress("moduleID", &pn_moduleID, &b_pn_moduleID);
0265     pn_tree->SetBranchAddress("PN", pn_PN, &b_pn_PN);
0266     pn_tree->SetBranchAddress("PNoPN", pn_PNoPN, &b_pn_PNoPN);
0267     pn_tree->SetBranchAddress("PNoPNA", pn_PNoPNA, &b_pn_PNoPNA);
0268     pn_tree->SetBranchAddress("PNoPNB", pn_PNoPNB, &b_pn_PNoPNB);
0269 
0270     if (mtq_ok) {
0271       mtq_tree = (TTree*)mtq_file->Get(mtq_tree_name);
0272       assert(mtq_tree != nullptr);
0273       mtq_tree->SetMakeClass(1);
0274       mtq_tree->SetBranchAddress("side", &mtq_side, &b_mtq_side);
0275 
0276       for (int ii = 0; ii < iSize_mtq; ii++) {
0277         mtq_tree->SetBranchAddress(mtq_varName[ii], &mtq_mtq[ii], &b_mtq_mtq[ii]);
0278       }
0279     }
0280   } else if (_type == ME::iTestPulse) {
0281     TString _TPAPD_fname = cur;
0282     _TPAPD_fname += "APDPN_TESTPULSE.root";
0283 
0284     bool tpapd_ok;
0285     tpapd_ok = false;
0286 
0287     FILE* test;
0288     test = fopen(_TPAPD_fname, "r");
0289     if (test) {
0290       tpapd_ok = true;
0291       fclose(test);
0292     }
0293 
0294     if (tpapd_ok)
0295       tpapd_file = TFile::Open(_TPAPD_fname);
0296 
0297     if (verbose_) {
0298       std::cout << _TPAPD_fname << " ok=" << tpapd_ok << std::endl;
0299     }
0300     if (!tpapd_ok)
0301       return;
0302 
0303     TString tpapd_tree_name;
0304     TString tppn_tree_name;
0305 
0306     tpapd_tree_name = "TPAPD";
0307     tppn_tree_name = "TPPN";
0308 
0309     tpapd_tree = (TTree*)tpapd_file->Get(tpapd_tree_name);
0310     assert(tpapd_tree != nullptr);
0311     tpapd_tree->SetMakeClass(1);
0312     tpapd_tree->SetBranchAddress("ieta", &tpapd_ieta, &b_tpapd_ieta);
0313     tpapd_tree->SetBranchAddress("iphi", &tpapd_iphi, &b_tpapd_iphi);
0314     tpapd_tree->SetBranchAddress("dccID", &tpapd_dccID, &b_tpapd_dccID);
0315     tpapd_tree->SetBranchAddress("side", &tpapd_side, &b_tpapd_side);
0316     tpapd_tree->SetBranchAddress("towerID", &tpapd_towerID, &b_tpapd_towerID);
0317     tpapd_tree->SetBranchAddress("channelID", &tpapd_channelID, &b_tpapd_channelID);
0318     tpapd_tree->SetBranchAddress("moduleID", &tpapd_moduleID, &b_tpapd_moduleID);
0319     tpapd_tree->SetBranchAddress("flag", &tpapd_flag, &b_tpapd_flag);
0320     tpapd_tree->SetBranchAddress("gain", &tpapd_gain, &b_tpapd_gain);
0321     tpapd_tree->SetBranchAddress("APD", tpapd_APD, &b_tpapd_APD);
0322 
0323     tppn_tree = (TTree*)tpapd_file->Get(tppn_tree_name);
0324     assert(tppn_tree != nullptr);
0325     tppn_tree->SetMakeClass(1);
0326     tppn_tree->SetBranchAddress("side", &tppn_side, &b_tppn_side);
0327     tppn_tree->SetBranchAddress("pnID", &tppn_pnID, &b_tppn_pnID);
0328     tppn_tree->SetBranchAddress("moduleID", &tppn_moduleID, &b_tppn_moduleID);
0329     tppn_tree->SetBranchAddress("gain", &tppn_gain, &b_tppn_gain);
0330     tppn_tree->SetBranchAddress("PN", tppn_PN, &b_tppn_PN);
0331   }
0332   init_ok = true;
0333 }
0334 
0335 void MELaserPrim::bookHistograms() {
0336   if (!init_ok)
0337     return;
0338   refresh();
0339 
0340   TString i_name, d_name;
0341 
0342   if (_isBarrel) {
0343     ixmin = 0;
0344     ixmax = 85;
0345     nx = ixmax - ixmin;
0346     iymin = 0;
0347     iymax = 20;
0348     ny = iymax - iymin;
0349 
0350     //       for( int ilmod=1; ilmod<=9; ilmod++ )
0351     //  {
0352     //    _pn[ilmod] = MEEBGeom::pn( ilmod );
0353     //  }
0354   } else  // fixme --- to be implemented
0355   {
0356     ixmin = 1;
0357     ixmax = 101;
0358     nx = ixmax - ixmin;
0359     iymin = 1;
0360     iymax = 101;
0361     ny = iymax - iymin;
0362     //       for( int ilmod=1; ilmod<=21; ilmod++ )  // modules 20 and 21 are fake...
0363     //  {
0364     //    _pn[ilmod] = MEEEGeom::pn( ilmod );
0365     //  }
0366     //      abort();
0367   }
0368 
0369   TString t_name;
0370 
0371   //
0372   // Laser Run
0373   //
0374   t_name = "LMF_RUN_DAT";
0375   addBranchI(t_name, "LOGIC_ID");
0376   addBranchI(t_name, "NEVENTS");
0377   addBranchI(t_name, "QUALITY_FLAG");
0378 
0379   //
0380   // Laser Run IOV
0381   //
0382   t_name = "LMF_RUN_IOV";
0383   addBranchI(t_name, "TAG_ID");
0384   addBranchI(t_name, "SUB_RUN_NUM");
0385   addBranchI(t_name, "SUB_RUN_START_LOW");
0386   addBranchI(t_name, "SUB_RUN_START_HIGH");
0387   addBranchI(t_name, "SUB_RUN_END_LOW");
0388   addBranchI(t_name, "SUB_RUN_END_HIGH");
0389   addBranchI(t_name, "DB_TIMESTAMP_LOW");
0390   addBranchI(t_name, "DB_TIMESTAMP_HIGH");
0391   addBranchC(t_name, "SUB_RUN_TYPE");
0392 
0393   if (_type == ME::iLaser) {
0394     //
0395     // Laser ADC Primitives
0396     //
0397     bookHistoI(_primStr, "LOGIC_ID");
0398     bookHistoI(_primStr, "FLAG");
0399     bookHistoF(_primStr, "MEAN");
0400     bookHistoF(_primStr, "RMS");
0401     bookHistoF(_primStr, "M3");
0402     bookHistoF(_primStr, "APD_OVER_PNA_MEAN");
0403     bookHistoF(_primStr, "APD_OVER_PNA_RMS");
0404     bookHistoF(_primStr, "APD_OVER_PNA_M3");
0405     bookHistoF(_primStr, "APD_OVER_PNB_MEAN");
0406     bookHistoF(_primStr, "APD_OVER_PNB_RMS");
0407     bookHistoF(_primStr, "APD_OVER_PNB_M3");
0408     bookHistoF(_primStr, "APD_OVER_PN_MEAN");
0409     bookHistoF(_primStr, "APD_OVER_PN_RMS");
0410     bookHistoF(_primStr, "APD_OVER_PN_M3");
0411     bookHistoF(_primStr, "APD_OVER_APDA_MEAN");
0412     bookHistoF(_primStr, "APD_OVER_APDA_RMS");
0413     bookHistoF(_primStr, "APD_OVER_APDA_M3");
0414     bookHistoF(_primStr, "APD_OVER_APDB_MEAN");
0415     bookHistoF(_primStr, "APD_OVER_APDB_RMS");
0416     bookHistoF(_primStr, "APD_OVER_APDB_M3");
0417     bookHistoF(_primStr, "SHAPE_COR");
0418     bookHistoF(_primStr, "ALPHA");
0419     bookHistoF(_primStr, "BETA");
0420     // NEW GHM 08/06 --> SCHEMA MODIFIED?
0421     bookHistoF(_primStr, "TIME_MEAN");
0422     bookHistoF(_primStr, "TIME_RMS");
0423     bookHistoF(_primStr, "TIME_M3");
0424     bookHistoF(_primStr, "TIME_NEVT");
0425 
0426     //
0427     // Laser PN Primitives
0428     //
0429     t_name = lmfLaserName(ME::iLmfLaserPnPrim, _type, _color);
0430     addBranchI(t_name, "LOGIC_ID");
0431     addBranchI(t_name, "FLAG");
0432     addBranchF(t_name, "MEAN");
0433     addBranchF(t_name, "RMS");
0434     addBranchF(t_name, "M3");
0435     addBranchF(t_name, "PNA_OVER_PNB_MEAN");
0436     addBranchF(t_name, "PNA_OVER_PNB_RMS");
0437     addBranchF(t_name, "PNA_OVER_PNB_M3");
0438 
0439     //
0440     // Laser Pulse
0441     //
0442     t_name = lmfLaserName(ME::iLmfLaserPulse, _type, _color);
0443     addBranchI(t_name, "LOGIC_ID");
0444     addBranchI(t_name, "FIT_METHOD");
0445     addBranchF(t_name, "MTQ_AMPL");
0446     addBranchF(t_name, "MTQ_TIME");
0447     addBranchF(t_name, "MTQ_RISE");
0448     addBranchF(t_name, "MTQ_FWHM");
0449     addBranchF(t_name, "MTQ_FW20");
0450     addBranchF(t_name, "MTQ_FW80");
0451     addBranchF(t_name, "MTQ_SLIDING");
0452 
0453     //
0454     // Laser Config
0455     //
0456     t_name = lmfLaserName(ME::iLmfLaserConfig, _type);
0457     addBranchI(t_name, "LOGIC_ID");
0458     addBranchI(t_name, "WAVELENGTH");
0459     addBranchI(t_name, "VFE_GAIN");
0460     addBranchI(t_name, "PN_GAIN");
0461     addBranchI(t_name, "LSR_POWER");
0462     addBranchI(t_name, "LSR_ATTENUATOR");
0463     addBranchI(t_name, "LSR_CURRENT");
0464     addBranchI(t_name, "LSR_DELAY_1");
0465     addBranchI(t_name, "LSR_DELAY_2");
0466 
0467     //
0468     // Laser LaserRun config dat
0469     //
0470     t_name = "RUN_LASERRUN_CONFIG_DAT";
0471     addBranchI(t_name, "LOGIC_ID");
0472     addBranchC(t_name, "LASER_SEQUENCE_TYPE");
0473     addBranchC(t_name, "LASER_SEQUENCE_COND");
0474 
0475   } else if (_type == ME::iTestPulse) {
0476     //
0477     // Test Pulse ADC Primitives
0478     //
0479     bookHistoI(_tpPrimStr, "LOGIC_ID");
0480     bookHistoI(_tpPrimStr, "FLAG");
0481     bookHistoF(_tpPrimStr, "MEAN");
0482     bookHistoF(_tpPrimStr, "RMS");
0483     bookHistoF(_tpPrimStr, "M3");
0484     bookHistoF(_tpPrimStr, "NEVT");
0485 
0486     //
0487     // Test Pulse PN Primitives
0488     //
0489     t_name = lmfLaserName(ME::iLmfTestPulsePnPrim, _type);
0490     addBranchI(t_name, "LOGIC_ID");
0491     addBranchI(t_name, "FLAG");
0492     addBranchI(t_name, "GAIN");
0493     addBranchF(t_name, "MEAN");
0494     addBranchF(t_name, "RMS");
0495     addBranchF(t_name, "M3");
0496 
0497     //
0498     // Test Pulse Config
0499     //
0500     t_name = lmfLaserName(ME::iLmfTestPulseConfig, _type);
0501     addBranchI(t_name, "LOGIC_ID");
0502     addBranchI(t_name, "VFE_GAIN");
0503     addBranchI(t_name, "PN_GAIN");
0504   }
0505 }
0506 
0507 void MELaserPrim::fillHistograms() {
0508   TString t_name;
0509 
0510   if (!init_ok)
0511     return;
0512 
0513   Long64_t nentries = 0;
0514   Long64_t ientry = 0;
0515 
0516   int channelView_(0);
0517   int id1_(0), id2_(0);
0518   int logic_id_(0);
0519 
0520   if (_type == ME::iLaser) {
0521     nentries = apdpn_tree->GetEntriesFast();
0522     for (Long64_t jentry = 0; jentry < nentries; jentry++) {
0523       ientry = apdpn_tree->LoadTree(jentry);
0524       assert(ientry >= 0);
0525       apdpn_tree->GetEntry(jentry);
0526 
0527       if (ab_tree) {
0528         ientry = ab_tree->LoadTree(jentry);
0529         assert(ientry >= 0);
0530         ab_tree->GetEntry(jentry);
0531       }
0532 
0533       if (apdpn_iphi < 0)
0534         continue;
0535 
0536       // fixme remove until coordinated are fine
0537       //if(ab_tree) assert( apdpn_ieta==ab_ieta && apdpn_iphi==ab_iphi );
0538 
0539       int ix(0);
0540       int iy(0);
0541       if (_isBarrel) {
0542         // Barrel, global coordinates
0543         id1_ = _sm;
0544         if (apdpn_side != _side)
0545           continue;
0546         int ieta = apdpn_ieta;
0547         int iphi = apdpn_iphi;
0548         MEEBGeom::XYCoord xy_ = MEEBGeom::localCoord(ieta, iphi);
0549         ix = xy_.first;
0550         iy = xy_.second;
0551         id2_ = MEEBGeom::crystal_channel(ix, iy);
0552         channelView_ = iEB_crystal_number;
0553       } else {
0554         // EndCaps, global coordinates
0555         id1_ = apdpn_iphi;
0556         id2_ = apdpn_ieta;
0557         ix = id1_;
0558         iy = id2_;
0559         channelView_ = iEE_crystal_number;
0560       }
0561 
0562       logic_id_ = logicId(channelView_, id1_, id2_);
0563 
0564       int flag = apdpn_flag;
0565 
0566       setInt("LOGIC_ID", ix, iy, logic_id_);
0567       setInt("FLAG", ix, iy, flag);
0568       setVal("MEAN", ix, iy, apdpn_apdpn[iAPD][iMean]);
0569       setVal("RMS", ix, iy, apdpn_apdpn[iAPD][iRMS]);
0570       setVal("M3", ix, iy, apdpn_apdpn[iAPD][iM3]);  // fixme --- peak?
0571       setVal("APD_OVER_PNA_MEAN", ix, iy, apdpn_apdpn[iAPDoPNA][iMean]);
0572       setVal("APD_OVER_PNA_RMS", ix, iy, apdpn_apdpn[iAPDoPNA][iRMS]);
0573       setVal("APD_OVER_PNA_M3", ix, iy, apdpn_apdpn[iAPDoPNA][iM3]);  // fixme
0574       setVal("APD_OVER_PNB_MEAN", ix, iy, apdpn_apdpn[iAPDoPNB][iMean]);
0575       setVal("APD_OVER_PNB_RMS", ix, iy, apdpn_apdpn[iAPDoPNB][iRMS]);
0576       setVal("APD_OVER_PNB_M3", ix, iy, apdpn_apdpn[iAPDoPNB][iM3]);  // fixme
0577       setVal("APD_OVER_PN_MEAN", ix, iy, apdpn_apdpn[iAPDoPN][iMean]);
0578       setVal("APD_OVER_PN_RMS", ix, iy, apdpn_apdpn[iAPDoPN][iRMS]);
0579       setVal("APD_OVER_PN_M3", ix, iy, apdpn_apdpn[iAPDoPN][iM3]);  // fixme
0580       // JM
0581       setVal("APD_OVER_APD_MEAN", ix, iy, apdpn_apdpn[iAPDoAPDA][iMean]);
0582       setVal("APD_OVER_APD_RMS", ix, iy, apdpn_apdpn[iAPDoAPDA][iRMS]);
0583       setVal("APD_OVER_APD_M3", ix, iy, apdpn_apdpn[iAPDoAPDA][iM3]);  // fixme
0584       setVal("APD_OVER_APDA_MEAN", ix, iy, apdpn_apdpn[iAPDoAPDA][iMean]);
0585       setVal("APD_OVER_APDA_RMS", ix, iy, apdpn_apdpn[iAPDoAPDA][iRMS]);
0586       setVal("APD_OVER_APDA_M3", ix, iy, apdpn_apdpn[iAPDoAPDA][iM3]);  // fixme
0587       setVal("APD_OVER_APDB_MEAN", ix, iy, apdpn_apdpn[iAPDoAPDB][iMean]);
0588       setVal("APD_OVER_APDB_RMS", ix, iy, apdpn_apdpn[iAPDoAPDB][iRMS]);
0589       setVal("APD_OVER_APDB_M3", ix, iy, apdpn_apdpn[iAPDoAPDB][iM3]);  // fixme
0590       // JM
0591       setVal("SHAPE_COR", ix, iy, apdpn_ShapeCor);
0592       if (ab_tree) {
0593         setVal("ALPHA", ix, iy, ab_ab[iAlpha]);
0594         setVal("BETA", ix, iy, ab_ab[iBeta]);
0595       } else {
0596         setVal("ALPHA", ix, iy, 0.);
0597         setVal("BETA", ix, iy, 0.);
0598       }
0599       // NEW GHM 08/06
0600       setVal("TIME_MEAN", ix, iy, apdpn_apdpn[iTime][iMean]);
0601       setVal("TIME_RMS", ix, iy, apdpn_apdpn[iTime][iRMS]);
0602       setVal("TIME_M3", ix, iy, apdpn_apdpn[iTime][iM3]);
0603       setVal("TIME_NEVT", ix, iy, apdpn_apdpn[iTime][iNevt]);
0604     }
0605 
0606     //
0607     // PN primitives
0608     //
0609     t_name = lmfLaserName(ME::iLmfLaserPnPrim, _type, _color);
0610 
0611     nentries = pn_tree->GetEntriesFast();
0612     assert(nentries % 2 == 0);
0613     int module_(0);
0614     id1_ = _sm;
0615     id2_ = 0;
0616 
0617     Long64_t jentry = 0;
0618 
0619     while (jentry < nentries) {
0620       for (int jj = 0; jj < 2; jj++) {
0621         // jj=0 --> PNA
0622         // jj=1 --> PNB
0623 
0624         int zentry = jentry + jj;
0625         assert(zentry < nentries);
0626 
0627         ientry = pn_tree->LoadTree(zentry);
0628         assert(ientry >= 0);
0629         pn_tree->GetEntry(zentry);
0630 
0631         if (_side != pn_side)
0632           break;
0633 
0634         if (jj == 1)
0635           assert(pn_moduleID == module_);
0636         module_ = pn_moduleID;
0637         assert(pn_pnID == jj);
0638 
0639         // get the PN number
0640         std::pair<int, int> memPn_ = ME::pn(_lmr, module_, (ME::PN)jj);
0641         if (_isBarrel) {
0642           //          assert( memPn_.first%600==_dcc%600 );
0643           id1_ = _sm;
0644           id2_ = memPn_.second;
0645         } else {
0646           int dee_ = MEEEGeom::dee(_lmr);
0647           //          int mem_ = memPn_.first%600;
0648           //          if(      dee_==1 )
0649           //            {
0650           //              if( jj==ME::iPNA ) assert( mem_==50 );
0651           //              else               assert( mem_==51 );
0652           //            }
0653           //          else if( dee_==2 )
0654           //            {
0655           //              if( jj==ME::iPNA ) assert( mem_==47 );  // warning !
0656           //              else               assert( mem_==46 );
0657           //              //              assert( mem_==46 || mem_==47 );
0658           //            }
0659           //          else if( dee_==3 )
0660           //            {
0661           //              if( jj==ME::iPNA ) assert( mem_==1 );
0662           //              else               assert( mem_==2 );
0663           //            }
0664           //          else if( dee_==4 )
0665           //            {
0666           //              if( jj==ME::iPNA ) assert( mem_==5 );
0667           //              else               assert( mem_==6 );
0668           //            }
0669           id1_ = dee_;
0670           id2_ = (jj + 1) * 100 + memPn_.second;
0671         }
0672 
0673         if (_isBarrel) {
0674           channelView_ = iEB_LM_PN;
0675         } else {
0676           channelView_ = iEE_LM_PN;
0677         }
0678         logic_id_ = logicId(channelView_, id1_, id2_);
0679 
0680         i_t[t_name + separator + "LOGIC_ID"] = logic_id_;
0681         f_t[t_name + separator + "MEAN"] = pn_PN[iMean];
0682         f_t[t_name + separator + "RMS"] = pn_PN[iRMS];
0683         f_t[t_name + separator + "M3"] = pn_PN[iM3];
0684         f_t[t_name + separator + "PNA_OVER_PNB_MEAN"] = (jj == 0) ? pn_PNoPNB[iMean] : pn_PNoPNA[iMean];
0685         f_t[t_name + separator + "PNA_OVER_PNB_RMS"] = (jj == 0) ? pn_PNoPNB[iRMS] : pn_PNoPNA[iRMS];
0686         f_t[t_name + separator + "PNA_OVER_PNB_M3"] = (jj == 0) ? pn_PNoPNB[iM3] : pn_PNoPNA[iM3];
0687 
0688         t_t[t_name]->Fill();
0689       }
0690       //      std::cout << "Module=" << module_ << "\tPNA=" << pn_[0] << "\tPNB=" << pn_[1] << std::endl;
0691 
0692       //      if( _isBarrel )
0693       //        jentry += 4;
0694       //      else
0695       //        jentry += 2;
0696       jentry += 2;
0697     }
0698 
0699     logic_id_ = logicId(iECAL_LMR, _lmr);
0700 
0701     //
0702     // MATACQ primitives
0703     //
0704 
0705     if (mtq_tree) {
0706       t_name = lmfLaserName(ME::iLmfLaserPulse, _type, _color);
0707 
0708       nentries = mtq_tree->GetEntriesFast();
0709       assert(nentries == 2);
0710       for (Long64_t jentry = 0; jentry < nentries; jentry++) {
0711         ientry = mtq_tree->LoadTree(jentry);
0712         assert(ientry >= 0);
0713         mtq_tree->GetEntry(jentry);
0714 
0715         if (mtq_side != _side)
0716           continue;
0717 
0718         i_t[t_name + separator + "LOGIC_ID"] = logic_id_;
0719         i_t[t_name + separator + "FIT_METHOD"] = 0;  // fixme  --- what's this? ? ?
0720         f_t[t_name + separator + "MTQ_AMPL"] = mtq_mtq[iAmpl];
0721         f_t[t_name + separator + "MTQ_TIME"] = mtq_mtq[iPeak];
0722         f_t[t_name + separator + "MTQ_RISE"] = mtq_mtq[iTrise];
0723         f_t[t_name + separator + "MTQ_FWHM"] = mtq_mtq[iFwhm];
0724         f_t[t_name + separator + "MTQ_FW20"] = mtq_mtq[iFw20];
0725         f_t[t_name + separator + "MTQ_FW80"] = mtq_mtq[iFw80];
0726         f_t[t_name + separator + "MTQ_SLIDING"] =
0727             mtq_mtq[iSlide];  // fixme  --- sliding: max of average in sliding window
0728 
0729         t_t[t_name]->Fill();
0730       }
0731     } else {
0732       t_name = lmfLaserName(ME::iLmfLaserPulse, _type, _color);
0733 
0734       i_t[t_name + separator + "LOGIC_ID"] = logic_id_;
0735       i_t[t_name + separator + "FIT_METHOD"] = 0;  // fixme
0736       f_t[t_name + separator + "MTQ_AMPL"] = 0.0;
0737       f_t[t_name + separator + "MTQ_TIME"] = 0.0;
0738       f_t[t_name + separator + "MTQ_RISE"] = 0.0;
0739       f_t[t_name + separator + "MTQ_FWHM"] = 0.0;
0740       f_t[t_name + separator + "MTQ_FW20"] = 0.0;
0741       f_t[t_name + separator + "MTQ_FW80"] = 0.0;
0742       f_t[t_name + separator + "MTQ_SLIDING"] = 0.0;
0743 
0744       t_t[t_name]->Fill();
0745     }
0746 
0747     //
0748     // Laser Run
0749     //
0750     t_name = lmfLaserName(ME::iLmfLaserRun, _type);
0751     //std::cout << "Fill "<< t_name << std::endl;
0752     i_t[t_name + separator + "LOGIC_ID"] = logic_id_;
0753     i_t[t_name + separator + "NEVENTS"] = _events;
0754     i_t[t_name + separator + "QUALITY_FLAG"] = 1;  // fixme
0755     t_t[t_name]->Fill();
0756 
0757     //
0758     // Laser Config
0759     //
0760     t_name = lmfLaserName(ME::iLmfLaserConfig, _type);
0761     //std::cout << "Fill "<< t_name << std::endl;
0762     i_t[t_name + separator + "LOGIC_ID"] = logic_id_;
0763     i_t[t_name + separator + "WAVELENGTH"] = _color;
0764     i_t[t_name + separator + "VFE_GAIN"] = _mgpagain;      // fixme
0765     i_t[t_name + separator + "PN_GAIN"] = _memgain;        // fixme
0766     i_t[t_name + separator + "LSR_POWER"] = _power;        // will be available from MATACQ data
0767     i_t[t_name + separator + "LSR_ATTENUATOR"] = _filter;  // idem
0768     i_t[t_name + separator + "LSR_CURRENT"] = 0;           // idem
0769     i_t[t_name + separator + "LSR_DELAY_1"] = _delay;      // idem
0770     i_t[t_name + separator + "LSR_DELAY_2"] = 0;           // idem
0771     t_t[t_name]->Fill();
0772 
0773   } else if (_type == ME::iTestPulse) {
0774     //       nentries = tpapd_tree->GetEntriesFast();
0775     //       for( Long64_t jentry=0; jentry<nentries; jentry++ )
0776     //  {
0777     //    ientry = tpapd_tree->LoadTree( jentry );
0778     //    assert( ientry>=0 );
0779     //    nb     = tpapd_tree->GetEntry( jentry );
0780 
0781     //    if( tpapd_iphi<0 ) continue;
0782 
0783     //    const bool new_= true;
0784 
0785     //    int ix;
0786     //    int iy;
0787     //    int id2;
0788     //    if( new_ )
0789     //      {
0790     //        // for Cruzet3 , global coordinates
0791     //        if ( tpapd_side != _side ) continue;
0792     //        int ieta=tpapd_ieta;
0793     //        int iphi=tpapd_iphi;
0794     //        MEEBGeom::XYCoord xy_ = MEEBGeom::localCoord( ieta, iphi );
0795     //        ix = xy_.first;
0796     //        iy = xy_.second;
0797     //        id2 = ix*20 + iy;  // !!! TO BE CHECKED !!!
0798     //      }
0799     //    else
0800     //      {
0801     //        // for Cruzet2 , local coordinates
0802     //        ix = tpapd_ieta;
0803     //        iy = 19-tpapd_iphi;
0804     //        id2 = ix*20 + (20 - iy);  // !!! TO BE CHECKED !!!
0805     //      }
0806     //    //
0807 
0808     //    int id1 = _sm;   // fixme --- this is for barrel
0809     //    int logic_id_ = 1011000000;    // fixme
0810     //    logic_id_ += 10000*id1 + id2;
0811 
0812     //    int flag = tpapd_flag;
0813 
0814     //    setInt( "LOGIC_ID",           ix, iy,  logic_id_ );
0815     //    setInt( "FLAG",               ix, iy,  flag );
0816     //    setVal( "MEAN",               ix, iy,  tpapd_APD[iMean] );
0817     //    setVal( "RMS",                ix, iy,  tpapd_APD[iRMS] );
0818     //    setVal( "M3",                 ix, iy,  tpapd_APD[iM3] );
0819     //    setVal( "NEVT",               ix, iy,  tpapd_APD[iNevt] );
0820     //  }
0821 
0822     //       //
0823     //       // PN primitives
0824     //       //
0825     //       t_name = lmfLaserName( ME::iLmfTestPulsePnPrim, _type );
0826 
0827     //       nentries = tppn_tree->GetEntriesFast();
0828     //       assert( nentries%2==0 );
0829     //       int module_, pn_[2];
0830     //       int id1_(_sm), id2_(0);
0831     //       int logic_id_(0);
0832 
0833     //       Long64_t jentry=0;
0834     //       if( _side==1 ) jentry+=2;  // fixme : true also for endcaps?
0835     //       while( jentry<nentries )
0836     //  {
0837     //    for( int jj=0; jj<2; jj++ )
0838     //      {
0839     //        // jj=0 --> PNA
0840     //        // jj=1 --> PNB
0841 
0842     //        int zentry = jentry+jj;
0843     //        assert( zentry<nentries );
0844 
0845     //        ientry = tppn_tree->LoadTree( zentry );
0846     //        assert( ientry>=0 );
0847     //        nb     = tppn_tree->GetEntry( zentry );
0848 
0849     //        if( jj==1 ) assert( tppn_moduleID==module_ );
0850     //        module_ = tppn_moduleID;
0851     //        assert( tppn_pnID==jj );
0852 
0853     //        pn_[jj] = ( jj==0 ) ? _pn[module_].first : _pn[module_].second;
0854     //        id2_ = pn_[jj];
0855     //        logic_id_ = 1131000000 ;
0856     //        //      logic_id_ = 0;    // fixme
0857     //        logic_id_ += 10000*id1_ + id2_;
0858 
0859     //        i_t[t_name+separator+"LOGIC_ID"] = logic_id_;
0860     //        i_t[t_name+separator+"GAIN"]  = tppn_gain;
0861     //        f_t[t_name+separator+"MEAN"]  = tppn_PN[iMean];
0862     //        f_t[t_name+separator+"RMS"]   = tppn_PN[iRMS];
0863     //        f_t[t_name+separator+"M3"]    = tppn_PN[iM3];     // fixme --- peak?
0864 
0865     //        t_t[t_name]->Fill();
0866 
0867     //      }
0868 
0869     //    //      std::cout << "Module=" << module_ << "\tPNA=" << pn_[0] << "\tPNB=" << pn_[1] << std::endl;
0870 
0871     //    jentry += 4;
0872     //  }
0873 
0874     //       logic_id_ = 1041000000;
0875     //       logic_id_ += 10000*id1_;
0876     //       logic_id_ += id1_;
0877 
0878     //       //
0879     //       // Test Pulse Run
0880     //       //
0881     //       t_name = lmfLaserName( ME::iLmfTestPulseRun, _type );
0882     //       //std::cout << "Fill "<< t_name << std::endl;
0883     //       i_t[t_name+separator+"LOGIC_ID"]       = logic_id_;  // fixme --- is there a channelview for this?
0884     //       i_t[t_name+separator+"NEVENTS"]        =  _events;
0885     //       i_t[t_name+separator+"QUALITY_FLAG"]   = 1;                // fixme
0886     //       t_t[t_name]->Fill();
0887 
0888     //       //
0889     //       // Test Pulse Config
0890     //       //
0891     //       t_name = lmfLaserName( ME::iLmfTestPulseConfig, _type );
0892     //       //std::cout << "Fill "<< t_name << std::endl;
0893     //       i_t[t_name+separator+"LOGIC_ID"]        = logic_id_;  // fixme
0894     //       i_t[t_name+separator+"VFE_GAIN"]        = _mgpagain; // fixme
0895     //       i_t[t_name+separator+"PN_GAIN"]         = _memgain;  // fixme
0896     //       t_t[t_name]->Fill();
0897   }
0898 
0899   //
0900   // Laser Run IOV
0901   //
0902   t_name = "LMF_RUN_IOV";
0903   //std::cout << "Fill "<< t_name << std::endl;
0904   i_t[t_name + separator + "TAG_ID"] = 0;          // fixme
0905   i_t[t_name + separator + "SUB_RUN_NUM"] = _run;  // fixme
0906   i_t[t_name + separator + "SUB_RUN_START_LOW"] = ME::time_low(_ts_beg);
0907   i_t[t_name + separator + "SUB_RUN_START_HIGH"] = ME::time_high(_ts_beg);
0908   i_t[t_name + separator + "SUB_RUN_END_LOW"] = ME::time_low(_ts_end);
0909   i_t[t_name + separator + "SUB_RUN_END_HIGH"] = ME::time_high(_ts_end);
0910   i_t[t_name + separator + "DB_TIMESTAMP_LOW"] = ME::time_low(_ts);
0911   i_t[t_name + separator + "DB_TIMESTAMP_HIGH"] = ME::time_high(_ts);
0912   c_t[t_name + separator + "SUB_RUN_TYPE"] = "LASER TEST CRUZET";  //fixme
0913   t_t[t_name]->Fill();
0914 }
0915 
0916 void MELaserPrim::writeHistograms() {
0917   if (!init_ok)
0918     return;
0919 
0920   out_file = new TFile(_outfile, "RECREATE");
0921   //  out_file->cd();
0922 
0923   std::map<TString, TH2*>::iterator it;
0924 
0925   for (it = i_h.begin(); it != i_h.end(); ++it) {
0926     it->second->Write();
0927     delete it->second;
0928   }
0929 
0930   for (it = f_h.begin(); it != f_h.end(); ++it) {
0931     it->second->Write();
0932     delete it->second;
0933   }
0934 
0935   std::map<TString, TTree*>::iterator it_t;
0936   for (it_t = t_t.begin(); it_t != t_t.end(); ++it_t) {
0937     it_t->second->Write();
0938     delete it_t->second;
0939   }
0940 
0941   //  std::cout << "Closing " << _outfile << std::endl;
0942   out_file->Close();
0943   delete out_file;
0944   out_file = nullptr;
0945 }
0946 
0947 MELaserPrim::~MELaserPrim() {
0948   delete apdpn_tree;
0949   delete ab_tree;
0950   delete pn_tree;
0951   delete mtq_tree;
0952   delete tpapd_tree;
0953   delete tppn_tree;
0954   if (apdpn_file != nullptr) {
0955     //      std::cout << "Closing apdpn_file " << std::endl;
0956     apdpn_file->Close();
0957     delete apdpn_file;
0958     apdpn_file = nullptr;
0959   }
0960   if (ab_file != nullptr) {
0961     //      std::cout << "Closing ab_file " << std::endl;
0962     ab_file->Close();
0963     delete ab_file;
0964     ab_file = nullptr;
0965   }
0966   if (mtq_file != nullptr) {
0967     //      std::cout << "Closing mtq_file " << std::endl;
0968     mtq_file->Close();
0969     delete mtq_file;
0970     mtq_file = nullptr;
0971   }
0972   if (tpapd_file != nullptr) {
0973     //      std::cout << "Closing tpapd_file " << std::endl;
0974     tpapd_file->Close();
0975     delete tpapd_file;
0976     tpapd_file = nullptr;
0977   }
0978 }
0979 
0980 void MELaserPrim::print(std::ostream& o) {
0981   o << "DCC/SM/side/type/color/run/ts " << _dcc << "/" << _sm << "/" << _side << "/" << _type << "/" << _color << "/"
0982     << _run << "/" << _ts << std::endl;
0983 
0984   //   for( int ix=ixmin; ix<ixmax; ix++ )
0985   //     {
0986   //       for( int iy=iymin; iy<iymax; iy++ )
0987   //    {
0988   //      int flag     = getInt( _primStr+"FLAG", ix, iy );
0989   //      if( flag==0 ) continue;
0990   //      int logic_id = getInt( _primStr+"LOGIC_ID", ix, iy );
0991   //      float   apd = getVal( _primStr+"MEAN", ix, iy );
0992   //      o << "Crystal ieta=" << ix << "\tiphi=" << iy << "\tlogic_id=" << logic_id << "\tAPD=" << apd <<  endl;
0993   //    }
0994   //     }
0995 }
0996 
0997 TString MELaserPrim::lmfLaserName(int table, int type, int color) {
0998   TString str("LMF_ERROR");
0999   if (table < 0 || table >= ME::iSizeLmf)
1000     return str;
1001   if (color < 0 || color >= ME::iSizeC)
1002     return str;
1003 
1004   if (type == ME::iLaser) {
1005     TString colstr;
1006     switch (color) {
1007       case ME::iBlue:
1008         colstr = "_BLUE";
1009         break;
1010       case ME::iGreen:
1011         colstr = "_GREEN";
1012         break;
1013       case ME::iRed:
1014         colstr = "_RED";
1015         break;
1016       case ME::iIRed:
1017         colstr = "_IRED";
1018         break;
1019       default:
1020         abort();
1021     }
1022     str = "LMF_LASER";
1023     switch (table) {
1024       case ME::iLmfLaserRun:
1025         str = "LMF_RUN";
1026         break;
1027       case ME::iLmfLaserConfig:
1028         str += "_CONFIG";
1029         break;
1030       case ME::iLmfLaserPulse:
1031         str += colstr;
1032         str += "_PULSE";
1033         break;
1034       case ME::iLmfLaserPrim:
1035         str += colstr;
1036         str += "_PRIM";
1037         break;
1038       case ME::iLmfLaserPnPrim:
1039         str += colstr;
1040         str += "_PN_PRIM";
1041         break;
1042       default:
1043         abort();
1044     }
1045   } else if (type == ME::iTestPulse) {
1046     str = "LMF_TEST_PULSE";
1047     switch (table) {
1048       case ME::iLmfTestPulseRun:
1049         str = "LMF_RUN";
1050         break;
1051       case ME::iLmfTestPulseConfig:
1052         str += "_CONFIG";
1053         break;
1054       case ME::iLmfTestPulsePrim:
1055         str += "_PRIM";
1056         break;
1057       case ME::iLmfTestPulsePnPrim:
1058         str += "_PN_PRIM";
1059         break;
1060       default:
1061         abort();
1062     }
1063   }
1064   str += "_DAT";
1065   return str;
1066 }
1067 
1068 void MELaserPrim::addBranchI(const char* t_name_, const char* v_name_) {
1069   TString slashI("/i");  // Warning: always unsigned
1070   TString t_name(t_name_);
1071   TString v_name(v_name_);
1072   if (t_t.count(t_name) == 0)
1073     t_t[t_name] = new TTree(t_name, t_name);
1074   t_t[t_name]->Branch(v_name, &i_t[t_name + separator + v_name], v_name + slashI);
1075 }
1076 
1077 void MELaserPrim::addBranchF(const char* t_name_, const char* v_name_) {
1078   TString slashF("/F");
1079   TString t_name(t_name_);
1080   TString v_name(v_name_);
1081   if (t_t.count(t_name) == 0)
1082     t_t[t_name] = new TTree(t_name, t_name);
1083   t_t[t_name]->Branch(v_name, &f_t[t_name + separator + v_name], v_name + slashF);
1084 }
1085 
1086 void MELaserPrim::addBranchC(const char* t_name_, const char* v_name_) {
1087   TString slashC("/C");
1088   TString t_name(t_name_);
1089   TString v_name(v_name_);
1090   if (t_t.count(t_name) == 0)
1091     t_t[t_name] = new TTree(t_name, t_name);
1092   t_t[t_name]->Branch(v_name, &c_t[t_name + separator + v_name], v_name + slashC);
1093 }
1094 
1095 void MELaserPrim::bookHistoI(const char* h_name_, const char* v_name_) {
1096   TString i_name = TString(h_name_) + TString(v_name_);
1097   TH2* h_ = new TH2I(i_name, i_name, nx, ixmin, ixmax, ny, iymin, iymax);
1098   setHistoStyle(h_);
1099   i_h[i_name] = h_;
1100 }
1101 
1102 void MELaserPrim::bookHistoF(const char* h_name_, const char* v_name_) {
1103   TString d_name = TString(h_name_) + TString(v_name_);
1104   TH2* h_ = new TH2F(d_name, d_name, nx, ixmin, ixmax, ny, iymin, iymax);
1105   setHistoStyle(h_);
1106   f_h[d_name] = h_;
1107 }
1108 
1109 bool MELaserPrim::setInt(const char* name, int ix, int iy, int ival) {
1110   TString name_;
1111   if (_type == ME::iLaser)
1112     name_ = _primStr + name;
1113   else if (_type == ME::iTestPulse)
1114     name_ = _tpPrimStr + name;
1115 
1116   int _ival = getInt(name_, ix, iy);
1117   assert(_ival != -99);
1118   if (_ival != 0)
1119     return false;
1120 
1121   TH2I* h_ = (TH2I*)i_h[name_];
1122   assert(h_ != nullptr);
1123   h_->Fill(ix + 0.5, iy + 0.5, ival);
1124 
1125   return true;
1126 }
1127 
1128 bool MELaserPrim::setVal(const char* name, int ix, int iy, float val) {
1129   TString name_;
1130   if (_type == ME::iLaser)
1131     name_ = _primStr + name;
1132   else if (_type == ME::iTestPulse)
1133     name_ = _tpPrimStr + name;
1134 
1135   float _val = getVal(name_, ix, iy);
1136   assert(_val != -99);
1137   if (_val != 0)
1138     return false;
1139 
1140   TH2F* h_ = (TH2F*)f_h[name_];
1141   assert(h_ != nullptr);
1142 
1143   h_->Fill(ix + 0.5, iy + 0.5, val);
1144 
1145   return true;
1146 }
1147 
1148 Int_t MELaserPrim::getInt(const char* name, int ix, int iy) {
1149   Int_t ival = -99;
1150   if (i_h.count(name) == 1) {
1151     TH2I* h_ = (TH2I*)i_h[name];
1152     assert(h_ != nullptr);
1153     int binx = h_->GetXaxis()->FindBin(ix + 0.5);
1154     int biny = h_->GetYaxis()->FindBin(iy + 0.5);
1155     ival = (Int_t)h_->GetCellContent(binx, biny);
1156   }
1157   return ival;
1158 }
1159 
1160 Float_t MELaserPrim::getVal(const char* name, int ix, int iy) {
1161   Float_t val = -99.;
1162   if (f_h.count(name) == 1) {
1163     TH2F* h_ = (TH2F*)f_h[name];
1164     assert(h_ != nullptr);
1165     int binx = h_->GetXaxis()->FindBin(ix + 0.5);
1166     int biny = h_->GetYaxis()->FindBin(iy + 0.5);
1167     val = h_->GetCellContent(binx, biny);
1168   }
1169   return val;
1170 }
1171 
1172 bool MELaserPrim::setInt(const char* tname, const char* vname, int ival) {
1173   TString key_(tname);
1174   key_ += separator;
1175   key_ += vname;
1176   assert(i_t.count(key_) == 1);
1177   i_t[key_] = ival;
1178   return true;
1179 }
1180 
1181 bool MELaserPrim::setVal(const char* tname, const char* vname, float val) {
1182   TString key_(tname);
1183   key_ += separator;
1184   key_ += vname;
1185 
1186   // ghm
1187   if (f_t.count(key_) != 1) {
1188     std::cout << key_ << std::endl;
1189   }
1190   assert(f_t.count(key_) == 1);
1191   f_t[key_] = val;
1192   return true;
1193 }
1194 
1195 bool MELaserPrim::fill(const char* tname) {
1196   TString key_(tname);
1197   assert(t_t.count(key_) == 1);
1198   t_t[key_]->Fill();
1199   return true;
1200 }
1201 
1202 void MELaserPrim::setHistoStyle(TH1* h) {
1203   if (h == nullptr)
1204     return;
1205 
1206   float _scale = 1;
1207 
1208   h->SetLineColor(4);
1209   h->SetLineWidth(1);
1210   h->SetFillColor(38);
1211   TAxis* axis[3];
1212   axis[0] = h->GetXaxis();
1213   axis[1] = h->GetYaxis();
1214   axis[2] = h->GetZaxis();
1215   for (int ii = 0; ii < 3; ii++) {
1216     TAxis* a = axis[ii];
1217     if (!a)
1218       continue;
1219     a->SetLabelFont(132);
1220     a->SetLabelOffset(_scale * 0.005);
1221     a->SetLabelSize(_scale * 0.04);
1222     a->SetTitleFont(132);
1223     a->SetTitleOffset(_scale * 1);
1224     a->SetTitleSize(_scale * 0.04);
1225   }
1226   h->SetStats(kTRUE);
1227 }
1228 
1229 void MELaserPrim::refresh() {
1230   std::map<TString, TH2*>::iterator it;
1231 
1232   for (it = i_h.begin(); it != i_h.end(); ++it) {
1233     delete it->second;
1234     it->second = nullptr;
1235   }
1236   i_h.clear();
1237 
1238   for (it = f_h.begin(); it != f_h.end(); ++it) {
1239     delete it->second;
1240     it->second = nullptr;
1241   }
1242   f_h.clear();
1243 
1244   std::map<TString, TTree*>::iterator it_t;
1245   for (it_t = t_t.begin(); it_t != t_t.end(); ++it_t) {
1246     delete it_t->second;
1247     it->second = nullptr;
1248   }
1249   t_t.clear();
1250 }