Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:20:26

0001 #include "LL1HEADER.h"
0002 
0003 #include <onlmon/OnlMon.h>
0004 
0005 // Event library
0006 
0007 #include <Event/Event.h>
0008 #include <Event/EventTypes.h>
0009 #include <Event/msg_profile.h>
0010 
0011 // ROOT
0012 
0013 #include <TFile.h>
0014 #include <TTree.h>
0015 #include <TSystem.h>
0016 
0017 #include <cstdlib>
0018 #include <cmath>
0019 #include <iostream>
0020 #include <fstream>
0021 #include <sstream>
0022 
0023 void process_jet(Packet *p, LL1HEADER *&ll1h);
0024 void process_emcal(int pid, Packet *p, LL1HEADER *&ll1h);
0025 void ll1setup(Event *evt, LL1HEADER *&ll1h)
0026 {
0027 
0028   for ( Int_t ixmit = 0; ixmit < NXMIT; ixmit++ ) 
0029   {
0030 
0031     Packet* p;
0032     Int_t pid= PACKET[ixmit];
0033     p = evt->getPacket( pid );
0034     if ( p ) 
0035     {
0036       ll1h->runnumber = evt->getRunNumber();
0037       ll1h->evtnr = p->iValue(0,"EVTNR");
0038       ll1h->clock = p->lValue(0,"CLOCK");
0039       ll1h->monitor = p->iValue(0,"MONITOR");
0040       ll1h->nsamples = p->iValue(0,"SAMPLES");
0041 
0042       if (pid == 13002)
0043       {
0044     process_jet(p, ll1h);
0045       }
0046       else if (pid >= 13010)
0047       {
0048     process_emcal(pid, p, ll1h);
0049       }            
0050       else if (pid == 13001)
0051       {
0052     int indx1=0; int inhit1=0; 
0053     int indx2=0; int inhit2=0;
0054     for(int is = 0; is < NRSAM; is++)
0055       {
0056         ll1h->nhit_s1[is] =  p->iValue(is,NHITCHANNEL+NADCSH*0);
0057         ll1h->nhit_s2[is] =  p->iValue(is,NHITCHANNEL+NADCSH*1);
0058         ll1h->nhit_n1[is] =  p->iValue(is,NHITCHANNEL+NADCSH*2);
0059         ll1h->nhit_n2[is] =  p->iValue(is,NHITCHANNEL+NADCSH*3);
0060         
0061         ll1h->nhit_n[is] = ll1h->nhit_n1[is] + ll1h->nhit_n2[is];
0062         ll1h->nhit_s[is] = ll1h->nhit_s1[is] + ll1h->nhit_s2[is];
0063 
0064         ll1h->timesum_s1[is] =  p->iValue(is,NHITCHANNEL+1+NADCSH*0);
0065         ll1h->timesum_s2[is] =  p->iValue(is,NHITCHANNEL+2+NADCSH*1);
0066         ll1h->timesum_n1[is] =  p->iValue(is,NHITCHANNEL+3+NADCSH*2);
0067         ll1h->timesum_n2[is] =  p->iValue(is,NHITCHANNEL+4+NADCSH*3);
0068 
0069         ll1h->timesum_s[is] = ll1h->timesum_n1[is] + ll1h->timesum_n2[is];
0070         ll1h->timesum_n[is] = ll1h->timesum_s1[is] + ll1h->timesum_s2[is];
0071 
0072         if( ll1h->nhit_n[is]>0)
0073         {
0074           if(inhit1 <  ll1h->nhit_n[is])
0075           {
0076         indx1 = is;
0077         inhit1 = ll1h->nhit_n[is];
0078           } 
0079         }
0080 
0081         if( ll1h->nhit_s[is]>0)
0082         {
0083           if(inhit2 <  ll1h->nhit_s[is])
0084         {
0085           indx2 = is;
0086           inhit2 = ll1h->nhit_s[is];
0087         } 
0088         }
0089         
0090         ll1h->chargesum_n1[is] = 0;
0091         ll1h->chargesum_n2[is] = 0;
0092         ll1h->chargesum_s1[is] = 0;
0093         ll1h->chargesum_s2[is] = 0;
0094         ll1h->chargesum_s[is] = 0;
0095         ll1h->chargesum_n[is] = 0;
0096 
0097         for(int ic = 0; ic <NCH; ic++) 
0098           {
0099         ll1h->channel[ic][is] = p->iValue(is,ic);
0100         if(ic<NHITCHANNEL){ ll1h->chargesum_s1[is] += p->iValue(is,ic);
0101         }
0102         else if(ic>=NADCSH*1 && ic<NHITCHANNEL+NADCSH*1){
0103           ll1h->chargesum_s2[is] += p->iValue(is,ic);
0104         }
0105         else if(ic>=NADCSH*2 && ic<NHITCHANNEL+NADCSH*2){
0106           ll1h->chargesum_n1[is] += p->iValue(is,ic);
0107         }
0108         else if(ic>=NADCSH*3 && ic<NHITCHANNEL+NADCSH*3){
0109           ll1h->chargesum_n2[is] += p->iValue(is,ic);
0110         }
0111         ll1h->chargesum_n[is] = ll1h->chargesum_n1[is] + ll1h->chargesum_n2[is];
0112         ll1h->chargesum_s[is] = ll1h->chargesum_s1[is] + ll1h->chargesum_s2[is];
0113           }
0114         for(int it=NCH; it<(NCH+NTRIGWORDS); it++)
0115           {
0116         ll1h->triggerwords[it - NCH][is] =  p->iValue(is,it);
0117           }
0118       } 
0119     ll1h->idxhitn = indx1;
0120       ll1h->idxhits = indx2;
0121 
0122       ll1h->idxsample = (ll1h->idxhitn==ll1h->idxhits) ? ll1h->idxhitn : -1;
0123 
0124       delete p;
0125     }
0126     }
0127     }
0128 }
0129 
0130 void process_jet(Packet *p, LL1HEADER *&ll1h)
0131 {
0132   
0133   // go through all input fibers
0134   static int jet_event_counter = 0;
0135   if (jet_event_counter == 0)
0136     {
0137       ll1h->hit_format_jet = p->getHitFormat();
0138     }
0139   for (int i = 0; i < 16; i++)
0140     {
0141       ll1h->jet_sample[i] = -1;
0142       for (int j =0; j < 24; j++)
0143     {
0144       for (int is = 0; is < p->iValue(0, "SAMPLES"); is++)
0145         {
0146           int ieta = j%12;
0147           int iphi = j/12 + i*2;
0148 
0149           int value = p->iValue(is, iphi + 32*ieta);
0150           if (value)
0151         {
0152           ll1h->jet_sample[i] = is;
0153           ll1h->jet_input[ieta][iphi] = value;
0154         }
0155         }
0156     }
0157     }
0158 
0159   for (int i = 0; i < p->iValue(0, "TRIGGERWORDS"); i++)
0160     {
0161       for (int is = 0; is < p->iValue(0,"SAMPLES"); is++)
0162     {
0163       int value = p->iValue(is, 16*24 + i);
0164       if (value)
0165       {
0166         ll1h->jet_output[i/32][i%32] = value;
0167       }
0168     }
0169     }
0170 
0171   return;
0172 }
0173 
0174 void process_emcal(int pid, Packet *p, LL1HEADER *&ll1h)
0175 {
0176   
0177   // go through all input fibers
0178 
0179   int emcal_board = pid - 13010;
0180 
0181   ll1h->emcal_sample[emcal_board] = -1;
0182 
0183   for (int i = 0; i < p->iValue(0, "CHANNELS"); i++)
0184     {
0185       for (int is = 0; is < p->iValue(0, "SAMPLES"); is++)
0186     {
0187       int value = p->iValue(is, i);
0188 
0189       ll1h->emcal_2x2_map[((i/16)%12) * 4 + (i%4)][emcal_board*8 + (i%16)/4] = value;
0190       
0191     }
0192     }
0193 
0194   for (int i = 0; i < p->iValue(0, "TRIGGERWORDS"); i++)
0195     {
0196       for (int is = 0; is < p->iValue(0, "SAMPLES"); is++)
0197     {
0198       int value = p->iValue(is, 16*24 + i);
0199       if (value)
0200         {
0201           ll1h->emcal_sample[emcal_board] = is;
0202           ll1h->emcal_8x8_map[i%12][emcal_board*2 + i/12] = value;
0203         }
0204     }
0205     }
0206   return;
0207 }