CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CSCDQM_EventProcessor.h
Go to the documentation of this file.
1 /*
2  * =====================================================================================
3  *
4  * Filename: EventProcessor.h
5  *
6  * Description: Object which processes Event and provides Hits to
7  * HitHandler object.
8  *
9  * Version: 1.0
10  * Created: 10/03/2008 10:26:04 AM
11  * Revision: none
12  * Compiler: gcc
13  *
14  * Author: Valdas Rapsevicius, valdas.rapsevicius@cern.ch
15  * Company: CERN, CH
16  *
17  * =====================================================================================
18  */
19 
20 #ifndef CSCDQM_EventProcessor_H
21 #define CSCDQM_EventProcessor_H
22 
23 #include <set>
24 #include <string>
25 #include <math.h>
26 
27 #include <TString.h>
28 
29 #ifdef DQMGLOBAL
30 
37 
38 #endif
39 
40 #include "CSCDQM_Logger.h"
41 #include "CSCDQM_Summary.h"
43 #include "CSCDQM_Configuration.h"
44 #include "CSCDQM_Configuration.h"
45 
50 
53 
54 namespace cscdqm {
55 
59  struct HWStandbyType {
60 
61  // if standby flag should be considered at all?
62  // at the start it will be false, thus good for last value ;)
63  bool process;
64 
65  // ME+
66  bool MeP;
67 
68  // ME-
69  bool MeM;
70 
72  process = MeP = MeM = false;
73  }
74 
75  void applyMeP(bool ready) {
76  MeP = MeP || !ready;
77  }
78 
79  void applyMeM(bool ready) {
80  MeM = MeM || !ready;
81  }
82 
83  bool fullStandby() const {
84  return (MeM && MeP);
85  }
86 
87  bool operator==(const HWStandbyType& t) const {
88  return (t.MeP == MeP && t.MeM == MeM && t.process == process);
89  }
90 
91  bool operator!=(const HWStandbyType& t) const {
92  return !(*this == t);
93  }
94 
96  MeP = t.MeP;
97  MeM = t.MeM;
98  process = t.process;
99  return *this;
100  }
101 
102  };
103 
104  typedef std::map<CSCIdType, ExaminerStatusType> CSCExaminerMapType;
105  typedef std::vector<DDUIdType> DDUExaminerVectorType;
106  // typedef std::map<int, long> CSCExaminerMapType;
107  // typedef std::vector<int> DDUExaminerVectorType;
108 
114 
115 // ===================================================================================================
116 // General stuff
117 // ===================================================================================================
118 
119  public:
120 
121  EventProcessor(Configuration* const p_config);
122 
123 #ifdef DQMGLOBAL
124 
125  EventProcessor(Configuration* const p_config, const edm::InputTag& itag, edm::ConsumesCollector& coco);
126 
127 #endif
128 
133 
134  void init();
135 
136  void updateFractionHistos();
137  void updateEfficiencyHistos();
138  void standbyEfficiencyHistos(HWStandbyType& standby);
139  void writeShifterHistograms();
140 
141  unsigned int maskHWElements(std::vector<std::string>& tokens);
142 
143  private:
144 
145  bool processExaminer(const CSCDCCExaminer& binChecker);
146  bool processExaminer(const CSCDCCExaminer& binChecker, const CSCDCCFormatStatusDigi& digi);
147  void processDDU(const CSCDDUEventData& data, const CSCDCCExaminer& binChecker);
148  void processCSC(const CSCEventData& data, const int dduID, const CSCDCCExaminer& binChecker);
149 
150  void calcEMUFractionHisto(const HistoId& result, const HistoId& set, const HistoId& subset);
151 
152  const bool getEMUHisto(const HistoId& histo, MonitorObject*& me);
153  const bool getFEDHisto(const HistoId& histo, const HwId& fedID, MonitorObject*& me);
154  const bool getDDUHisto(const HistoId& histo, const HwId& dduID, MonitorObject*& me);
155  const bool getCSCHisto(const HistoId& histo, const HwId& crateID, const HwId& dmbSlot, MonitorObject*& me);
156  const bool getCSCHisto(const HistoId& histo, const HwId& crateID, const HwId& dmbSlot, const HwId& adId, MonitorObject*& me);
157  const bool getParHisto(const HistoId& histo, MonitorObject*& me);
158  void preProcessEvent();
159 
160  const bool getCSCFromMap(const unsigned int& crateId, const unsigned int& dmbId, unsigned int& cscType, unsigned int& cscPosition) const;
161  void setEmuEventDisplayBit(MonitorObject*& mo, const unsigned int x, const unsigned int y, const unsigned int bit);
162  void resetEmuEventDisplays();
163 
166 
169 
170  std::map<uint32_t, uint32_t> L1ANumbers;
171  std::map<uint32_t, bool> fNotFirstEvent;
172  uint32_t L1ANumber;
173  uint32_t BXN;
174  uint32_t cntDMBs;
175  uint32_t cntCFEBs;
176  uint32_t cntALCTs;
177  uint32_t cntTMBs;
178 
179 
180  uint16_t theFormatVersion;
181 
182 
183  // bool fFirstEvent;
184  bool fCloseL1As; // Close L1A bit from DDU Trailer
186 
187  // token for input
189 
190 // ===================================================================================================
191 // Local ONLY stuff
192 // ===================================================================================================
193 
194 #ifdef DQMLOCAL
195 
196  public:
197 
198  void processEvent(const char* data, const int32_t dataSize, const uint32_t errorStat, const int32_t nodeNumber);
199 
200 #endif
201 
202 // ===================================================================================================
203 // Global ONLY stuff
204 // ===================================================================================================
205 
206 #ifdef DQMGLOBAL
207 
208  private:
209 
210  bool bCSCEventCounted;
211 
212  public:
213 
214  void processEvent(const edm::Event& e, const edm::InputTag& inputTag );
215 
216 #endif
217 
218  };
219 
220 }
221 
222 #endif
Object used to process Events and compute statistics.
const bool getCSCFromMap(const unsigned int &crateId, const unsigned int &dmbId, unsigned int &cscType, unsigned int &cscPosition) const
Get CSC type and position from crate and dmb identifiers.
uint32_t cntCFEBs
Total Number of DMBs per event from DDU Header DAV.
unsigned int HwId
bool operator!=(const HWStandbyType &t) const
void preProcessEvent()
Common Local and Global DQM function to be called before processing Event.
const bool getDDUHisto(const HistoId &histo, const HwId &dduID, MonitorObject *&me)
Get DDU Level Monitoring Object.
std::map< CSCIdType, ExaminerStatusType > CSCExaminerMapType
Monitoring Object interface used to cover Root object and provide common interface to EventProcessor ...
void updateEfficiencyHistos()
Update Efficiency MOs.
dispatcher processEvent(e, inputTag, standby)
Hardware and Physics Efficiency data structures and routines.
void processCSC(const CSCEventData &data, const int dduID, const CSCDCCExaminer &binChecker)
Process Chamber Data and fill MOs.
std::map< uint32_t, uint32_t > L1ANumbers
void init()
Initialize EventProcessor: reading out config information.
bool operator==(const HWStandbyType &t) const
void setEmuEventDisplayBit(MonitorObject *&mo, const unsigned int x, const unsigned int y, const unsigned int bit)
Set a single bit in the 3D Histogram (aka EMU level event display). Checks if mo and x != null...
const bool getCSCHisto(const HistoId &histo, const HwId &crateID, const HwId &dmbSlot, MonitorObject *&me)
Get CSC (Chamber) Level Monitoring Object.
uint32_t cntTMBs
Total Number of ALCTs per event from DMB DAV.
CSCDQM Framework Global Configuration.
unsigned int HistoId
EventProcessor(Configuration *const p_config)
Constructor.
uint32_t cntALCTs
Total Number of CFEBs per event from DMB DAV.
tuple result
Definition: query.py:137
const bool getParHisto(const HistoId &histo, MonitorObject *&me)
Get Parameter Monitoring Object.
void processDDU(const CSCDDUEventData &data, const CSCDCCExaminer &binChecker)
Process DDU output and fill MOs.
bool processExaminer(const CSCDCCExaminer &binChecker)
unsigned int maskHWElements(std::vector< std::string > &tokens)
Mask HW elements from the efficiency calculations. Can be applied on runtime!
const HWStandbyType & operator=(const HWStandbyType &t)
const bool getFEDHisto(const HistoId &histo, const HwId &fedID, MonitorObject *&me)
Get FED Level Monitoring Object.
CSC Format Status Object.
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
edm::EDGetTokenT< FEDRawDataCollection > frdtoken
void calcEMUFractionHisto(const HistoId &result, const HistoId &set, const HistoId &subset)
Calculate fractional histogram.
std::vector< DDUIdType > DDUExaminerVectorType
void resetEmuEventDisplays()
Reset Emu level EventDisplay histograms once per event.
bool fCloseL1As
Data Format version (2005, 2013)
uint16_t theFormatVersion
Total Number of TMBs per event from DMB DAV.
std::map< uint32_t, bool > fNotFirstEvent
void updateFractionHistos()
Update Fractional MOs.
void standbyEfficiencyHistos(HWStandbyType &standby)
apply standby flags/parameters
const bool getEMUHisto(const HistoId &histo, MonitorObject *&me)
Get EMU (Top Level) Monitoring Object.