CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Alignment/LaserAlignment/plugins/LaserAlignmentEventFilter.cc

Go to the documentation of this file.
00001 #include "LaserAlignmentEventFilter.h"
00002 
00003 #include <FWCore/Framework/interface/Event.h> 
00004 #include <FWCore/Framework/interface/EventSetup.h> 
00005 #include <FWCore/Framework/interface/MakerMacros.h>
00006 #include <FWCore/Framework/interface/EventSetupRecordKey.h>
00007 #include <FWCore/Framework/interface/ESHandle.h>
00008 
00009 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00010 #include "CondFormats/SiStripObjects/interface/SiStripFedCabling.h"
00011 #include "CondFormats/DataRecord/interface/SiStripFedCablingRcd.h"
00012 #include "EventFilter/SiStripRawToDigi/interface/SiStripFEDBuffer.h"
00013 
00014 #include <algorithm>
00015 
00019 LaserAlignmentEventFilter::LaserAlignmentEventFilter( const edm::ParameterSet& iConfig ) :
00020   FED_collection(iConfig.getParameter<edm::InputTag>("FedInputTag")),
00021   signal_filter(true),
00022   single_channel_thresh(11),
00023   channel_count_thresh(4),
00024   LAS_event_count(0),
00025   cabling(0),
00026   cacheId_(0)
00027 {
00028   // Read in Filter Lists
00029   std::vector<int> FED_IDs = iConfig.getParameter<std::vector<int> >("FED_IDs");
00030   set_las_fed_ids(FED_IDs);
00031 
00032   std::vector<int> SIGNAL_IDs = iConfig.getParameter<std::vector<int> >("SIGNAL_IDs");
00033   set_las_signal_ids(SIGNAL_IDs);
00034 
00035   // Read in Filter Flags
00036   signal_filter = iConfig.getParameter<bool>("SIGNAL_Filter");
00037 
00038   // Read in Filter Parameters
00039   single_channel_thresh =  iConfig.getParameter<unsigned>("SINGLE_CHANNEL_THRESH");
00040   channel_count_thresh = iConfig.getParameter<unsigned>("CHANNEL_COUNT_THRESH");
00041 
00042   edm::LogInfo("LasFilterConstructor") << "\n" 
00043                                        << "\nSIGNAL_Filter: " << (signal_filter ? "true" : "false")
00044                                        << ", SIGNLE_CHANNEL_THRESH: " << single_channel_thresh 
00045                                        << ", CHANNEL_COUNT_THRESH: " << channel_count_thresh;
00046 }
00047 
00051 LaserAlignmentEventFilter::~LaserAlignmentEventFilter() {
00052 }
00053 
00057 void LaserAlignmentEventFilter::beginRun( const edm::EventSetup& iSetup) {
00058   updateCabling( iSetup );
00059 }
00060 
00064 
00065 bool LaserAlignmentEventFilter::filter( edm::Event& iEvent, const edm::EventSetup& iSetup ) 
00066 {
00067   updateCabling( iSetup );
00068   unsigned int det_ctr=0; // Count how many modules are tested for signal
00069   unsigned int sig_ctr=0; // Count how many modules have signal
00070   unsigned long buffer_sum=0; // Sum of buffer sizes
00071 
00072   // Retrieve FED raw data (by label, which is "source" by default)
00073   edm::Handle<FEDRawDataCollection> buffers;
00074   iEvent.getByLabel( FED_collection, buffers ); 
00075 
00076 
00077   std::vector<uint16_t>::const_iterator ifed = las_fed_ids.begin();
00078   for ( ; ifed != las_fed_ids.end(); ifed++ ) {
00079     // Retrieve FED raw data for given FED 
00080     const FEDRawData& input = buffers->FEDData( static_cast<int>(*ifed) );
00081     LogDebug("LaserAlignmentEventFilter") << "Examining FED " << *ifed; 
00082 
00083      // Check on FEDRawData pointer
00084      if ( !input.data() ) {
00085        continue;
00086      }  
00087     
00088      // Check on FEDRawData size
00089      if ( !input.size() ) {
00090        continue;
00091      }
00092 
00093      // get the cabling connections for this FED
00094      const std::vector<FedChannelConnection>& conns = cabling->connections(*ifed);
00095 
00096      // construct FEDBuffer
00097      std::auto_ptr<sistrip::FEDBuffer> buffer;
00098      try {
00099        buffer.reset(new sistrip::FEDBuffer(input.data(),input.size()));
00100        if (!buffer->doChecks()) {
00101         throw cms::Exception("FEDBuffer") << "FED Buffer check fails for FED ID " << *ifed << ". (BW)";
00102         if(buffer->daqEventType() != sistrip::DAQ_EVENT_TYPE_CALIBRATION){
00103           edm::LogWarning("LaserAlignmentEventFilter") << "Event is not calibration type";
00104           return false;
00105         }
00106        }
00107      }
00108      catch (const cms::Exception& e) { 
00109        if ( edm::isDebugEnabled() ) {
00110         edm::LogWarning("LaserAlignmentEventFilter") << "Exception caught when creating FEDBuffer object for FED " << *ifed << ": " << e.what();
00111        }
00112        continue;
00113      }
00114 
00115     // Iterate through FED channels, extract payload and create Digis
00116     std::vector<FedChannelConnection>::const_iterator iconn = conns.begin();
00117     for ( ; iconn != conns.end(); iconn++ ) {
00118 
00119       if(signal_filter){
00120         if ( std::binary_search(las_signal_ids.begin(), las_signal_ids.end(), iconn->detId())){ 
00121           LogDebug("LaserAlignmentEventFilter")
00122             << " Found LAS signal module in FED " 
00123             << *ifed
00124             << "  DetId: " 
00125             << iconn->detId() << "\n"
00126             << "buffer->channel(iconn->fedCh()).size(): " << buffer->channel(iconn->fedCh()).length();
00127           buffer_size.push_back(buffer->channel(iconn->fedCh()).length());
00128           det_ctr ++;
00129           if(buffer->channel(iconn->fedCh()).length() > single_channel_thresh) sig_ctr++;
00130           buffer_sum += buffer->channel(iconn->fedCh()).length();
00131           if(sig_ctr > channel_count_thresh){
00132             LAS_event_count ++;
00133             LogDebug("LaserAlignmentEventFilter") << "Event identified as LAS";
00134             return true;
00135           }
00136         }
00137       }
00138     } // channel loop
00139   } // FED loop
00140   
00141 //   LogDebug("LaserAlignmentEventFilter") << det_ctr << " channels were tested for signal\n" 
00142 //                          <<sig_ctr << " channels have signal\n"
00143 //                          << "Sum of buffer sizes: " << buffer_sum;
00144 
00145 
00146   return false;
00147 }
00148 
00149 
00153 void LaserAlignmentEventFilter::endJob() {
00154   //edm::LogInfo("LaserAlignmentEventFilter") << "found " << LAS_event_count << " LAS events";
00155 }
00156 
00157 // Create the table of FEDs that contain LAS Modules
00158 void LaserAlignmentEventFilter::set_las_fed_ids(const std::vector<int>& las_feds)
00159 {
00160   // Convert the std::vector to a std::set
00161   las_fed_ids = std::vector<uint16_t>(las_feds.begin(), las_feds.end());
00162   // Sort for binary search
00163   std::sort(las_fed_ids.begin(), las_fed_ids.end());
00164 }
00165 
00166 // Create table of modules to be tested for signals
00167 void LaserAlignmentEventFilter::set_las_signal_ids(const std::vector<int>& las_signal)
00168 {
00169   // Convert the std::vector to a std::set
00170   las_signal_ids = std::vector<uint32_t>(las_signal.begin(), las_signal.end());
00171   // Sort for binary search
00172   std::sort(las_signal_ids.begin(), las_signal_ids.end());
00173 }
00174 
00175 
00176 // This method was copied from EventFilter/SiStripRawToDigi/plugin/SiStripRawToDigiModule
00177 void LaserAlignmentEventFilter::updateCabling( const edm::EventSetup& setup ) {
00178 
00179   uint32_t cache_id = setup.get<SiStripFedCablingRcd>().cacheIdentifier();
00180 
00181   if ( cacheId_ != cache_id ) {
00182    edm::ESHandle<SiStripFedCabling> c;
00183    setup.get<SiStripFedCablingRcd>().get( c );
00184    cabling = c.product();
00185    cacheId_ = cache_id;
00186   }
00187 }
00188 
00189 
00190 //define this as a plug-in
00191 DEFINE_FWK_MODULE(LaserAlignmentEventFilter);