Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:14

0001 // -*- C++ -*-
0002 //
0003 // Package:    GctFibreAnalyzer
0004 // Class:      GctFibreAnalyzer
0005 //
0006 /**\class GctFibreAnalyzer GctFibreAnalyzer.cc L1Trigger/L1GctAnalzyer/src/GctFibreAnalyzer.cc
0007 
0008 Description: Analyzer individual fibre channels from the source card.
0009 
0010 */
0011 //
0012 // Original Author:  Alex Tapper
0013 //         Created:  Thu Jul 12 14:21:06 CEST 2007
0014 //
0015 //
0016 
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"  // Logger
0018 #include "FWCore/Utilities/interface/Exception.h"          // Exceptions
0019 
0020 // Include file
0021 #include "L1Trigger/L1GctAnalyzer/interface/GctFibreAnalyzer.h"
0022 
0023 // Data format
0024 #include "DataFormats/L1GlobalCaloTrigger/interface/L1GctCollections.h"
0025 
0026 GctFibreAnalyzer::GctFibreAnalyzer(const edm::ParameterSet& iConfig)
0027     : m_fibreSource(iConfig.getUntrackedParameter<edm::InputTag>("FibreSource")),
0028       m_doLogicalID(iConfig.getUntrackedParameter<bool>("doLogicalID")),
0029       m_doCounter(iConfig.getUntrackedParameter<bool>("doCounter")),
0030       m_numZeroEvents(0),
0031       m_numInconsistentPayloadEvents(0),
0032       m_numConsistentEvents(0) {}
0033 
0034 GctFibreAnalyzer::~GctFibreAnalyzer() {
0035   edm::LogInfo("Zero Fibreword events") << "Total number of events with zero fibrewords: " << m_numZeroEvents;
0036   edm::LogInfo("Inconsistent Payload events")
0037       << "Total number of events with inconsistent payloads: " << m_numInconsistentPayloadEvents;
0038   if (m_doCounter) {
0039     edm::LogInfo("Successful events") << "Total number of Successful events: " << m_numConsistentEvents;
0040   }
0041 }
0042 
0043 void GctFibreAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0044   using namespace edm;
0045 
0046   Handle<L1GctFibreCollection> fibre;
0047   iEvent.getByLabel(m_fibreSource, fibre);
0048 
0049   bool bc0 = false;
0050   int flag_for_zeroes = 0;
0051   int flag_for_inconsistent_events = 0;
0052   unsigned int flag_for_consistency = 0;
0053   int flag_for_consistent_events = 0;
0054 
0055   for (L1GctFibreCollection::const_iterator f = fibre->begin(); f != fibre->end(); f++) {
0056     if (f->data() != 0) {
0057       if (m_doCounter) {
0058         if (f == fibre->begin()) {
0059           flag_for_consistency = f->data();
0060         } else if (flag_for_consistency == f->data()) {
0061           flag_for_consistent_events++;
0062         }
0063       }
0064 
0065       // Check for corrupt fibre data
0066       if (!CheckFibreWord(*f)) {
0067         edm::LogInfo("GCT fibre data error") << "Missing phase bit (clock) in fibre data " << (*f);
0068       }
0069 
0070       // Check for BC0
0071       if (CheckForBC0(*f) && (f == fibre->begin())) {
0072         bc0 = true;
0073       }
0074 
0075       // Check for mismatch between fibres
0076       if ((bc0 && !CheckForBC0(*f)) || (!bc0 && CheckForBC0(*f))) {
0077         edm::LogInfo("GCT fibre data error") << "BC0 mismatch in fibre data " << (*f);
0078       }
0079 
0080       // Check logical ID pattern
0081       if (m_doLogicalID)
0082         CheckLogicalID(*f);
0083 
0084       // Check counter pattern
0085       if (m_doCounter)
0086         CheckCounter(*f);
0087 
0088       flag_for_zeroes = 1;
0089     } else {
0090       //this basically checks that all the data is 0s by the time it gets to the last iteration
0091       if (flag_for_zeroes == 0 && f == (fibre->end() - 1)) {
0092         m_numZeroEvents++;
0093       }
0094       //if the zero flag is set to 1 and it managed to find its way into here then something is wrong!
0095       if (flag_for_zeroes == 1) {
0096         flag_for_inconsistent_events++;
0097       }
0098     }
0099   }
0100   //check for inconsistent events i.e. those with one(or more) zeroes, and the rest data
0101   if (flag_for_inconsistent_events != 0) {
0102     m_numInconsistentPayloadEvents++;
0103   }
0104   //check for consistent events with the counter
0105   if (m_doCounter && flag_for_consistent_events != 0) {
0106     m_numConsistentEvents++;
0107   }
0108 }
0109 
0110 bool GctFibreAnalyzer::CheckForBC0(const L1GctFibreWord fibre) {
0111   // Check for BC0 on this event
0112   if (fibre.data() & 0x8000) {
0113     return true;
0114   } else {
0115     return false;
0116   }
0117 }
0118 
0119 bool GctFibreAnalyzer::CheckFibreWord(const L1GctFibreWord fibre) {
0120   // Check that the phase or clock bit (MSB bit on cycle 1) is set as it should be
0121   if (fibre.data() & 0x80000000) {
0122     return true;
0123   } else {
0124     return false;
0125   }
0126 }
0127 
0128 void GctFibreAnalyzer::CheckCounter(const L1GctFibreWord fibre) {
0129   // Remove MSB from both cycles
0130   int cycle0Data, cycle1Data;
0131 
0132   cycle0Data = fibre.data() & 0x7FFF;
0133   cycle1Data = (fibre.data() >> 16) & 0x7FFF;
0134 
0135   // Check to see if fibre numbers are consistent
0136   if ((cycle0Data + 1) != cycle1Data) {
0137     edm::LogInfo("GCT fibre data error") << "Fibre data not incrementing in cycles 0 and 1 "
0138                                          << " Cycle 0 data=" << cycle0Data << " Cycle 1 data=" << cycle1Data << " "
0139                                          << fibre;
0140   }
0141 
0142   // For now just write out the data
0143   edm::LogInfo("GCT fibre counter data") << " Fibre data: cycle0=" << cycle0Data << " cycle1=" << cycle1Data << " "
0144                                          << fibre;
0145 }
0146 
0147 void GctFibreAnalyzer::CheckLogicalID(const L1GctFibreWord fibre) {
0148   //added by Jad Marrouche, May 08
0149 
0150   unsigned ref_jf_link[] = {1, 2, 3, 4, 1, 2, 3, 4};
0151   //this array lists the link number ordering we expect from the 3 JFs in positive eta
0152   //below, we modify indices 2 and 3 from 3,4 to 1,2 to represent negative eta
0153 
0154   unsigned ref_eta0_link[] = {3, 4, 3, 4, 3, 4};
0155   //this array lists the link number ordering we expect from the ETA0
0156 
0157   unsigned ref_jf_type[] = {2, 2, 3, 3, 1, 1, 1, 1};
0158   //this array lists the SC_type ordering we expect for the JFs
0159 
0160   unsigned ref_eta0_type[] = {2, 2, 2, 2, 2, 2};
0161   //this array lists the SC_type ordering we expect for the ETA0 (for consistency)
0162 
0163   int eta_region, rct_phi_region, leaf_phi_region, jf_type, elec_type, local_source_card_id, source_card_id_read,
0164       source_card_id_expected;
0165 
0166   // Check that data in cycle 0 and cycle 1 are equal
0167   if ((fibre.data() & 0x7FFF) != ((fibre.data() & 0x7FFF0000) >> 16)) {
0168     edm::LogInfo("GCT fibre data error") << "Fibre data different on cycles 0 and 1 " << fibre;
0169   }
0170 
0171   //fibre.block() gives 0x90c etc
0172   //fibre.index() runs from 0 to 6/8
0173 
0174   if ((fibre.block() >> 10) & 0x1) {
0175     eta_region = 0;      //negative eta region
0176     ref_jf_link[2] = 1;  //modify indices to represent neg_eta fibre mapping
0177     ref_jf_link[3] = 2;
0178   } else
0179     eta_region = 1;  //positive eta region
0180 
0181   if (((fibre.block() >> 8) & 0x7) == 0 || ((fibre.block() >> 8) & 0x7) == 4)  //i.e. electron leaf cards
0182   {
0183     if ((fibre.block() & 0xFF) == 0x04)
0184       elec_type = 1;
0185     else if ((fibre.block() & 0xFF) == 0x84)
0186       elec_type = 0;
0187     else
0188       throw cms::Exception("Unknown GCT fibre data block ") << fibre.block();  //else something screwed up
0189 
0190     rct_phi_region = (fibre.index() / 3) + (4 * elec_type);
0191 
0192     local_source_card_id = (4 * eta_region);
0193 
0194     source_card_id_expected = (8 * rct_phi_region) + local_source_card_id;
0195 
0196     source_card_id_read = (fibre.data() >> 8) & 0x7F;
0197 
0198     if (source_card_id_expected != source_card_id_read) {
0199       edm::LogInfo("GCT fibre data error")
0200           << "Electron Source Card IDs do not match "
0201           << "Expected ID = " << source_card_id_expected << " ID read from data = " << source_card_id_read << " "
0202           << fibre;  //screwed up
0203     }
0204 
0205     if ((fibre.data() & 0xFF) != (unsigned int)(2 + fibre.index() % 3)) {
0206       edm::LogInfo("GCT fibre data error")
0207           << "Electron Fibres do not match "
0208           << "Expected Fibre = " << (2 + fibre.index() % 3) << " Fibre read from data = " << (fibre.data() & 0xFF)
0209           << " " << fibre;  //screwed up
0210     }
0211 
0212   } else  //i.e. jet leaf cards
0213   {
0214     //the reason we use these values for eta_region is so it is easy to add 4 to the local source card ID
0215     //remember that 0x9.. 0xA.. and 0xB.. are +ve eta block headers
0216     //whereas 0xD.., 0xE.. and 0xF.. are -ve eta
0217     //can distinguish between them using the above mask and shift
0218 
0219     if ((fibre.block() & 0xFF) == 0x04)
0220       jf_type = 1;  //JF2
0221     else if ((fibre.block() & 0xFF) == 0x0C)
0222       jf_type = 2;  //JF3
0223     else if ((fibre.block() & 0xFF) == 0x84)
0224       jf_type = -1;  //ETA0
0225     else if ((fibre.block() & 0xFF) == 0x8C)
0226       jf_type = 0;  //JF1
0227     else
0228       throw cms::Exception("Unknown GCT fibre data block ") << fibre.block();  //else something screwed up
0229 
0230     //edm::LogInfo("JF Type Info") << "JF TYPE = " << jf_type << " block = " << fibre;  //INFO ONLY
0231 
0232     leaf_phi_region = ((fibre.block() >> 8) & 0x7) - 1;  //0,1,2,3,4,5 for leaf cards
0233     if (eta_region == 0)
0234       leaf_phi_region--;  //need to do this because block index goes 9.. A.. B.. D.. E.. F.. - 8 and C are reserved for electron leafs which are dealt with above
0235     if (leaf_phi_region < 0 || leaf_phi_region > 5)
0236       throw cms::Exception("Unknown Leaf Card ") << fibre.block();
0237     //throw exception if number is outside 0-5 which means something screwed up
0238     //if(leaf_phi_region <0 || leaf_phi_region >5) edm::LogInfo("GCT fibre data error") << "Unknown leaf card " << fibre;
0239     //write to logger if number is outside 0-5 which means something screwed up
0240 
0241     if (jf_type == -1) {
0242       //in this case fibre.index() runs from 0-5
0243       //JF1 comes first, followed by JF2 and JF3
0244 
0245       if (fibre.index() <= 5)  //the compiler warning is because fibre.index() is unsigned int and hence is always >=0
0246       {
0247         rct_phi_region = ((8 + ((leaf_phi_region % 3) * 3) + (fibre.index() / 2)) % 9);
0248         //fibre.index()/2 will give 0 for 0,1 1 for 2,3 and 2 for 4,5
0249 
0250         //local_source_card_id = ref_eta0_type[ fibre.index() ] + (4 * (1 % eta_region));
0251         //take the ones complement of the eta_region because this is the shared part (i.e. other eta0 region)
0252         //this is done by (1 % eta_region) since 1%0 = 1 and 1%1=0
0253         //FLAWED - since 1%0 is a floating point exception you idiot!
0254 
0255         local_source_card_id = ref_eta0_type[fibre.index()] + (4 + (eta_region * -4));
0256         //this gives what you want - adds 4 when eta_region = 0 (neg) and adds 0 when eta_region = 1 (pos)
0257 
0258         source_card_id_expected = (8 * rct_phi_region) + local_source_card_id;
0259         //from GCT_refdoc_v2_2.pdf
0260 
0261         source_card_id_read = (fibre.data() >> 8) & 0x7F;
0262 
0263         if (source_card_id_expected != source_card_id_read) {
0264           edm::LogInfo("GCT fibre data error")
0265               << "ETA0 Source Card IDs do not match "
0266               << "Expected ID = " << source_card_id_expected << " ID read from data = " << source_card_id_read << " "
0267               << fibre;  //screwed up
0268         }
0269 
0270         if ((fibre.data() & 0xFF) != ref_eta0_link[fibre.index()]) {
0271           edm::LogInfo("GCT fibre data error")
0272               << "ETA0 Fibres do not match "
0273               << "Expected Fibre = " << ref_eta0_link[fibre.index()]
0274               << " Fibre read from data = " << (fibre.data() & 0xFF) << " " << fibre;  //screwed up
0275         }
0276       } else
0277         edm::LogInfo("GCT fibre data error") << "ETA0 Fibre index out of bounds " << fibre;
0278       //edm::LogInfo("Fibre Index Info") << "ETA0 Fibre index = " << fibre.index();
0279     }
0280 
0281     if (jf_type >= 0) {
0282       if (fibre.index() <= 7) {
0283         rct_phi_region = ((8 + ((leaf_phi_region % 3) * 3) + jf_type) % 9);  //see table below
0284 
0285         /*
0286         Leaf Card   |   RCT crate   |   Jet Finder
0287         ___________________________________________
0288         LC3 |   LC0 |   17  |   8   |   JF1
0289                 |       |   9   |   0   |   JF2
0290                 |       |   10  |   1   |   JF3
0291         ___________________________________________
0292         LC4 |   LC1 |   11  |   2   |   JF1
0293                 |       |   12  |   3   |   JF2
0294                 |       |   13  |   4   |   JF3
0295         ___________________________________________
0296         LC5 |   LC2 |   14  |   5   |   JF1
0297                 |       |   15  |   6   |   JF2
0298                 |       |   16  |   7   |   JF3
0299         ___________________________________________
0300         The phase results in the 17/8 being at the top
0301         This can be adjusted as necessary by changing
0302         the number 8 added before modulo 9 operation
0303               */
0304 
0305         local_source_card_id = ref_jf_type[fibre.index()] + (4 * eta_region);
0306 
0307         //since the SC sharing scheme renumbers SC 7 as SC3:
0308         if (local_source_card_id == 7)
0309           local_source_card_id = 3;
0310         //there is probably a more elegant way to do this
0311 
0312         source_card_id_expected = (8 * rct_phi_region) + local_source_card_id;
0313 
0314         source_card_id_read = (fibre.data() >> 8) & 0x7F;
0315 
0316         if (source_card_id_expected != source_card_id_read) {
0317           edm::LogInfo("GCT fibre data error")
0318               << "Source Card IDs do not match "
0319               << "Expected ID = " << source_card_id_expected << " ID read from data = " << source_card_id_read << " "
0320               << fibre;  //screwed up
0321         }
0322 
0323         if ((fibre.data() & 0xFF) != ref_jf_link[fibre.index()]) {
0324           edm::LogInfo("GCT fibre data error")
0325               << "Fibres do not match "
0326               << "Expected Fibre = " << ref_jf_link[fibre.index()]
0327               << " Fibre read from data = " << (fibre.data() & 0xFF) << " " << fibre;  //screwed up
0328         }
0329       } else
0330         edm::LogInfo("GCT fibre data error") << "Fibre index out of bounds " << fibre;
0331     }
0332   }
0333 }