Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:18:51

0001 #include "TriggerValid.h"
0002 
0003 // Trigger includes
0004 #include <calobase/TowerInfo.h>
0005 #include <calobase/TowerInfoContainer.h>
0006 #include <calobase/TowerInfoDefs.h>
0007 
0008 #include <calotrigger/LL1Out.h>
0009 #include <calotrigger/LL1Outv1.h>
0010 #include <calotrigger/TriggerDefs.h>
0011 #include <calotrigger/TriggerPrimitiveContainerv1.h>
0012 #include <calotrigger/TriggerPrimitivev1.h>
0013 
0014 #include <qautils/QAHistManagerDef.h>
0015 
0016 #include <fun4all/Fun4AllHistoManager.h>
0017 #include <fun4all/Fun4AllReturnCodes.h>
0018 
0019 #include <phool/getClass.h>
0020 #include <phool/phool.h>  // for PHWHERE
0021 
0022 #include <TFile.h>
0023 #include <TH1.h>
0024 #include <TH2.h>
0025 #include <TLorentzVector.h>
0026 #include <TNtuple.h>
0027 #include <TProfile2D.h>
0028 #include <TSystem.h>
0029 #include <TTree.h>
0030 
0031 #include <boost/format.hpp>
0032 
0033 #include <cmath>     // for log10, pow, sqrt, abs, M_PI
0034 #include <iostream>  // for operator<<, endl, basic_...
0035 #include <limits>
0036 #include <map>  // for operator!=, _Rb_tree_con...
0037 #include <string>
0038 #include <utility>  // for pair
0039 
0040 TriggerValid::TriggerValid(const std::string& name)
0041   : SubsysReco(name)
0042 {
0043 }
0044 
0045 int TriggerValid::Init(PHCompositeNode* /*unused*/)
0046 {
0047   auto hm = QAHistManagerDef::getHistoManager();
0048   assert(hm);
0049   // create and register your histos (all types) here
0050 
0051   if (m_debug)
0052   {
0053     std::cout << "In TriggerValid::Init" << std::endl;
0054   }
0055   auto h_gl1_triggers = new TH1D("h_gl1_triggers", ";Trigger Name;#Triggers", 64, -0.5, 63.5);
0056   hm->registerHisto(h_gl1_triggers);
0057   for (int i = 0; i < 8; i++)
0058   {
0059     auto h_gl1_photon_energy = new TH1D((boost::format("h_gl1_photon_energy_%d") % i).str().c_str(), ";8x8 EMCAL energy [GeV]; counts", 100, 0, 50);
0060     auto h_gl1_jetpatch_energy = new TH1D((boost::format("h_gl1_jetpatch_energy_%d") % i).str().c_str(), ";.8 x .8 EMCAL+HCAL energy [GeV]", 50, 0, 50);
0061     hm->registerHisto(h_gl1_photon_energy);
0062     hm->registerHisto(h_gl1_jetpatch_energy);
0063   }
0064 
0065   auto h_emu_emcal_2x2_frequency = new TH2F("h_emu_emcal_2x2_frequency", ";#eta;#phi", 48, 0, 96, 128, 0, 256);
0066   hm->registerHisto(h_emu_emcal_2x2_frequency);
0067   auto h_emu_emcal_8x8_frequency = new TH2F("h_emu_emcal_8x8_frequency", ";#eta;#phi", 12, 0, 96, 32, 0, 256);
0068   hm->registerHisto(h_emu_emcal_8x8_frequency);
0069   auto h_emu_ihcal_2x2_frequency = new TH2F("h_emu_ihcal_2x2_frequency", ";#eta;#phi", 12, 0, 24, 32, 0, 64);
0070   hm->registerHisto(h_emu_ihcal_2x2_frequency);
0071   auto h_emu_ohcal_2x2_frequency = new TH2F("h_emu_ohcal_2x2_frequency", ";#eta;#phi", 12, 0, 24, 32, 0, 64);
0072   hm->registerHisto(h_emu_ohcal_2x2_frequency);
0073   auto h_emu_hcal_2x2_frequency = new TH2F("h_emu_hcal_2x2_frequency", ";#eta;#phi", 12, 0, 24, 32, 0, 64);
0074   hm->registerHisto(h_emu_hcal_2x2_frequency);
0075   for (int i = 0; i < 4; i++)
0076   {
0077     auto h_emu_jet_frequency_trig = new TH2F((boost::format("h_emu_jet_frequency_%d") % i).str().c_str(), ";#eta;#phi", 9, 0, 9, 32, 0, 32);
0078     hm->registerHisto(h_emu_jet_frequency_trig);
0079     auto h_emu_photon_frequency_trig = new TH2F((boost::format("h_emu_photon_frequency_%d") % i).str().c_str(), ";#eta;#phi", 12, 0, 12, 32, 0, 32);
0080     hm->registerHisto(h_emu_photon_frequency_trig);
0081   }
0082 
0083   auto h_emu_emcal_2x2_avg_out = new TProfile2D("h_emu_emcal_2x2_avg_out", ";#eta;#phi", 48, 0, 96, 128, 0, 256);
0084   hm->registerHisto(h_emu_emcal_2x2_avg_out);
0085   auto h_emu_emcal_8x8_avg_out = new TProfile2D("h_emu_emcal_8x8_avg_out", ";#eta;#phi", 12, 0, 96, 32, 0, 256);
0086   hm->registerHisto(h_emu_emcal_8x8_avg_out);
0087   auto h_emu_ihcal_2x2_avg_out = new TProfile2D("h_emu_ihcal_2x2_avg_out", ";#eta;#phi", 12, 0, 24, 32, 0, 64);
0088   hm->registerHisto(h_emu_ihcal_2x2_avg_out);
0089   auto h_emu_ohcal_2x2_avg_out = new TProfile2D("h_emu_ohcal_2x2_avg_out", ";#eta;#phi", 12, 0, 24, 32, 0, 64);
0090   hm->registerHisto(h_emu_ohcal_2x2_avg_out);
0091   auto h_emu_hcal_2x2_avg_out = new TProfile2D("h_emu_hcal_2x2_avg_out", ";#eta;#phi", 12, 0, 24, 64, 0, 64);
0092   hm->registerHisto(h_emu_hcal_2x2_avg_out);
0093   auto h_emcal_2x2_energy_lutsum = new TH2F("h_emcal_2x2_energy_lutsum", ";LUT output; Energy [GeV]", 256, -0.5, 255.5, 200, 0, 20);
0094   hm->registerHisto(h_emcal_2x2_energy_lutsum);
0095   auto h_emcal_8x8_energy_lutsum = new TH2F("h_emcal_8x8_energy_lutsum", ";LUT output; Energy [GeV]", 256, -0.5, 255.5, 200, 0, 20);
0096   hm->registerHisto(h_emcal_8x8_energy_lutsum);
0097   auto h_hcal_2x2_energy_lutsum = new TH2F("h_hcal_2x2_energy_lutsum", ";LUT output; Energy [GeV]", 256, -0.5, 255.5, 200, 0, 20);
0098   hm->registerHisto(h_hcal_2x2_energy_lutsum);
0099   auto h_hcalin_2x2_energy_lutsum = new TH2F("h_hcalin_2x2_energy_lutsum", ";LUT output; Energy [GeV]", 256, -0.5, 255.5, 200, 0, 20);
0100   hm->registerHisto(h_hcalin_2x2_energy_lutsum);
0101   auto h_hcalout_2x2_energy_lutsum = new TH2F("h_hcalout_2x2_energy_lutsum", ";LUT output; Energy [GeV]", 256, -0.5, 255.5, 200, 0, 20);
0102   hm->registerHisto(h_hcalout_2x2_energy_lutsum);
0103 
0104   auto h_jet_energy_lutsum = new TH2F("h_jet_energy_lutsum", ";LUT output; Energy [GeV]", 4096, -0.5, 4095.5, 200, 0, 20);
0105   hm->registerHisto(h_jet_energy_lutsum);
0106 
0107   if (m_debug)
0108   {
0109     std::cout << "Leaving TriggerValid::Init" << std::endl;
0110   }
0111   return 0;
0112 }
0113 
0114 int TriggerValid::process_event(PHCompositeNode* topNode)
0115 {
0116   _eventcounter++;
0117 
0118   process_towers(topNode);
0119 
0120   return Fun4AllReturnCodes::EVENT_OK;
0121 }
0122 
0123 int TriggerValid::process_towers(PHCompositeNode* topNode)
0124 {
0125   if (m_debug)
0126   {
0127     std::cout << _eventcounter << std::endl;
0128   }
0129 
0130   TowerInfoContainer* towers_emcal = findNode::getClass<TowerInfoContainer>(topNode, "TOWERINFO_CALIB_CEMC");
0131 
0132   TowerInfoContainer* towers_hcalin = findNode::getClass<TowerInfoContainer>(topNode, "TOWERINFO_CALIB_HCALIN");
0133 
0134   TowerInfoContainer* towers_hcalout = findNode::getClass<TowerInfoContainer>(topNode, "TOWERINFO_CALIB_HCALOUT");
0135 
0136   LL1Out* ll1out_jet = findNode::getClass<LL1Out>(topNode, "LL1OUT_JET");
0137 
0138   LL1Out* ll1out_photon = findNode::getClass<LL1Out>(topNode, "LL1OUT_PHOTON");
0139 
0140   // y  LL1Out* ll1out_pair = findNode::getClass<LL1Out>(topNode, "LL1OUT_PAIR");
0141 
0142   TriggerPrimitiveContainer* trigger_primitives_emcal = findNode::getClass<TriggerPrimitiveContainer>(topNode, "TRIGGERPRIMITIVES_EMCAL");
0143 
0144   TriggerPrimitiveContainer* trigger_primitives_emcal_ll1 = findNode::getClass<TriggerPrimitiveContainer>(topNode, "TRIGGERPRIMITIVES_EMCAL_LL1");
0145 
0146   TriggerPrimitiveContainer* trigger_primitives_hcalin = findNode::getClass<TriggerPrimitiveContainer>(topNode, "TRIGGERPRIMITIVES_HCALIN");
0147 
0148   TriggerPrimitiveContainer* trigger_primitives_hcalout = findNode::getClass<TriggerPrimitiveContainer>(topNode, "TRIGGERPRIMITIVES_HCALOUT");
0149 
0150   TriggerPrimitiveContainer* trigger_primitives_hcal_ll1 = findNode::getClass<TriggerPrimitiveContainer>(topNode, "TRIGGERPRIMITIVES_HCAL_LL1");
0151 
0152   Gl1Packet* gl1_packet = findNode::getClass<Gl1Packet>(topNode, "GL1Packet");
0153 
0154   auto hm = QAHistManagerDef::getHistoManager();
0155   assert(hm);
0156 
0157   auto h_emu_emcal_2x2_frequency = dynamic_cast<TH2*>(hm->getHisto("h_emu_emcal_2x2_frequency"));
0158   auto h_emu_emcal_8x8_frequency = dynamic_cast<TH2*>(hm->getHisto("h_emu_emcal_8x8_frequency"));
0159   auto h_emu_ihcal_2x2_frequency = dynamic_cast<TH2*>(hm->getHisto("h_emu_ihcal_2x2_frequency"));
0160   auto h_emu_ohcal_2x2_frequency = dynamic_cast<TH2*>(hm->getHisto("h_emu_ohcal_2x2_frequency"));
0161   auto h_emu_hcal_2x2_frequency = dynamic_cast<TH2*>(hm->getHisto("h_emu_hcal_2x2_frequency"));
0162 
0163   TH2* h_emu_jet_frequency_trig[4];
0164   TH2* h_emu_photon_frequency_trig[4];
0165   for (int i = 0; i < 4; i++)
0166   {
0167     h_emu_jet_frequency_trig[i] = dynamic_cast<TH2*>(hm->getHisto((boost::format("h_emu_jet_frequency_%d") % i).str().c_str()));
0168     h_emu_photon_frequency_trig[i] = dynamic_cast<TH2*>(hm->getHisto((boost::format("h_emu_photon_frequency_%d") % i).str().c_str()));
0169   }
0170 
0171   auto h_gl1_triggers = dynamic_cast<TH1*>(hm->getHisto("h_gl1_triggers"));
0172   TH1* h_gl1_jetpatch_energy[8];
0173   TH1* h_gl1_photon_energy[8];
0174 
0175   for (int i = 0; i < 8; i++)
0176   {
0177     h_gl1_jetpatch_energy[i] = dynamic_cast<TH1*>(hm->getHisto((boost::format("h_gl1_jetpatch_energy_%d") % i).str().c_str()));
0178     h_gl1_photon_energy[i] = dynamic_cast<TH1*>(hm->getHisto((boost::format("h_gl1_photon_energy_%d") % i).str().c_str()));
0179   }
0180 
0181   auto h_emu_emcal_2x2_avg_out = dynamic_cast<TProfile2D*>(hm->getHisto("h_emu_emcal_2x2_avg_out"));
0182   auto h_emu_emcal_8x8_avg_out = dynamic_cast<TProfile2D*>(hm->getHisto("h_emu_emcal_8x8_avg_out"));
0183   auto h_emu_ihcal_2x2_avg_out = dynamic_cast<TProfile2D*>(hm->getHisto("h_emu_ihcal_2x2_avg_out"));
0184   auto h_emu_ohcal_2x2_avg_out = dynamic_cast<TProfile2D*>(hm->getHisto("h_emu_ohcal_2x2_avg_out"));
0185   auto h_emu_hcal_2x2_avg_out = dynamic_cast<TProfile2D*>(hm->getHisto("h_emu_hcal_2x2_avg_out"));
0186   auto h_emcal_2x2_energy_lutsum = dynamic_cast<TH2*>(hm->getHisto("h_emcal_2x2_energy_lutsum"));
0187   auto h_emcal_8x8_energy_lutsum = dynamic_cast<TH2*>(hm->getHisto("h_emcal_8x8_energy_lutsum"));
0188   auto h_hcal_2x2_energy_lutsum = dynamic_cast<TH2*>(hm->getHisto("h_hcal_2x2_energy_lutsum"));
0189   auto h_hcalin_2x2_energy_lutsum = dynamic_cast<TH2*>(hm->getHisto("h_hcalin_2x2_energy_lutsum"));
0190   auto h_hcalout_2x2_energy_lutsum = dynamic_cast<TH2*>(hm->getHisto("h_hcalout_2x2_energy_lutsum"));
0191   auto h_jet_energy_lutsum = dynamic_cast<TH2*>(hm->getHisto("h_jet_energy_lutsum"));
0192 
0193   std::map<TriggerDefs::TriggerSumKey, unsigned int> v_emcal_emu_2x2 = {};
0194 
0195   std::map<TriggerDefs::TriggerSumKey, unsigned int> v_emcal_emu_8x8 = {};
0196 
0197   std::map<TriggerDefs::TriggerSumKey, unsigned int> v_hcal_emu_2x2 = {};
0198   std::map<TriggerDefs::TriggerSumKey, unsigned int> v_hcalin_emu_2x2 = {};
0199   std::map<TriggerDefs::TriggerSumKey, unsigned int> v_hcalout_emu_2x2 = {};
0200 
0201   std::map<TriggerDefs::TriggerSumKey, unsigned int> v_jet_emu = {};
0202 
0203   std::vector<int> trig_bits{};
0204   if (gl1_packet)
0205   {
0206     ULong64_t gl1_scaledvec = gl1_packet->lValue(0, "ScaledVector");
0207     for (unsigned int bit = 0; bit < 64; bit++)
0208     {
0209       if (((gl1_scaledvec >> bit) & 0x1U) == 0x1U)
0210       {
0211         trig_bits.push_back(bit);
0212         h_gl1_triggers->Fill(bit);
0213       }
0214     }
0215   }
0216 
0217   if (trigger_primitives_emcal)
0218   {
0219     TriggerPrimitiveContainerv1::Range range = trigger_primitives_emcal->getTriggerPrimitives();
0220     for (TriggerPrimitiveContainerv1::Iter iter = range.first; iter != range.second; ++iter)
0221     {
0222       TriggerPrimitive* trigger_primitive = (*iter).second;
0223       TriggerPrimitivev1::Range srange = trigger_primitive->getSums();
0224       for (TriggerPrimitive::Iter siter = srange.first; siter != srange.second; ++siter)
0225       {
0226         std::vector<unsigned int>* sum = (*siter).second;
0227         std::vector<unsigned int>::iterator it = max_element(sum->begin(), sum->end());
0228         unsigned int summ = 0;
0229         unsigned int sumk = (*siter).first;
0230 
0231         if (it != sum->end())
0232         {
0233           summ = (*it);
0234         }
0235 
0236         uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 4 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0237         uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk) + 4 * TriggerDefs::getPrimitiveEtaId_from_TriggerSumKey(sumk);
0238 
0239         if (summ)
0240         {
0241           float i2x2x = h_emu_emcal_2x2_frequency->GetXaxis()->GetBinCenter(sum_eta + 1);
0242           float i2x2y = h_emu_emcal_2x2_frequency->GetYaxis()->GetBinCenter(sum_phi + 1);
0243           h_emu_emcal_2x2_frequency->Fill(i2x2x, i2x2y, 1);
0244           h_emu_emcal_2x2_avg_out->Fill(i2x2x, i2x2y, summ);
0245         }
0246         v_emcal_emu_2x2[sumk] = summ;
0247       }
0248     }
0249   }
0250 
0251   if (trigger_primitives_emcal_ll1)
0252   {
0253     TriggerPrimitiveContainerv1::Range range = trigger_primitives_emcal_ll1->getTriggerPrimitives();
0254     for (TriggerPrimitiveContainerv1::Iter iter = range.first; iter != range.second; ++iter)
0255     {
0256       TriggerPrimitive* trigger_primitive = (*iter).second;
0257       TriggerPrimitivev1::Range srange = trigger_primitive->getSums();
0258       for (TriggerPrimitive::Iter siter = srange.first; siter != srange.second; ++siter)
0259       {
0260         std::vector<unsigned int>* sum = (*siter).second;
0261         std::vector<unsigned int>::iterator it = max_element(sum->begin(), sum->end());
0262         unsigned int summ = 0;
0263         unsigned int sumk = (*siter).first;
0264 
0265         if (it != sum->end())
0266         {
0267           summ = (*it);
0268         }
0269 
0270         uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 2 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0271         uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk);
0272 
0273         if (summ)
0274         {
0275           float i2x2x = h_emu_emcal_8x8_frequency->GetXaxis()->GetBinCenter(sum_eta + 1);
0276           float i2x2y = h_emu_emcal_8x8_frequency->GetYaxis()->GetBinCenter(sum_phi + 1);
0277 
0278           h_emu_emcal_8x8_frequency->Fill(i2x2x, i2x2y, 1);
0279           h_emu_emcal_8x8_avg_out->Fill(i2x2x, i2x2y, summ);
0280         }
0281         v_emcal_emu_8x8[sumk] = summ;
0282       }
0283     }
0284   }
0285 
0286   if (trigger_primitives_hcalin)
0287   {
0288     TriggerPrimitiveContainerv1::Range range = trigger_primitives_hcalin->getTriggerPrimitives();
0289     for (TriggerPrimitiveContainerv1::Iter iter = range.first; iter != range.second; ++iter)
0290     {
0291       TriggerPrimitive* trigger_primitive = (*iter).second;
0292       TriggerPrimitivev1::Range srange = trigger_primitive->getSums();
0293       for (TriggerPrimitive::Iter siter = srange.first; siter != srange.second; ++siter)
0294       {
0295         std::vector<unsigned int>* sum = (*siter).second;
0296         std::vector<unsigned int>::iterator it = max_element(sum->begin(), sum->end());
0297         unsigned int summ = 0;
0298         unsigned int sumk = (*siter).first;
0299 
0300         if (it != sum->end())
0301         {
0302           summ = (*it);
0303         }
0304 
0305         uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 4 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0306         uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk) + 4 * TriggerDefs::getPrimitiveEtaId_from_TriggerSumKey(sumk);
0307 
0308         if (summ)
0309         {
0310           float i2x2x = h_emu_ihcal_2x2_frequency->GetXaxis()->GetBinCenter(sum_eta + 1);
0311           float i2x2y = h_emu_ihcal_2x2_frequency->GetYaxis()->GetBinCenter(sum_phi + 1);
0312           h_emu_ihcal_2x2_frequency->Fill(i2x2x, i2x2y, 1);
0313           h_emu_ihcal_2x2_avg_out->Fill(i2x2x, i2x2y, summ);
0314         }
0315         v_hcalin_emu_2x2[sumk] = summ;
0316       }
0317     }
0318   }
0319 
0320   if (trigger_primitives_hcalout)
0321   {
0322     TriggerPrimitiveContainerv1::Range range = trigger_primitives_hcalout->getTriggerPrimitives();
0323     for (TriggerPrimitiveContainerv1::Iter iter = range.first; iter != range.second; ++iter)
0324     {
0325       TriggerPrimitive* trigger_primitive = (*iter).second;
0326       TriggerPrimitivev1::Range srange = trigger_primitive->getSums();
0327       for (TriggerPrimitive::Iter siter = srange.first; siter != srange.second; ++siter)
0328       {
0329         std::vector<unsigned int>* sum = (*siter).second;
0330         std::vector<unsigned int>::iterator it = max_element(sum->begin(), sum->end());
0331         unsigned int summ = 0;
0332         unsigned int sumk = (*siter).first;
0333 
0334         if (it != sum->end())
0335         {
0336           summ = (*it);
0337         }
0338         uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 4 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0339         uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk) + 4 * TriggerDefs::getPrimitiveEtaId_from_TriggerSumKey(sumk);
0340 
0341         if (summ)
0342         {
0343           float i2x2x = h_emu_ohcal_2x2_frequency->GetXaxis()->GetBinCenter(sum_eta + 1);
0344           float i2x2y = h_emu_ohcal_2x2_frequency->GetYaxis()->GetBinCenter(sum_phi + 1);
0345           h_emu_ohcal_2x2_frequency->Fill(i2x2x, i2x2y, 1);
0346           h_emu_ohcal_2x2_avg_out->Fill(i2x2x, i2x2y, summ);
0347         }
0348         v_hcalout_emu_2x2[sumk] = summ;
0349       }
0350     }
0351   }
0352 
0353   if (trigger_primitives_hcal_ll1)
0354   {
0355     TriggerPrimitiveContainerv1::Range range = trigger_primitives_hcal_ll1->getTriggerPrimitives();
0356     for (TriggerPrimitiveContainerv1::Iter iter = range.first; iter != range.second; ++iter)
0357     {
0358       TriggerPrimitive* trigger_primitive = (*iter).second;
0359       TriggerPrimitivev1::Range srange = trigger_primitive->getSums();
0360       for (TriggerPrimitive::Iter siter = srange.first; siter != srange.second; ++siter)
0361       {
0362         std::vector<unsigned int>* sum = (*siter).second;
0363         std::vector<unsigned int>::iterator it = max_element(sum->begin(), sum->end());
0364         unsigned int summ = 0;
0365         unsigned int sumk = (*siter).first;
0366 
0367         if (it != sum->end())
0368         {
0369           summ = (*it);
0370         }
0371 
0372         uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 2 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0373         uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk);
0374 
0375         if (summ)
0376         {
0377           float i2x2x = h_emu_hcal_2x2_frequency->GetXaxis()->GetBinCenter(sum_eta + 1);
0378           float i2x2y = h_emu_hcal_2x2_frequency->GetYaxis()->GetBinCenter(sum_phi + 1);
0379           h_emu_hcal_2x2_frequency->Fill(i2x2x, i2x2y, 1);
0380           h_emu_hcal_2x2_avg_out->Fill(i2x2x, i2x2y, summ);
0381         }
0382         v_hcal_emu_2x2[sumk] = summ;
0383       }
0384     }
0385   }
0386 
0387   if (ll1out_photon)
0388   {
0389     for (int i = 0; i < 4; i++)
0390     {
0391       auto triggered_sums = ll1out_photon->getTriggeredSumKeys(i + 1);
0392       for (auto& key : triggered_sums)
0393       {
0394         unsigned int phi = TriggerDefs::getSumPhiId(key) + TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(key) * 2;
0395         unsigned int eta = TriggerDefs::getSumEtaId(key);
0396 
0397         h_emu_photon_frequency_trig[i]->Fill(eta, phi);
0398       }
0399     }
0400   }
0401 
0402   if (ll1out_jet)
0403   {
0404     for (int i = 0; i < 4; i++)
0405     {
0406       auto triggered_sums = ll1out_jet->getTriggeredSumKeys(i + 1);
0407       for (auto& key : triggered_sums)
0408       {
0409         int eta = (key >> 16U) & 0xffffU;
0410         int phi = (key & 0xffffU);
0411         h_emu_jet_frequency_trig[i]->Fill(eta, phi);
0412       }
0413     }
0414   }
0415 
0416   // h_emcal_2x2_energy_lutsum = new TH2F("h_emcal_2x2_energy_lutsum",";Energy [GeV];LUT output", 100, 0, 10, 64, 0, 256);
0417   // h_emcal_8x8_energy_lutsum = new TH2F("h_emcal_8x8_energy_lutsum",";Energy [GeV];LUT output", 100, 0, 10, 64, 0, 256);
0418   // h_hcal_2x2_energy_lutsum = new TH2F("h_hcal_2x2_energy_lutsum",";Energy [GeV];LUT output", 100, 0, 10, 64, 0, 256);
0419   // h_jet_energy_lutsum = new TH2F("h_jet_energy_lutsum",";Energy [GeV];LUT output", 100, 0, 10, 64, 0, 256*16);
0420 
0421   float emcal_energies[12][35]{};
0422   float hcal_energies[12][35]{};
0423   float max_energy_emcal = 0.0;
0424   if (towers_emcal)
0425   {
0426     // go through the emulated 2x2 map for emcal
0427     for (auto& it : v_emcal_emu_2x2)
0428     {
0429       unsigned int sumk = it.first;
0430       uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 4 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0431       uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk) + 4 * TriggerDefs::getPrimitiveEtaId_from_TriggerSumKey(sumk);
0432 
0433       float energy_sum = 0.0;
0434       for (int itower = 0; itower < 4; itower++)
0435       {
0436         int ieta = sum_eta * 2 + itower % 2;
0437         int iphi = sum_phi * 2 + itower / 2;
0438         TowerInfo* tower = towers_emcal->get_tower_at_key(TowerInfoDefs::encode_emcal(ieta, iphi));
0439         float offlineenergy = tower->get_energy();
0440         if (!tower->get_isGood())
0441         {
0442           continue;
0443         }
0444         energy_sum += offlineenergy;
0445       }
0446       h_emcal_2x2_energy_lutsum->Fill(it.second, energy_sum);
0447     }
0448 
0449     // now the 8x8
0450 
0451     for (auto& it : v_emcal_emu_8x8)
0452     {
0453       unsigned int sumk = it.first;
0454       uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 2 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0455       uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk);
0456 
0457       float energy_sum = 0.0;
0458       for (int itower = 0; itower < 64; itower++)
0459       {
0460         int ieta = sum_eta * 8 + itower % 8;
0461         int iphi = sum_phi * 8 + itower / 8;
0462         TowerInfo* tower = towers_emcal->get_tower_at_key(TowerInfoDefs::encode_emcal(ieta, iphi));
0463         float offlineenergy = tower->get_energy();
0464         if (!tower->get_isGood())
0465         {
0466           continue;
0467         }
0468         energy_sum += offlineenergy;
0469       }
0470       if (energy_sum > max_energy_emcal)
0471       {
0472         max_energy_emcal = energy_sum;
0473       }
0474       emcal_energies[sum_eta][sum_phi] = energy_sum;
0475       h_emcal_8x8_energy_lutsum->Fill(it.second, energy_sum);
0476     }
0477   }
0478 
0479   if (towers_hcalin || towers_hcalout)
0480   {
0481     // go through the emulated 2x2 map for emcal
0482     for (auto& it : v_hcal_emu_2x2)
0483     {
0484       unsigned int sumk = it.first;
0485       uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 2 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0486       uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk);
0487 
0488       float energy_sum = 0.0;
0489       for (int itower = 0; itower < 4; itower++)
0490       {
0491         int ieta = sum_eta * 2 + itower % 2;
0492         int iphi = sum_phi * 2 + itower / 2;
0493         if (towers_hcalin)
0494         {
0495           TowerInfo* tower = towers_hcalin->get_tower_at_key(TowerInfoDefs::encode_hcal(ieta, iphi));
0496           float offlineenergy = tower->get_energy();
0497 
0498           if (!tower->get_isGood())
0499           {
0500             continue;
0501           }
0502           energy_sum += offlineenergy;
0503         }
0504         if (towers_hcalin)
0505         {
0506           TowerInfo* tower = towers_hcalout->get_tower_at_key(TowerInfoDefs::encode_hcal(ieta, iphi));
0507           float offlineenergy = tower->get_energy();
0508           if (!tower->get_isGood())
0509           {
0510             continue;
0511           }
0512           energy_sum += offlineenergy;
0513         }
0514       }
0515       h_hcal_2x2_energy_lutsum->Fill(it.second, energy_sum);
0516       hcal_energies[sum_eta][sum_phi] = energy_sum;
0517     }
0518 
0519     // go through the emulated 2x2 map for emcal
0520     for (auto& it : v_hcalin_emu_2x2)
0521     {
0522       unsigned int sumk = it.first;
0523       uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 2 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0524       uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk);
0525 
0526       float energy_sum = 0.0;
0527       for (int itower = 0; itower < 4; itower++)
0528       {
0529         int ieta = sum_eta * 2 + itower % 2;
0530         int iphi = sum_phi * 2 + itower / 2;
0531         if (towers_hcalin)
0532         {
0533           TowerInfo* tower = towers_hcalin->get_tower_at_key(TowerInfoDefs::encode_hcal(ieta, iphi));
0534           float offlineenergy = tower->get_energy();
0535 
0536           if (!tower->get_isGood())
0537           {
0538             continue;
0539           }
0540           energy_sum += offlineenergy;
0541         }
0542       }
0543       h_hcalin_2x2_energy_lutsum->Fill(it.second, energy_sum);
0544     }
0545     // go through the emulated 2x2 map for emcal
0546     for (auto& it : v_hcalout_emu_2x2)
0547     {
0548       unsigned int sumk = it.first;
0549       uint16_t sum_phi = TriggerDefs::getSumPhiId(sumk) + 2 * TriggerDefs::getPrimitivePhiId_from_TriggerSumKey(sumk);
0550       uint16_t sum_eta = TriggerDefs::getSumEtaId(sumk);
0551 
0552       float energy_sum = 0.0;
0553       for (int itower = 0; itower < 4; itower++)
0554       {
0555         int ieta = sum_eta * 2 + itower % 2;
0556         int iphi = sum_phi * 2 + itower / 2;
0557 
0558         if (towers_hcalout)
0559         {
0560           TowerInfo* tower = towers_hcalout->get_tower_at_key(TowerInfoDefs::encode_hcal(ieta, iphi));
0561           float offlineenergy = tower->get_energy();
0562           if (!tower->get_isGood())
0563           {
0564             continue;
0565           }
0566           energy_sum += offlineenergy;
0567         }
0568       }
0569       h_hcalout_2x2_energy_lutsum->Fill(it.second, energy_sum);
0570     }
0571   }
0572 
0573   float jet_energies[9][32]{};
0574   float max_energy_jetpatch = 0.0;
0575   for (int i = 0; i < 3; i++)
0576   {
0577     for (int j = 0; j < 12; j++)
0578     {
0579       emcal_energies[j][i + 32] = emcal_energies[j][i];
0580       hcal_energies[j][i + 32] = hcal_energies[j][i];
0581     }
0582   }
0583   for (int i = 0; i < 9; i++)
0584   {
0585     for (int j = 0; j < 32; j++)
0586     {
0587       for (int k = 0; k < 16; k++)
0588       {
0589         jet_energies[i][j] += emcal_energies[i + k % 4][j + k / 4];
0590         jet_energies[i][j] += hcal_energies[i + k % 4][j + k / 4];
0591       }
0592     }
0593   }
0594 
0595   for (auto& it : v_jet_emu)
0596   {
0597     unsigned int sumk = it.first;
0598     uint16_t sum_phi = sumk & 0xffffU;
0599     uint16_t sum_eta = (sumk >> 16U) & 0xffffU;
0600 
0601     if (max_energy_jetpatch < jet_energies[sum_eta][sum_phi])
0602     {
0603       max_energy_jetpatch = jet_energies[sum_eta][sum_phi];
0604     }
0605     h_jet_energy_lutsum->Fill(it.second, jet_energies[sum_eta][sum_phi]);
0606   }
0607 
0608   for (auto& j : trig_bits)
0609   {
0610     if (j >= 16 && j < 24)
0611     {
0612       h_gl1_jetpatch_energy[j - 16]->Fill(max_energy_jetpatch);
0613     }
0614     else if (j >= 24 && j < 32)
0615     {
0616       h_gl1_photon_energy[j - 24]->Fill(max_energy_emcal);
0617     }
0618   }
0619   return Fun4AllReturnCodes::EVENT_OK;
0620 }