Back to home page

sPhenix code displayed by LXR

 
 

    


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

0001 #include "MBDTriggerEmulator.h"
0002 
0003 #include "LL1Outv1.h"
0004 #include "TriggerPrimitive.h"
0005 #include "TriggerPrimitiveContainerv1.h"
0006 
0007 #include <ffarawobjects/CaloPacket.h>
0008 #include <ffarawobjects/CaloPacketContainer.h>
0009 
0010 #include <ffamodules/CDBInterface.h>
0011 
0012 #include <cdbobjects/CDBHistos.h>  // for CDBHistos
0013 
0014 #include <fun4all/Fun4AllReturnCodes.h>
0015 
0016 #include <phool/PHCompositeNode.h>
0017 #include <phool/PHIODataNode.h>
0018 #include <phool/PHNode.h>
0019 #include <phool/PHNodeIterator.h>
0020 #include <phool/PHObject.h>
0021 #include <phool/getClass.h>
0022 #include <phool/phool.h>
0023 
0024 #include <Event/Event.h>
0025 #include <Event/EventTypes.h>
0026 #include <Event/packet.h>
0027 
0028 #include <TH1.h>
0029 
0030 #include <algorithm>
0031 #include <bitset>
0032 #include <cstdint>
0033 #include <cstdlib>
0034 #include <iostream>
0035 #include <string>
0036 #include <utility>
0037 
0038 // constructor
0039 MBDTriggerEmulator::MBDTriggerEmulator(const std::string &name)
0040   : SubsysReco(name)
0041   , m_nevent(0)
0042   , m_mbd_npassed(0)
0043   , m_isdata(false)
0044   , m_nsamples(16)
0045   , m_idx(8)
0046 {
0047   // initialize all important counters
0048 
0049   // default nsamples is 16 for mbd, 12 for calos
0050 
0051   // reset variables
0052   for (int j = 0; j < 8; j++)
0053   {
0054     m_trig_charge[j] = 0;
0055     for (auto &k : m2_trig_charge)
0056     {
0057       k[j] = 0;
0058     }
0059   }
0060   m_trig_nhit = 0;
0061   for (unsigned int &k : m2_trig_nhit)
0062   {
0063     k = 0;
0064   }
0065   for (int j = 0; j < 4; j++)
0066   {
0067     m_trig_time[j] = 0;
0068     for (int k = 0; k < 4; k++)
0069     {
0070       m2_trig_time[j][k] = 0;
0071     }
0072   }
0073 
0074   for (int i = 0; i < 2; i++)
0075   {
0076     m_out_tsum[i] = 0;
0077     m_out_nhit[i] = 0;
0078     m_out_tavg[i] = 0;
0079     m_out_trem[i] = 0;
0080   }
0081   m_out_vtx_sub = 0;
0082   m_out_vtx_add = 0;
0083 
0084   for (uint16_t i = 0; i < 128; i++)
0085   {
0086     uint16_t key = (i & 0x7fffU);
0087     h_mbd_charge_lut[key] = nullptr;
0088     key = (0x1U << 0xfU) + (i & 0x7fffU);
0089     h_mbd_time_lut[key] = nullptr;
0090     h_mbd_slewing_lut[key] = nullptr;
0091   }
0092 
0093   m_ll1out_mbd = nullptr;
0094 
0095   m_n_primitives = 4;
0096   m_n_sums = 13;
0097   m_trig_sample = -1;
0098   m_trig_sub_delay = 4;
0099   m_threshold = 1;
0100 
0101   m_nhit1 = 1;
0102   m_nhit2 = 2;
0103   m_timediff1 = 10;
0104   m_timediff2 = 20;
0105   m_timediff3 = 30;
0106 }
0107 
0108 // make file and histomanager (but nothing goes in the file at the moment)
0109 int MBDTriggerEmulator::Init(PHCompositeNode * /*topNode*/)
0110 {
0111   return 0;
0112 }
0113 
0114 // at the beginning of the run
0115 int MBDTriggerEmulator::InitRun(PHCompositeNode *topNode)
0116 {
0117   // Get the detectors that are used for a given trigger.
0118 
0119   m_ll1_nodename = "LL1OUT_MBD";
0120   m_prim_nodename = "TRIGGERPRIMITIVES_MBD";
0121 
0122   // Get the calibrations and produce the lookup tables;
0123 
0124   if (Download_Calibrations())
0125   {
0126     std::cout << "here" << std::endl;
0127     return Fun4AllReturnCodes::ABORTRUN;
0128   }
0129 
0130   CreateNodes(topNode);
0131 
0132   return 0;
0133 }
0134 
0135 int MBDTriggerEmulator::Download_Calibrations()
0136 {
0137   if (!m_mbd_charge_lutname.empty())
0138   {
0139     cdbttree_mbd_charge = new CDBHistos(m_mbd_charge_lutname);
0140   }
0141   else
0142   {
0143     std::string calibdir = CDBInterface::instance()->getUrl("mbd_charge_trigger_lut");
0144     if (calibdir.empty())
0145     {
0146       m_default_lut_mbd = true;
0147       std::cout << "Could not find and load histograms for MBD LUTs! defaulting to the identity table!" << std::endl;
0148     }
0149     else
0150     {
0151       cdbttree_mbd_charge = new CDBHistos(calibdir);
0152     }
0153   }
0154   if (cdbttree_mbd_charge)
0155   {
0156     cdbttree_mbd_charge->LoadCalibrations();
0157 
0158     for (uint16_t i = 0; i < 128; i++)
0159     {
0160       std::string histoname = "h_mbd_charge_lut_" + std::to_string(i);
0161       uint16_t key = (i & 0x7fffU);
0162       h_mbd_charge_lut[key] = (TH1I *) cdbttree_mbd_charge->getHisto(histoname);
0163     }
0164   }
0165   else
0166   {
0167     return Fun4AllReturnCodes::ABORTRUN;
0168   }
0169   // now time channels
0170   if (!m_mbd_time_lutname.empty())
0171   {
0172     cdbttree_mbd_time = new CDBHistos(m_mbd_time_lutname);
0173   }
0174   else
0175   {
0176     std::string calibdir = CDBInterface::instance()->getUrl("mbd_time_trigger_lut");
0177     if (calibdir.empty())
0178     {
0179       m_default_lut_mbd = true;
0180       std::cout << "Could not find and load histograms for MBD LUTs! defaulting to the identity table!" << std::endl;
0181     }
0182     else
0183     {
0184       cdbttree_mbd_time = new CDBHistos(calibdir);
0185     }
0186   }
0187   if (cdbttree_mbd_time)
0188   {
0189     cdbttree_mbd_time->LoadCalibrations();
0190 
0191     for (uint16_t i = 0; i < 128; i++)
0192     {
0193       std::string histoname = "h_mbd_time_lut_" + std::to_string(i);
0194       uint16_t key = (i & 0x7fffU) + (0x1U << 0xfU);
0195       h_mbd_time_lut[key] = (TH1I *) cdbttree_mbd_time->getHisto(histoname);
0196     }
0197   }
0198   else
0199   {
0200     return Fun4AllReturnCodes::ABORTRUN;
0201   }
0202 
0203   if (!m_mbd_slewing_lutname.empty())
0204   {
0205     cdbttree_mbd_slewing = new CDBHistos(m_mbd_slewing_lutname);
0206   }
0207   else
0208   {
0209     std::string calibdir = CDBInterface::instance()->getUrl("mbd_slewing_trigger_lut");
0210     if (calibdir.empty())
0211     {
0212       m_default_lut_mbd = true;
0213       std::cout << "Could not find and load histograms for MBD LUTs! defaulting to the identity table!" << std::endl;
0214     }
0215     else
0216     {
0217       cdbttree_mbd_slewing = new CDBHistos(calibdir);
0218     }
0219   }
0220   if (cdbttree_mbd_slewing)
0221   {
0222     cdbttree_mbd_slewing->LoadCalibrations();
0223 
0224     for (uint16_t i = 0; i < 128; i++)
0225     {
0226       std::string histoname = "h_mbd_slewing_lut_" + std::to_string(i);
0227       uint16_t key = (i & 0x7fffU) + (0x1U << 0xfU);
0228       h_mbd_slewing_lut[key] = (TH1I *) cdbttree_mbd_slewing->getHisto(histoname);
0229     }
0230   }
0231   else
0232   {
0233     return Fun4AllReturnCodes::ABORTRUN;
0234   }
0235   return 0;
0236 }
0237 // process event procedure
0238 int MBDTriggerEmulator::process_event(PHCompositeNode *topNode)
0239 {
0240   if (Verbosity() >= 1)
0241   {
0242     std::cout << __FUNCTION__ << ": event " << m_nevent << std::endl;
0243   }
0244 
0245   // Get all nodes needed fo
0246   GetNodes(topNode);
0247 
0248   // process waveforms from the waveform container into primitives
0249   if (process_waveforms())
0250   {
0251     return Fun4AllReturnCodes::EVENT_OK;
0252   }
0253 
0254   if (Verbosity())
0255   {
0256     std::cout << __FILE__ << __FUNCTION__ << __LINE__ << "::"
0257               << "done with waveforms" << std::endl;
0258   }
0259   // process all the primitives into sums.
0260   process_primitives();
0261 
0262   // calculate the true LL1 trigger at emcal and hcal.
0263 
0264   // calculate the true LL1 trigger algorithm.
0265   if (process_trigger())
0266   {
0267     return Fun4AllReturnCodes::EVENT_OK;
0268   }
0269 
0270   m_nevent++;
0271 
0272   if (Verbosity() >= 2)
0273   {
0274     identify();
0275   }
0276 
0277   return Fun4AllReturnCodes::EVENT_OK;
0278 }
0279 
0280 // RESET event procedure that takes all variables to 0 and clears the primitives.
0281 int MBDTriggerEmulator::ResetEvent(PHCompositeNode * /*topNode*/)
0282 {
0283   m_peak_sub_ped_mbd.clear();
0284   return 0;
0285 }
0286 
0287 int MBDTriggerEmulator::process_waveforms()
0288 {
0289   // Get range of waveforms
0290   if (Verbosity())
0291   {
0292     std::cout << __FILE__ << "::" << __FUNCTION__ << ":: Processing waveforms" << std::endl;
0293   }
0294 
0295   if (!m_waveforms_mbd)
0296   {
0297     return process_raw();
0298   }
0299   if (Verbosity())
0300   {
0301     std::cout << __LINE__ << std::endl;
0302   }
0303   int sample_start = 1;
0304   int sample_end = m_nsamples;
0305   if (m_trig_sample > 0)
0306   {
0307     sample_start = m_trig_sample;
0308     sample_end = m_trig_sample + 1;
0309   }
0310   if (Verbosity())
0311   {
0312     std::cout << __LINE__ << std::endl;
0313   }
0314 
0315   // for each waveform, clauclate the peak - pedestal given the sub-delay setting
0316   CaloPacket *dstp[2]{nullptr};
0317 
0318   for (int ipkt = 0; ipkt < 2; ipkt++)
0319   {
0320     if (Verbosity())
0321     {
0322       std::cout << __LINE__ << std::endl;
0323     }
0324 
0325     int pktid = 1001 + ipkt;  // packet id
0326 
0327     dstp[ipkt] = m_waveforms_mbd->getPacketbyId(pktid);
0328 
0329     if (dstp[ipkt])
0330     {
0331       for (int ich = 0; ich < 128; ich++)
0332       {
0333         //        int feech = ipkt * 128 + ich;
0334         uint16_t isTime = (ich % 16 < 8 ? 1 : 0);
0335 
0336         uint16_t ipmt = (ipkt * 64) + ((ich / 16) * 8) + (ich % 8);
0337         uint16_t key = ((isTime & 0x1U) << 0xfU) + (ipmt & 0x7fffU);
0338 
0339         std::vector<unsigned int> v_peak_sub_ped;
0340         unsigned int peak_sub_ped = 0;
0341 
0342         for (int i = sample_start; i < sample_end; i++)
0343         {
0344           int16_t maxim = dstp[ipkt]->iValue(i, ich);  // > dstp[ipkt]->iValue(i+1, ich) ? dstp[ipkt]->iValue(i, ich) : dstp[ipkt]->iValue(i+1, ich));
0345           //          maxim = (maxim > dstp[ipkt]->iValue(i+2, ich) ? maxim : dstp[ipkt]->iValue(i+2, ich));
0346 
0347           if (Verbosity())
0348           {
0349             std::cout << __LINE__ << std::endl;
0350           }
0351 
0352           int subtraction = maxim - dstp[ipkt]->iValue(((i - m_trig_sub_delay > 0 ? i - m_trig_sub_delay : 0)), ich);
0353 
0354           subtraction = std::max(subtraction, 0);
0355 
0356           peak_sub_ped = (((unsigned int) subtraction) & 0x3fffU);
0357           v_peak_sub_ped.push_back(peak_sub_ped);
0358         }
0359         // save in global.
0360         m_peak_sub_ped_mbd[key] = v_peak_sub_ped;
0361       }
0362       if (Verbosity())
0363       {
0364         std::cout << __LINE__ << std::endl;
0365       }
0366       delete dstp[ipkt];
0367     }
0368 
0369     if (Verbosity())
0370     {
0371       std::cout << __LINE__ << std::endl;
0372     }
0373   }
0374   if (Verbosity())
0375   {
0376     std::cout << __LINE__ << std::endl;
0377   }
0378 
0379   return Fun4AllReturnCodes::EVENT_OK;
0380 }
0381 int MBDTriggerEmulator::process_raw()
0382 {
0383   // Get range of waveforms
0384   if (Verbosity())
0385   {
0386     std::cout << __FILE__ << "::" << __FUNCTION__ << ":: Processing waveforms" << std::endl;
0387   }
0388 
0389   if (m_useoffline)
0390   {
0391     return Fun4AllReturnCodes::ABORTEVENT;
0392   }
0393 
0394   if (m_event == nullptr)
0395   {
0396     std::cout << PHWHERE << " Event not found" << std::endl;
0397     return Fun4AllReturnCodes::ABORTEVENT;
0398   }
0399 
0400   if (m_event->getEvtType() != DATAEVENT)
0401   {
0402     return Fun4AllReturnCodes::ABORTEVENT;
0403   }
0404 
0405   int sample_start = 1;
0406   int sample_end = m_nsamples;
0407   if (m_trig_sample > 0)
0408   {
0409     sample_start = m_trig_sample;
0410     sample_end = m_trig_sample + 1;
0411   }
0412 
0413   // for each waveform, clauclate the peak - pedestal given the sub-delay setting
0414   Packet *dstp[2]{nullptr};
0415 
0416   for (int ipkt = 0; ipkt < 2; ipkt++)
0417   {
0418     int pktid = 1001 + ipkt;  // packet id
0419 
0420     dstp[ipkt] = m_event->getPacket(pktid);
0421 
0422     if (dstp[ipkt])
0423     {
0424       for (int ich = 0; ich < 128; ich++)
0425       {
0426         //        int feech = ipkt * 128 + ich;
0427         uint16_t isTime = (ich % 16 < 8 ? 1 : 0);
0428 
0429         uint16_t ipmt = (ipkt * 64) + ((ich / 16) * 8) + (ich % 8);
0430         uint16_t key = ((isTime & 0x1U) << 0xfU) + (ipmt & 0x7fffU);
0431 
0432         std::vector<unsigned int> v_peak_sub_ped;
0433         unsigned int peak_sub_ped = 0;
0434 
0435         for (int i = sample_start; i < sample_end; i++)
0436         {
0437           int16_t maxim = dstp[ipkt]->iValue(i, ich);  // > dstp[ipkt]->iValue(i+1, ich) ? dstp[ipkt]->iValue(i, ich) : dstp[ipkt]->iValue(i+1, ich));
0438           //          maxim = (maxim > dstp[ipkt]->iValue(i+2, ich) ? maxim : dstp[ipkt]->iValue(i+2, ich));
0439 
0440           int subtraction = maxim - dstp[ipkt]->iValue(((i - m_trig_sub_delay > 0 ? i - m_trig_sub_delay : 0)), ich);
0441 
0442           subtraction = std::max(subtraction, 0);
0443 
0444           peak_sub_ped = (((unsigned int) subtraction) & 0x3fffU);
0445           v_peak_sub_ped.push_back(peak_sub_ped);
0446         }
0447         // save in global.
0448         m_peak_sub_ped_mbd[key] = v_peak_sub_ped;
0449       }
0450     }
0451     delete dstp[ipkt];
0452   }
0453 
0454   return Fun4AllReturnCodes::EVENT_OK;
0455 }
0456 
0457 // procedure to process the peak - pedestal into primitives.
0458 int MBDTriggerEmulator::process_primitives()
0459 {
0460   int i = 0;
0461   int nsample = m_nsamples - 1;
0462   if (m_trig_sample > 0)
0463   {
0464     nsample = 1;
0465   }
0466 
0467   if (Verbosity())
0468   {
0469     std::cout << __FILE__ << "::" << __FUNCTION__ << ":: Processing primitives" << std::endl;
0470   }
0471 
0472   for (i = 0; i < m_n_primitives; i++)
0473   {
0474     // make primitive key
0475     TriggerDefs::TriggerPrimKey primkey = TriggerDefs::getTriggerPrimKey(m_triggerid, TriggerDefs::GetDetectorId("MBD"), TriggerDefs::GetPrimitiveId("MBD"), i);
0476 
0477     // make primitive and check mask;
0478     TriggerPrimitive *primitive = m_primitives_mbd->get_primitive_at_key(primkey);
0479 
0480     // iterate through samples
0481     for (int is = 0; is < nsample; is++)
0482     {
0483       // reset variables
0484       for (unsigned int &j : m_trig_charge)
0485       {
0486         j = 0;
0487       }
0488       m_trig_nhit = 0;
0489       for (unsigned int &j : m_trig_time)
0490       {
0491         j = 0;
0492       }
0493 
0494       unsigned int tmp;
0495       unsigned int tmp2;
0496       unsigned int qadd[32];
0497 
0498       // for each section of the board (4 sections of 8 time and 8 charge
0499       // go through 8 charge channels
0500       for (int j = 0; j < 32; j++)
0501       {
0502         // pass upper 10 bits of charge to get 10 bit LUt outcome
0503         // tmp = m_l1_adc_table[m_peak_sub_ped_mbd[i * 64 + 8 + isec * 16 + j].at(is) >> 4U];
0504         uint16_t key = static_cast<uint16_t>((i * 32) + j) & 0x7fffU;
0505         unsigned int peak_minus_ped = m_peak_sub_ped_mbd[key].at(is) >> 0x4U;
0506         tmp = h_mbd_charge_lut[key]->GetBinContent(peak_minus_ped + 1);
0507 
0508         // put upper 3 bits of the 10 bits into slewing correction later
0509         qadd[j] = (tmp & 0x380U) >> 7U;
0510 
0511         // sum up to 11 bits.
0512         m_trig_charge[j / 4] += tmp & 0x7fU;
0513       }
0514 
0515       // 8 timing channels
0516       for (int j = 0; j < 32; j++)
0517       {
0518         // pass upper 10 bits of charge to get 10 bit LUt outcome
0519         // tmp = m_l1_adc_table[m_peak_sub_ped_mbd[i * 64 + 8 + isec * 16 + j].at(is) >> 4U];
0520         uint16_t key = (static_cast<uint16_t>((i * 32) + j) & 0x7fffU) + (0x1U << 0xfU);
0521         unsigned int peak_minus_ped = m_peak_sub_ped_mbd[key].at(is) >> 0x4U;
0522         tmp = h_mbd_time_lut[key]->GetBinContent(peak_minus_ped + 1);
0523 
0524         m_trig_nhit += (tmp & 0x200U) >> 9U;
0525         unsigned int slewkey = (qadd[j] << 9U) + (tmp & 0x01ffU);
0526         // get upper 3 bits of charge in the channel, and make it bits 9-11, the time of the chanel is the lower 9 bits from 0-8.
0527         tmp2 = h_mbd_slewing_lut[key]->GetBinContent(slewkey + 1);
0528         // attribute to the time sum
0529         m_trig_time[j / 8] += tmp2;
0530       }
0531       // ad in the charge sums
0532 
0533       for (int j = 0; j < 13; j++)
0534       {
0535         TriggerDefs::TriggerSumKey sumkey = TriggerDefs::getTriggerSumKey(m_triggerid, TriggerDefs::GetDetectorId("MBD"), TriggerDefs::GetPrimitiveId("MBD"), i, j);
0536         if (j < 8)
0537         {
0538           primitive->get_sum_at_key(sumkey)->push_back(m_trig_charge[j]);
0539         }
0540         else if (j == 8)
0541         {
0542           primitive->get_sum_at_key(sumkey)->push_back(m_trig_nhit);
0543         }
0544         else
0545         {
0546           primitive->get_sum_at_key(sumkey)->push_back(m_trig_time[j - 9]);
0547         }
0548       }
0549     }
0550   }
0551 
0552   return Fun4AllReturnCodes::EVENT_OK;
0553 }
0554 
0555 // Unless this is the MBD or HCAL Cosmics trigger, EMCAL and HCAL will go through here.
0556 // This creates the 8x8 non-overlapping sum and the 4x4 overlapping sum.
0557 
0558 // This is where the LL1 Jet/Pair/Cosmic algorithm is
0559 
0560 int MBDTriggerEmulator::process_trigger()
0561 {
0562   std::vector<unsigned int> bits;
0563 
0564   int nsample = m_nsamples - 1;
0565   // bits are to say whether the trigger has fired. this is what is sent to the GL1
0566   if (m_trig_sample > 0)
0567   {
0568     nsample = 1;
0569   }
0570   bits.reserve(nsample);
0571   for (int is = 0; is < nsample; is++)
0572   {
0573     bits.push_back(0);
0574   }
0575 
0576   std::vector<unsigned int> *trig_bits = m_ll1out_mbd->GetTriggerBits();
0577   if (!m_primitives_mbd)
0578   {
0579     std::cout << "There is no primitive container" << std::endl;
0580     return Fun4AllReturnCodes::EVENT_OK;
0581   }
0582 
0583   TriggerPrimitive::Range sumrange;
0584   int ip;
0585   int isum;
0586 
0587   TriggerPrimitiveContainer::Range range = m_primitives_mbd->getTriggerPrimitives();
0588 
0589   if (Verbosity() >= 2)
0590   {
0591     std::cout << __FUNCTION__ << " " << __LINE__ << " mbd primitives size: " << m_primitives_mbd->size() << std::endl;
0592   }
0593 
0594   std::vector<unsigned int> *word_mbd = nullptr;
0595 
0596   m_word_mbd.clear();
0597   for (int j = 0; j < 8; j++)
0598   {
0599     word_mbd = new std::vector<unsigned int>();
0600     m_word_mbd.push_back(word_mbd);
0601   }
0602 
0603   for (int is = 0; is < nsample; is++)
0604   {
0605     ip = 0;
0606     for (TriggerPrimitiveContainer::Iter iter = range.first; iter != range.second; ++iter, ip++)
0607     {
0608       TriggerPrimitive *primitive = (*iter).second;
0609       sumrange = primitive->getSums();
0610       isum = 0;
0611       for (TriggerPrimitive::Iter iter_sum = sumrange.first; iter_sum != sumrange.second; ++iter_sum, isum++)
0612       {
0613         if (isum < 8)
0614         {
0615           m2_trig_charge[ip][isum] = ((*iter_sum).second)->at(is);
0616         }
0617         else if (isum == 8)
0618         {
0619           m2_trig_nhit[ip] = ((*iter_sum).second)->at(is);
0620         }
0621         else
0622         {
0623           m2_trig_time[ip][isum - 9] = ((*iter_sum).second)->at(is);
0624         }
0625       }
0626     }
0627 
0628     if (Verbosity() && is == 11)
0629     {
0630       for (int q = 0; q < 8; q++)
0631       {
0632         std::cout << "Q" << std::dec << q << ": ";
0633         for (auto &ipp : m2_trig_charge)
0634         {
0635           std::cout << std::hex << ipp[q] << " ";
0636         }
0637         std::cout << " " << std::endl;
0638       }
0639       std::cout << "NH: ";
0640       for (unsigned int ipp : m2_trig_nhit)
0641       {
0642         std::cout << std::hex << ipp << " ";
0643       }
0644       std::cout << " " << std::endl;
0645 
0646       for (int q = 0; q < 4; q++)
0647       {
0648         std::cout << "T" << std::dec << q << ": ";
0649         for (auto &ipp : m2_trig_time)
0650         {
0651           std::cout << std::hex << ipp[q] << " ";
0652         }
0653         std::cout << " " << std::endl;
0654       }
0655     }
0656 
0657     m_out_tsum[0] = 0;
0658     m_out_tsum[1] = 0;
0659     m_out_nhit[0] = 0;
0660     m_out_nhit[1] = 0;
0661     m_out_tavg[0] = 0;
0662     m_out_tavg[1] = 0;
0663     m_out_trem[0] = 0;
0664     m_out_trem[1] = 0;
0665     m_out_vtx_sub = 0;
0666     m_out_vtx_add = 0;
0667 
0668     for (int i = 0; i < 2; i++)
0669     {
0670       for (int j = 0; j < 4; j++)
0671       {
0672         m_out_tsum[0] += m2_trig_time[i][j];
0673         m_out_tsum[1] += m2_trig_time[i + 2][j];
0674       }
0675       m_out_nhit[0] += m2_trig_nhit[i];
0676       m_out_nhit[1] += m2_trig_nhit[i + 2];
0677     }
0678 
0679     if (m_out_nhit[0] != 0)
0680     {
0681       m_out_tavg[0] = m_out_tsum[0] / m_out_nhit[0];
0682       m_out_trem[0] = m_out_tsum[0] % m_out_nhit[0];
0683     }
0684     if (m_out_nhit[1] != 0)
0685     {
0686       m_out_tavg[1] = m_out_tsum[1] / m_out_nhit[1];
0687       m_out_trem[1] = m_out_tsum[1] % m_out_nhit[1];
0688     }
0689 
0690     unsigned int max = m_out_tavg[0];
0691     unsigned int min = m_out_tavg[1];
0692     if (min > max)
0693     {
0694       max = m_out_tavg[1];
0695       min = m_out_tavg[0];
0696     }
0697 
0698     m_out_vtx_sub = (max - min) & 0x1ffU;
0699     m_out_vtx_add = (m_out_tavg[0] + m_out_tavg[1]) & 0x3ffU;
0700 
0701     m_word_mbd[0]->push_back(m_out_tavg[0]);
0702     m_word_mbd[1]->push_back(m_out_tavg[1]);
0703     m_word_mbd[2]->push_back(m_out_nhit[0]);
0704     m_word_mbd[3]->push_back(m_out_nhit[1]);
0705     m_word_mbd[4]->push_back(m_out_trem[0]);
0706     m_word_mbd[5]->push_back(m_out_trem[1]);
0707     m_word_mbd[6]->push_back(m_out_vtx_sub);
0708     m_word_mbd[7]->push_back(m_out_vtx_add);
0709 
0710     if (m_out_nhit[0] >= m_nhit1)
0711     {
0712       bits.at(is) ^= 1U << 0U;
0713     }
0714     if (m_out_nhit[1] >= m_nhit1)
0715     {
0716       bits.at(is) ^= 1U << 1U;
0717     }
0718     if (m_out_nhit[0] >= m_nhit2)
0719     {
0720       bits.at(is) ^= 1U << 2U;
0721     }
0722     if (m_out_nhit[1] >= m_nhit2)
0723     {
0724       bits.at(is) ^= 1U << 3U;
0725     }
0726 
0727     if (m_out_nhit[0] >= m_nhit1 && m_out_nhit[1] >= m_nhit1 && m_out_vtx_sub <= m_timediff1)
0728     {
0729       bits.at(is) ^= 1U << 4U;
0730     }
0731     if (m_out_nhit[0] >= m_nhit1 && m_out_nhit[1] >= m_nhit1 && m_out_vtx_sub <= m_timediff2)
0732     {
0733       bits.at(is) ^= 1U << 5U;
0734     }
0735     if (m_out_nhit[0] >= m_nhit1 && m_out_nhit[1] >= m_nhit1 && m_out_vtx_sub <= m_timediff3)
0736     {
0737       bits.at(is) ^= 1U << 6U;
0738     }
0739     if (m_out_nhit[0] >= m_nhit2 && m_out_nhit[1] >= m_nhit2 && m_out_vtx_sub <= m_timediff1)
0740     {
0741       bits.at(is) ^= 1U << 7U;
0742     }
0743     if (m_out_nhit[0] >= m_nhit2 && m_out_nhit[1] >= m_nhit2 && m_out_vtx_sub <= m_timediff2)
0744     {
0745       bits.at(is) ^= 1U << 8U;
0746     }
0747     if (m_out_nhit[0] >= m_nhit2 && m_out_nhit[1] >= m_nhit2 && m_out_vtx_sub <= m_timediff3)
0748     {
0749       bits.at(is) ^= 1U << 9U;
0750     }
0751 
0752     if (Verbosity())
0753     {
0754       std::cout << "Trigger Word : " << std::bitset<16>(bits.at(is)) << std::dec << std::endl;
0755     }
0756   }
0757   int pass = 0;
0758   for (int is = 0; is < nsample; is++)
0759   {
0760     trig_bits->push_back(bits.at(is));
0761     if (bits.at(is))
0762     {
0763       pass = 1;
0764     }
0765   }
0766   if (pass)
0767   {
0768     m_mbd_npassed++;
0769   }
0770   if (Verbosity() >= 2)
0771   {
0772     std::cout << " " << std::endl;
0773   }
0774 
0775   for (int iw = 0; iw < 8; iw++)
0776   {
0777     std::vector<unsigned int> *sum = m_ll1out_mbd->get_word(iw);
0778     for (int is = 0; is < nsample; is++)
0779     {
0780       sum->push_back(m_word_mbd[iw]->at(is));
0781     }
0782   }
0783 
0784   return Fun4AllReturnCodes::EVENT_OK;
0785 }
0786 
0787 void MBDTriggerEmulator::GetNodes(PHCompositeNode *topNode)
0788 {
0789   if (Verbosity() >= 2)
0790   {
0791     std::cout << __FUNCTION__ << std::endl;
0792   }
0793 
0794   m_event = findNode::getClass<Event>(topNode, "PRDF");
0795 
0796   m_ll1_nodename = "LL1OUT_MBD";
0797 
0798   m_ll1out_mbd = findNode::getClass<LL1Out>(topNode, m_ll1_nodename);
0799 
0800   if (!m_ll1out_mbd)
0801   {
0802     std::cout << "No LL1Out found... " << std::endl;
0803     exit(1);
0804   }
0805 
0806   m_prim_nodename = "TRIGGERPRIMITIVES_MBD";
0807   m_primitives_mbd = findNode::getClass<TriggerPrimitiveContainer>(topNode, m_prim_nodename);
0808 
0809   if (!m_primitives_mbd)
0810   {
0811     std::cout << "No TriggerPrimitives found... " << std::endl;
0812     exit(1);
0813   }
0814 
0815   m_waveforms_mbd = findNode::getClass<CaloPacketContainer>(topNode, "MBDPackets");
0816   if (m_waveforms_mbd)
0817   {
0818     m_useoffline = true;
0819   }
0820   // if (!m_waveforms_mbd)
0821   //   {
0822   //     std::cout << "No MBD Waveforms found... " << std::endl;
0823   //     exit(1);
0824   //   }
0825 
0826   return;
0827 }
0828 void MBDTriggerEmulator::CreateNodes(PHCompositeNode *topNode)
0829 {
0830   PHNodeIterator iter(topNode);
0831   PHCompositeNode *dstNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "DST"));
0832   if (!dstNode)
0833   {
0834     std::cout << PHWHERE << "DST Node missing doing nothing" << std::endl;
0835   }
0836 
0837   PHCompositeNode *ll1Node = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "LL1"));
0838   if (!ll1Node)
0839   {
0840     ll1Node = new PHCompositeNode("LL1");
0841     dstNode->addNode(ll1Node);
0842   }
0843 
0844   m_ll1_nodename = "LL1OUT_MBD";
0845   LL1Out *ll1out = findNode::getClass<LL1Out>(ll1Node, m_ll1_nodename);
0846   if (!ll1out)
0847   {
0848     ll1out = new LL1Outv1("MBD", "NONE");
0849     PHIODataNode<PHObject> *LL1OutNode = new PHIODataNode<PHObject>(ll1out, m_ll1_nodename, "PHObject");
0850     ll1Node->addNode(LL1OutNode);
0851   }
0852   m_prim_nodename = "TRIGGERPRIMITIVES_MBD";
0853   TriggerPrimitiveContainer *ll1out_prim = findNode::getClass<TriggerPrimitiveContainer>(ll1Node, m_prim_nodename);
0854   if (!ll1out_prim)
0855   {
0856     ll1out_prim = new TriggerPrimitiveContainerv1(TriggerDefs::TriggerId::mbdTId, TriggerDefs::DetectorId::mbdDId);
0857     PHIODataNode<PHObject> *LL1OutNode = new PHIODataNode<PHObject>(ll1out_prim, m_prim_nodename, "PHObject");
0858     ll1Node->addNode(LL1OutNode);
0859   }
0860   return;
0861 }
0862 
0863 int MBDTriggerEmulator::End(PHCompositeNode * /*topNode*/)
0864 {
0865   delete cdbttree_mbd_charge;
0866   delete cdbttree_mbd_time;
0867   delete cdbttree_mbd_slewing;
0868 
0869   std::cout << "------------------------" << std::endl;
0870   std::cout << "Total MBD passed: " << m_mbd_npassed << "/" << m_nevent << std::endl;
0871   std::cout << "------------------------" << std::endl;
0872 
0873   return 0;
0874 }
0875 
0876 void MBDTriggerEmulator::identify()
0877 {
0878   std::cout << " LL1Out: " << std::endl;
0879   m_ll1out_mbd->identify();
0880   std::cout << " TriggerPrimitives: " << std::endl;
0881   m_primitives_mbd->identify();
0882   std::cout << " Processed " << m_nevent << std::endl;
0883 }