CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTHcalCalibTypeFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HLTHcalCalibTypeFilter
4 // Class: HLTHcalCalibTypeFilter
5 //
13 //
14 // Original Author: Bryan DAHMES
15 // Created: Tue Jan 22 13:55:00 CET 2008
16 //
17 //
18 
19 
20 // system include files
21 #include <string>
22 #include <iostream>
23 #include <memory>
24 
25 // user include files
27 
35 
39 
41 
42 //
43 // constructors and destructor
44 //
46  DataInputToken_( consumes<FEDRawDataCollection>( config.getParameter<edm::InputTag>("InputTag") ) ),
47  CalibTypes_( config.getParameter< std::vector<int> >("CalibTypes") ),
48  Summary_( config.getUntrackedParameter<bool>("FilterSummary", false) ),
49  eventsByType_()
50 {
51  for (auto & i : eventsByType_) i = 0;
52 }
53 
54 
56 {
57 
58  // do anything here that needs to be done at desctruction time
59  // (e.g. close files, deallocate resources etc.)
60 
61 }
62 
63 void
66  desc.add<edm::InputTag>("InputTag",edm::InputTag("source"));
67  std::vector<int> temp; for (int i=1; i<=5; i++) temp.push_back(i);
68  desc.add<std::vector<int> >("CalibTypes", temp);
69  desc.addUntracked<bool>("FilterSummary",false);
70  descriptions.add("hltHcalCalibTypeFilter",desc);
71 }
72 
73 //
74 // member functions
75 //
76 
77 // ------------ method called on each new Event ------------
78 bool
80 {
81  using namespace edm;
82 
85 
86  // checking FEDs for calibration information
87  int calibType = -1 ; int numEmptyFEDs = 0 ;
88  std::vector<int> calibTypeCounter(8,0) ;
91  const FEDRawData& fedData = rawdata->FEDData(i) ;
92  if ( fedData.size() < 24 ) numEmptyFEDs++ ;
93  if ( fedData.size() < 24 ) continue ;
94  int value = ((const HcalDCCHeader*)(fedData.data()))->getCalibType() ;
95  calibTypeCounter.at(value)++ ; // increment the counter for this calib type
96  }
97  int maxCount = 0 ;
98  int numberOfFEDIds = FEDNumbering::MAXHCALFEDID - FEDNumbering::MINHCALFEDID + 1 ;
99  for (unsigned int i=0; i<calibTypeCounter.size(); i++) {
100  if ( calibTypeCounter.at(i) > maxCount ) { calibType = i ; maxCount = calibTypeCounter.at(i) ; }
101  if ( maxCount == numberOfFEDIds ) break ;
102  }
103  if ( calibType < 0 ) return false ; // No HCAL FEDs, thus no calibration type
104  if ( maxCount != (numberOfFEDIds-numEmptyFEDs) )
105  edm::LogWarning("HLTHcalCalibTypeFilter") << "Conflicting calibration types found. Assigning type "
106  << calibType ;
107  LogDebug("HLTHcalCalibTypeFilter") << "Calibration type is: " << calibType ;
108  if (Summary_)
109  ++eventsByType_.at(calibType);
110 
111  for (unsigned int i=0; i<CalibTypes_.size(); i++)
112  if ( calibType == CalibTypes_.at(i) ) return true ;
113  return false ;
114 }
115 
116 // ------------ method called once each job just after ending the event loop ------------
117 void
119  if ( Summary_ )
120  edm::LogWarning("HLTHcalCalibTypeFilter") << "Summary of filter decisions: "
121  << eventsByType_.at(hc_Null) << "(No Calib), "
122  << eventsByType_.at(hc_Pedestal) << "(Pedestal), "
123  << eventsByType_.at(hc_RADDAM) << "(RADDAM), "
124  << eventsByType_.at(hc_HBHEHPD) << "(HBHE/HPD), "
125  << eventsByType_.at(hc_HOHPD) << "(HO/HPD), "
126  << eventsByType_.at(hc_HFPMT) << "(HF/PMT)" ;
127 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
virtual void endJob(void) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
tuple calibType
Definition: diJetCalib.py:20
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const std::vector< int > CalibTypes_
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
const edm::EDGetTokenT< FEDRawDataCollection > DataInputToken_
int iEvent
Definition: GenABIO.cc:230
virtual bool filter(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::array< std::atomic< int >, 8 > eventsByType_
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
dictionary rawdata
Definition: lumiPlot.py:393
volatile std::atomic< bool > shutdown_flag false
HLTHcalCalibTypeFilter(const edm::ParameterSet &)