CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LaserAlignmentEventFilter.cc
Go to the documentation of this file.
2 
10 
14 
15 #include <algorithm>
16 
17 // make a sorted copy of a vector, optionally of a different type
18 template <typename T, typename U>
19 std::vector<T> copy_and_sort_vector(std::vector<U> const & input)
20 {
21  std::vector<T> copy(input.begin(), input.end());
22  std::sort(copy.begin(), copy.end());
23  return copy;
24 }
25 
30  FED_collection_token( consumes<FEDRawDataCollection>(iConfig.getParameter<edm::InputTag>("FedInputTag")) ),
31  las_fed_ids( copy_and_sort_vector<uint16_t>(iConfig.getParameter<std::vector<int>>("FED_IDs")) ),
32  las_signal_ids( copy_and_sort_vector<uint32_t>(iConfig.getParameter<std::vector<int>>("SIGNAL_IDs")) ),
33  single_channel_thresh(iConfig.getParameter<unsigned>("SINGLE_CHANNEL_THRESH") ),
34  channel_count_thresh( iConfig.getParameter<unsigned>("CHANNEL_COUNT_THRESH") )
35 {
36 }
37 
42 }
43 
47 {
48  //unsigned int det_ctr = 0; // count how many modules are tested for signal
49  unsigned int sig_ctr = 0; // count how many modules have signal
50  //unsigned long buffer_sum = 0; // sum of buffer sizes
51 
52  // retrieve FED raw data (by label, which is "source" by default)
54  iEvent.getByToken( FED_collection_token, buffers );
55 
56  // read the cabling map from the EventSetup
58  iSetup.get<SiStripFedCablingRcd>().get( cabling );
59 
60  std::vector<uint16_t>::const_iterator ifed = las_fed_ids.begin();
61  for ( ; ifed != las_fed_ids.end(); ifed++ ) {
62  // Retrieve FED raw data for given FED
63  const FEDRawData& input = buffers->FEDData( static_cast<int>(*ifed) );
64  LogDebug("LaserAlignmentEventFilter") << "Examining FED " << *ifed;
65 
66  // check on FEDRawData pointer
67  if (not input.data())
68  continue;
69 
70  // check on FEDRawData size
71  if (not input.size())
72  continue;
73 
74  // construct FEDBuffer
75  sistrip::FEDBuffer buffer(input.data(), input.size());
76  if (not buffer.doChecks()) {
77  edm::LogWarning("LaserAlignmentEventFilter") << "FED Buffer check fails for FED ID " << *ifed << ".";
78  continue;
79  }
80 
81  // get the cabling connections for this FED
82  auto const & conns = cabling->fedConnections(*ifed);
83 
84  // Iterate through FED channels, extract payload and create Digis
85  std::vector<FedChannelConnection>::const_iterator iconn = conns.begin();
86  for ( ; iconn != conns.end(); iconn++ ) {
87  if (std::binary_search(las_signal_ids.begin(), las_signal_ids.end(), iconn->detId())) {
88  LogDebug("LaserAlignmentEventFilter")
89  << " Found LAS signal module in FED "
90  << *ifed
91  << " DetId: "
92  << iconn->detId() << "\n"
93  << "buffer.channel(iconn->fedCh()).size(): " << buffer.channel(iconn->fedCh()).length();
94  //++det_ctr;
95  if (buffer.channel(iconn->fedCh()).length() > single_channel_thresh)
96  ++sig_ctr;
97  //buffer_sum += buffer.channel(iconn->fedCh()).length();
98  if (sig_ctr > channel_count_thresh){
99  LogDebug("LaserAlignmentEventFilter") << "Event identified as LAS";
100  return true;
101  }
102  }
103  } // channel loop
104  } // FED loop
105 
106 // LogDebug("LaserAlignmentEventFilter") << det_ctr << " channels were tested for signal\n"
107 // <<sig_ctr << " channels have signal\n"
108 // << "Sum of buffer sizes: " << buffer_sum;
109 
110 
111  return false;
112 }
113 
114 //define this as a plug-in
#define LogDebug(id)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const std::vector< uint16_t > las_fed_ids
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
static std::string const input
Definition: EdmProvDump.cc:43
int iEvent
Definition: GenABIO.cc:230
std::vector< T > copy_and_sort_vector(std::vector< U > const &input)
const edm::EDGetTokenT< FEDRawDataCollection > FED_collection_token
const std::vector< uint32_t > las_signal_ids
virtual bool filter(edm::StreamID, edm::Event &, edm::EventSetup const &) const override
LaserAlignmentEventFilter(const edm::ParameterSet &)
const T & get() const
Definition: EventSetup.h:56
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28