CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/L1Trigger/L1GctAnalyzer/src/GctFibreAnalyzer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    GctFibreAnalyzer
00004 // Class:      GctFibreAnalyzer
00005 // 
00011 //
00012 // Original Author:  Alex Tapper
00013 //         Created:  Thu Jul 12 14:21:06 CEST 2007
00014 // $Id: GctFibreAnalyzer.cc,v 1.17 2011/01/19 07:32:18 elmer Exp $
00015 //
00016 //
00017 
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h" // Logger
00019 #include "FWCore/Utilities/interface/Exception.h" // Exceptions
00020 
00021 // Include file
00022 #include "L1Trigger/L1GctAnalyzer/interface/GctFibreAnalyzer.h"
00023 
00024 // Data format
00025 #include "DataFormats/L1GlobalCaloTrigger/interface/L1GctCollections.h"
00026 
00027 GctFibreAnalyzer::GctFibreAnalyzer(const edm::ParameterSet& iConfig):
00028   m_fibreSource(iConfig.getUntrackedParameter<edm::InputTag>("FibreSource")),
00029   m_doLogicalID(iConfig.getUntrackedParameter<bool>("doLogicalID")),
00030   m_doCounter(iConfig.getUntrackedParameter<bool>("doCounter")),
00031   m_numZeroEvents(0),
00032   m_numInconsistentPayloadEvents(0),
00033   m_numConsistentEvents(0)
00034 {
00035 }
00036 
00037 GctFibreAnalyzer::~GctFibreAnalyzer()
00038 {
00039   edm::LogInfo("Zero Fibreword events") << "Total number of events with zero fibrewords: " << m_numZeroEvents;
00040   edm::LogInfo("Inconsistent Payload events") << "Total number of events with inconsistent payloads: " << m_numInconsistentPayloadEvents;
00041   if(m_doCounter){edm::LogInfo("Successful events") << "Total number of Successful events: " << m_numConsistentEvents;}
00042 }
00043 
00044 void GctFibreAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00045 {
00046   using namespace edm;
00047   
00048   Handle<L1GctFibreCollection> fibre;
00049   iEvent.getByLabel(m_fibreSource,fibre);
00050 
00051   bool bc0= false;
00052   int flag_for_zeroes = 0;
00053   int flag_for_inconsistent_events = 0;
00054   unsigned int flag_for_consistency = 0;
00055   int flag_for_consistent_events = 0;
00056 
00057   for (L1GctFibreCollection::const_iterator f=fibre->begin(); f!=fibre->end(); f++){
00058         
00059     if (f->data()!=0)
00060       {
00061         if(m_doCounter) 
00062           {
00063             if(f==fibre->begin()) {flag_for_consistency = f->data();}
00064             else if(flag_for_consistency == f->data()){flag_for_consistent_events++;}
00065           }
00066 
00067         // Check for corrupt fibre data
00068         if (!CheckFibreWord(*f)){
00069           edm::LogInfo("GCT fibre data error") << "Missing phase bit (clock) in fibre data " << (*f);
00070         }
00071     
00072         // Check for BC0
00073         if (CheckForBC0(*f) && (f==fibre->begin())) {
00074           bc0=true;
00075         }
00076 
00077         // Check for mismatch between fibres
00078         if ((bc0 && !CheckForBC0(*f)) ||
00079             (!bc0 && CheckForBC0(*f))){
00080           edm::LogInfo("GCT fibre data error") << "BC0 mismatch in fibre data " << (*f);
00081         }
00082 
00083         // Check logical ID pattern
00084         if (m_doLogicalID) CheckLogicalID(*f);
00085 
00086         // Check counter pattern
00087         if (m_doCounter) CheckCounter(*f);
00088 
00089         flag_for_zeroes = 1;    
00090       }
00091     else 
00092       {
00093         //this basically checks that all the data is 0s by the time it gets to the last iteration
00094         if(flag_for_zeroes == 0 && f==(fibre->end()-1)) {m_numZeroEvents++;}
00095         //if the zero flag is set to 1 and it managed to find its way into here then something is wrong!
00096         if(flag_for_zeroes == 1) {flag_for_inconsistent_events++;}
00097       }
00098 
00099   }
00100   //check for inconsistent events i.e. those with one(or more) zeroes, and the rest data
00101   if(flag_for_inconsistent_events != 0) {m_numInconsistentPayloadEvents++;}
00102   //check for consistent events with the counter
00103   if(m_doCounter && flag_for_consistent_events != 0) {m_numConsistentEvents++;}
00104 }
00105 
00106 bool GctFibreAnalyzer::CheckForBC0(const L1GctFibreWord fibre)
00107 {
00108   // Check for BC0 on this event
00109   if (fibre.data() & 0x8000){
00110     return true;
00111   } else {
00112     return false;
00113   }
00114 }
00115 
00116 bool GctFibreAnalyzer::CheckFibreWord(const L1GctFibreWord fibre)
00117 {
00118   // Check that the phase or clock bit (MSB bit on cycle 1) is set as it should be
00119   if (fibre.data() & 0x80000000){
00120     return true;
00121   } else {
00122     return false;
00123   }
00124 }
00125 
00126 void GctFibreAnalyzer::CheckCounter(const L1GctFibreWord fibre)
00127 {
00128   // Remove MSB from both cycles
00129   int cycle0Data, cycle1Data;
00130   
00131   cycle0Data = fibre.data() & 0x7FFF;
00132   cycle1Data = (fibre.data() >> 16) & 0x7FFF;
00133 
00134   // Check to see if fibre numbers are consistent
00135   if ((cycle0Data+1)!=cycle1Data){
00136     edm::LogInfo("GCT fibre data error") << "Fibre data not incrementing in cycles 0 and 1 "
00137                                          << " Cycle 0 data=" << cycle0Data
00138                                          << " Cycle 1 data=" << cycle1Data
00139                                          << " " << fibre;      
00140   }  
00141   
00142   // For now just write out the data
00143   edm::LogInfo("GCT fibre counter data") << " Fibre data: cycle0=" << cycle0Data 
00144                                          << " cycle1=" << cycle1Data
00145                                          << " " << fibre;
00146 }
00147 
00148 void GctFibreAnalyzer::CheckLogicalID(const L1GctFibreWord fibre)
00149 {
00150   //added by Jad Marrouche, May 08
00151 
00152   unsigned ref_jf_link[] = {1,2,3,4,1,2,3,4};
00153   //this array lists the link number ordering we expect from the 3 JFs in positive eta
00154   //below, we modify indices 2 and 3 from 3,4 to 1,2 to represent negative eta
00155 
00156   unsigned ref_eta0_link[] = {3,4,3,4,3,4};
00157   //this array lists the link number ordering we expect from the ETA0
00158 
00159   unsigned ref_jf_type[] = {2,2,3,3,1,1,1,1};
00160   //this array lists the SC_type ordering we expect for the JFs
00161 
00162   unsigned ref_eta0_type[] = {2,2,2,2,2,2};
00163   //this array lists the SC_type ordering we expect for the ETA0 (for consistency)
00164 
00165   int eta_region, rct_phi_region, leaf_phi_region, jf_type, elec_type, local_source_card_id, source_card_id_read, source_card_id_expected;
00166 
00167   // Check that data in cycle 0 and cycle 1 are equal
00168   if ((fibre.data()&0x7FFF)!=((fibre.data()&0x7FFF0000)>>16)){
00169     edm::LogInfo("GCT fibre data error") << "Fibre data different on cycles 0 and 1 " << fibre;
00170   }
00171 
00172   //fibre.block() gives 0x90c etc
00173   //fibre.index() runs from 0 to 6/8
00174 
00175   if((fibre.block() >> 10) & 0x1 ) 
00176     {
00177       eta_region = 0;           //negative eta region           
00178       ref_jf_link[2] = 1; //modify indices to represent neg_eta fibre mapping
00179       ref_jf_link[3] = 2;       
00180     }   
00181   else eta_region = 1;  //positive eta region 
00182 
00183   if(((fibre.block() >> 8) & 0x7) == 0 || ((fibre.block() >> 8) & 0x7) == 4)    //i.e. electron leaf cards
00184     {
00185         
00186       if((fibre.block() & 0xFF)==0x04)          elec_type=1;
00187       else if((fibre.block() & 0xFF)==0x84)     elec_type=0;
00188       else throw cms::Exception("Unknown GCT fibre data block ") << fibre.block(); //else something screwed up   
00189         
00190       rct_phi_region = (fibre.index() / 3) + (4*elec_type);
00191 
00192       local_source_card_id = (4*eta_region);
00193         
00194       source_card_id_expected = (8 * rct_phi_region) + local_source_card_id;
00195         
00196       source_card_id_read = (fibre.data() >> 8) & 0x7F;
00197                 
00198       if(source_card_id_expected != source_card_id_read ) 
00199         {
00200           edm::LogInfo("GCT fibre data error") << "Electron Source Card IDs do not match "  
00201                                                << "Expected ID = " << source_card_id_expected
00202                                                << " ID read from data = " << source_card_id_read
00203                                                << " " << fibre; //screwed up
00204         }
00205 
00206       if( (fibre.data() & 0xFF) != (unsigned int)(2 + fibre.index()%3))
00207         {
00208           edm::LogInfo("GCT fibre data error") << "Electron Fibres do not match "  
00209                                                << "Expected Fibre = " << (2 + fibre.index()%3)
00210                                                << " Fibre read from data = " << (fibre.data() & 0xFF)
00211                                                << " " << fibre; //screwed up
00212         }
00213 
00214 
00215     }
00216   else  //i.e. jet leaf cards
00217     {
00218       //the reason we use these values for eta_region is so it is easy to add 4 to the local source card ID
00219       //remember that 0x9.. 0xA.. and 0xB.. are +ve eta block headers
00220       //whereas 0xD.., 0xE.. and 0xF.. are -ve eta
00221       //can distinguish between them using the above mask and shift
00222 
00223       if((fibre.block() & 0xFF)==0x04)          jf_type=1;      //JF2
00224       else if((fibre.block() & 0xFF)==0x0C)     jf_type=2;      //JF3
00225       else if((fibre.block() & 0xFF)==0x84)     jf_type=-1;     //ETA0
00226       else if((fibre.block() & 0xFF)==0x8C)     jf_type=0;      //JF1
00227       else throw cms::Exception("Unknown GCT fibre data block ") << fibre.block(); //else something screwed up   
00228 
00229       //edm::LogInfo("JF Type Info") << "JF TYPE = " << jf_type << " block = " << fibre;        //INFO ONLY
00230 
00231       leaf_phi_region = ((fibre.block() >> 8) & 0x7)-1;         //0,1,2,3,4,5 for leaf cards
00232       if(eta_region == 0) 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
00233       if(leaf_phi_region <0 || leaf_phi_region >5) throw cms::Exception("Unknown Leaf Card ") << fibre.block(); 
00234       //throw exception if number is outside 0-5 which means something screwed up
00235       //if(leaf_phi_region <0 || leaf_phi_region >5) edm::LogInfo("GCT fibre data error") << "Unknown leaf card " << fibre;
00236       //write to logger if number is outside 0-5 which means something screwed up
00237 
00238       if(jf_type == -1)
00239         {
00240           //in this case fibre.index() runs from 0-5
00241           //JF1 comes first, followed by JF2 and JF3
00242         
00243           if(fibre.index() <=5 )        //the compiler warning is because fibre.index() is unsigned int and hence is always >=0
00244             {
00245               rct_phi_region = ( (8 + ((leaf_phi_region%3)*3) + (fibre.index() / 2) ) % 9);
00246               //fibre.index()/2 will give 0 for 0,1 1 for 2,3 and 2 for 4,5
00247                 
00248               //local_source_card_id = ref_eta0_type[ fibre.index() ] + (4 * (1 % eta_region));
00249               //take the ones complement of the eta_region because this is the shared part (i.e. other eta0 region)
00250               //this is done by (1 % eta_region) since 1%0 = 1 and 1%1=0
00251                           //FLAWED - since 1%0 is a floating point exception you idiot!
00252 
00253                           local_source_card_id = ref_eta0_type[ fibre.index() ] + (4 + (eta_region * -4));
00254                           //this gives what you want - adds 4 when eta_region = 0 (neg) and adds 0 when eta_region = 1 (pos)
00255         
00256               source_card_id_expected = (8 * rct_phi_region) + local_source_card_id;
00257               //from GCT_refdoc_v2_2.pdf
00258                 
00259               source_card_id_read = (fibre.data() >> 8) & 0x7F;
00260                 
00261               if(source_card_id_expected != source_card_id_read ) 
00262                 {
00263                   edm::LogInfo("GCT fibre data error") << "ETA0 Source Card IDs do not match "  
00264                                                        << "Expected ID = " << source_card_id_expected
00265                                                        << " ID read from data = " << source_card_id_read
00266                                                        << " " << fibre; //screwed up
00267                 }
00268 
00269               if( (fibre.data() & 0xFF) != ref_eta0_link[fibre.index()])
00270                 {
00271                   edm::LogInfo("GCT fibre data error") << "ETA0 Fibres do not match "  
00272                                                        << "Expected Fibre = " << ref_eta0_link[fibre.index()]
00273                                                        << " Fibre read from data = " << (fibre.data() & 0xFF)
00274                                                        << " " << fibre; //screwed up
00275                 }
00276             }
00277           else edm::LogInfo("GCT fibre data error") << "ETA0 Fibre index out of bounds " << fibre;
00278           //edm::LogInfo("Fibre Index Info") << "ETA0 Fibre index = " << fibre.index();
00279         }
00280 
00281 
00282       if(jf_type >=0) 
00283         {
00284           if(fibre.index() <=7 )
00285             {
00286               rct_phi_region = ( (8 + ((leaf_phi_region%3)*3) + jf_type ) % 9);         //see table below
00287 
00288               /*
00289                 Leaf Card       |       RCT crate       |       Jet Finder
00290                 ___________________________________________
00291                 LC3     |       LC0     |       17      |       8       |       JF1
00292                 |               |       9       |       0       |       JF2
00293                 |               |       10      |       1       |       JF3
00294                 ___________________________________________
00295                 LC4     |       LC1     |       11      |       2       |       JF1
00296                 |               |       12      |       3       |       JF2
00297                 |               |       13      |       4       |       JF3
00298                 ___________________________________________
00299                 LC5     |       LC2     |       14      |       5       |       JF1
00300                 |               |       15      |       6       |       JF2
00301                 |               |       16      |       7       |       JF3
00302                 ___________________________________________
00303                 The phase results in the 17/8 being at the top
00304                 This can be adjusted as necessary by changing
00305                 the number 8 added before modulo 9 operation
00306               */
00307 
00308               local_source_card_id = ref_jf_type[ fibre.index() ] + (4 * eta_region);
00309 
00310               //since the SC sharing scheme renumbers SC 7 as SC3:
00311               if(local_source_card_id == 7) local_source_card_id = 3;
00312               //there is probably a more elegant way to do this
00313 
00314               source_card_id_expected = (8 * rct_phi_region) + local_source_card_id;
00315 
00316               source_card_id_read = (fibre.data() >> 8) & 0x7F;
00317 
00318               if(source_card_id_expected != source_card_id_read ) 
00319                 {
00320                   edm::LogInfo("GCT fibre data error") << "Source Card IDs do not match "  
00321                                                        << "Expected ID = " << source_card_id_expected
00322                                                        << " ID read from data = " << source_card_id_read
00323                                                        << " " << fibre; //screwed up
00324                 }
00325 
00326               if( (fibre.data() & 0xFF) != ref_jf_link[fibre.index()])
00327                 {
00328                   edm::LogInfo("GCT fibre data error") << "Fibres do not match "  
00329                                                        << "Expected Fibre = " << ref_jf_link[fibre.index()]
00330                                                        << " Fibre read from data = " << (fibre.data() & 0xFF)
00331                                                        << " " << fibre; //screwed up
00332                 }
00333             }
00334           else edm::LogInfo("GCT fibre data error") << "Fibre index out of bounds " << fibre;
00335         }
00336 
00337     }
00338 
00339 
00340 }