Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:11:15

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