CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalCalibTypeFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HcalCalibTypeFilter
4 // Class: HcalCalibTypeFilter
5 //
13 //
14 // Original Author: Giovanni FRANZONI
15 // Created: Tue Jan 22 13:55:00 CET 2008
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
25 
31 
32 #include <string>
33 #include <iostream>
34 
40 
41 //
42 // class declaration
43 //
44 
46 public:
47  explicit HcalCalibTypeFilter(const edm::ParameterSet&);
48  virtual ~HcalCalibTypeFilter();
49 
50 private:
51  virtual void beginJob() override ;
52  virtual bool filter(edm::Event&, const edm::EventSetup&) override;
53  virtual void endJob() override ;
54 
55  // ----------member data ---------------------------
56 
58  bool Summary_ ;
59  std::vector<int> CalibTypes_ ;
60  std::vector<int> eventsByType ;
61 
62 };
63 
64 
65 //
66 // constructors and destructor
67 //
69 {
70  //now do what ever initialization is needed
71 
72  tok_data_ = consumes<FEDRawDataCollection>(edm::InputTag(iConfig.getParameter<std::string>("InputLabel") ));
73  Summary_ = iConfig.getUntrackedParameter<bool>("FilterSummary",false) ;
74  CalibTypes_ = iConfig.getParameter< std::vector<int> >("CalibTypes") ;
75 }
76 
77 
79 {
80 
81  // do anything here that needs to be done at desctruction time
82  // (e.g. close files, deallocate resources etc.)
83 
84 }
85 
86 
87 //
88 // member functions
89 //
90 
91 // ------------ method called on each new Event ------------
92 bool
94 {
95  using namespace edm;
96 
99 
100  if(!rawdata.isValid()){
101  return false;
102  }
103 
104 
105  // checking FEDs for calibration information
106  int calibType = -1 ; int numEmptyFEDs = 0 ;
107  std::vector<int> calibTypeCounter(8,0) ;
108  for (int i=FEDNumbering::MINHCALFEDID;
110  const FEDRawData& fedData = rawdata->FEDData(i) ;
111  if ( fedData.size() < 24 ) numEmptyFEDs++ ;
112  if ( fedData.size() < 24 ) continue ;
113  int value = ((const HcalDCCHeader*)(fedData.data()))->getCalibType() ;
114  calibTypeCounter.at(value)++ ; // increment the counter for this calib type
115  }
116 
117  int maxCount = 0 ;
118  int numberOfFEDIds = FEDNumbering::MAXHCALFEDID - FEDNumbering::MINHCALFEDID + 1 ;
119  for (unsigned int i=0; i<calibTypeCounter.size(); i++) {
120  if ( calibTypeCounter.at(i) > maxCount ) { calibType = i ; maxCount = calibTypeCounter.at(i) ; }
121  if ( maxCount == numberOfFEDIds ) break ;
122  }
123  if ( maxCount != (numberOfFEDIds-numEmptyFEDs) )
124  edm::LogWarning("HcalCalibTypeFilter") << "Conflicting calibration types found. Assigning type "
125  << calibType ;
126  LogDebug("HcalCalibTypeFilter") << "Calibration type is: " << calibType ;
127  eventsByType.at(calibType)++ ;
128  for (unsigned int i=0; i<CalibTypes_.size(); i++)
129  if ( calibType == CalibTypes_.at(i) ) return true ;
130  return false ;
131 }
132 
133 // ------------ method called once each job just before starting event loop ------------
134 void
136 {
137  eventsByType.clear() ;
138  eventsByType.resize(8,0) ;
139 }
140 
141 // ------------ method called once each job just after ending the event loop ------------
142 void
144  if ( Summary_ )
145  edm::LogWarning("HcalCalibTypeFilter") << "Summary of filter decisions: "
146  << eventsByType.at(hc_Null) << "(No Calib), "
147  << eventsByType.at(hc_Pedestal) << "(Pedestal), "
148  << eventsByType.at(hc_RADDAM) << "(RADDAM), "
149  << eventsByType.at(hc_HBHEHPD) << "(HBHE/HPD), "
150  << eventsByType.at(hc_HOHPD) << "(HO/HPD), "
151  << eventsByType.at(hc_HFPMT) << "(HF/PMT)" ;
152 }
153 
154 //define this as a plug-in
#define LogDebug(id)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
tuple calibType
Definition: diJetCalib.py:20
std::vector< int > CalibTypes_
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
int iEvent
Definition: GenABIO.cc:230
virtual bool filter(edm::Event &, const edm::EventSetup &) override
HcalCalibTypeFilter(const edm::ParameterSet &)
virtual void endJob() override
virtual void beginJob() override
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
edm::EDGetTokenT< FEDRawDataCollection > tok_data_
std::vector< int > eventsByType