CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 <cmath>
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  // if standby flag should be considered at all?
61  // at the start it will be false, thus good for last value ;)
62  bool process;
63 
64  // ME+
65  bool MeP;
66 
67  // ME-
68  bool MeM;
69 
70  HWStandbyType() { process = MeP = MeM = false; }
71 
72  void applyMeP(bool ready) { MeP = MeP || !ready; }
73 
74  void applyMeM(bool ready) { MeM = MeM || !ready; }
75 
76  bool fullStandby() const { return (MeM && MeP); }
77 
78  bool operator==(const HWStandbyType& t) const { return (t.MeP == MeP && t.MeM == MeM && t.process == process); }
79 
80  bool operator!=(const HWStandbyType& t) const { return !(*this == t); }
81 
83  MeP = t.MeP;
84  MeM = t.MeM;
85  process = t.process;
86  return *this;
87  }
88  };
89 
90  typedef std::map<CSCIdType, ExaminerStatusType> CSCExaminerMapType;
91  typedef std::vector<DDUIdType> DDUExaminerVectorType;
92  // typedef std::map<int, long> CSCExaminerMapType;
93  // typedef std::vector<int> DDUExaminerVectorType;
94 
100  // ===================================================================================================
101  // General stuff
102  // ===================================================================================================
103 
104  public:
105  EventProcessor(Configuration* const p_config);
106 
107 #ifdef DQMGLOBAL
108 
109  EventProcessor(Configuration* const p_config, const edm::InputTag& itag, edm::ConsumesCollector& coco);
110 
111 #endif
112 
117 
118  void init();
119 
120  void updateFractionHistos();
121  void updateEfficiencyHistos();
122  void standbyEfficiencyHistos(HWStandbyType& standby);
123  void writeShifterHistograms();
124 
125  unsigned int maskHWElements(std::vector<std::string>& tokens);
126 
127  private:
128  bool processExaminer(const CSCDCCExaminer& binChecker);
129  bool processExaminer(const CSCDCCExaminer& binChecker, const CSCDCCFormatStatusDigi& digi);
130  void processDDU(const CSCDDUEventData& data, const CSCDCCExaminer& binChecker);
131  void processCSC(const CSCEventData& data, const int dduID, const CSCDCCExaminer& binChecker);
132 
133  void calcEMUFractionHisto(const HistoId& result, const HistoId& set, const HistoId& subset);
134 
135  const bool getEMUHisto(const HistoId& histo, MonitorObject*& me);
136  const bool getFEDHisto(const HistoId& histo, const HwId& fedID, MonitorObject*& me);
137  const bool getDDUHisto(const HistoId& histo, const HwId& dduID, MonitorObject*& me);
138  const bool getCSCHisto(const HistoId& histo, const HwId& crateID, const HwId& dmbSlot, MonitorObject*& me);
139  const bool getCSCHisto(
140  const HistoId& histo, const HwId& crateID, const HwId& dmbSlot, const HwId& adId, MonitorObject*& me);
141  const bool getParHisto(const HistoId& histo, MonitorObject*& me);
142  void preProcessEvent();
143 
144  const bool getCSCFromMap(const unsigned int& crateId,
145  const unsigned int& dmbId,
146  unsigned int& cscType,
147  unsigned int& cscPosition) const;
148  void setEmuEventDisplayBit(MonitorObject*& mo, const unsigned int x, const unsigned int y, const unsigned int bit);
149  void resetEmuEventDisplays();
150 
153 
156 
157  std::map<uint32_t, uint32_t> L1ANumbers;
158  std::map<uint32_t, bool> fNotFirstEvent;
159  uint32_t L1ANumber;
160  uint32_t BXN;
161  uint32_t cntDMBs;
162  uint32_t cntCFEBs;
163  uint32_t cntALCTs;
164  uint32_t cntTMBs;
165 
166  uint16_t theFormatVersion;
167 
168  // bool fFirstEvent;
169  bool fCloseL1As; // Close L1A bit from DDU Trailer
171 
172  // token for input
174 
175  // ===================================================================================================
176  // Local ONLY stuff
177  // ===================================================================================================
178 
179 #ifdef DQMLOCAL
180 
181  public:
182  void processEvent(const char* data, const int32_t dataSize, const uint32_t errorStat, const int32_t nodeNumber);
183 
184 #endif
185 
186  // ===================================================================================================
187  // Global ONLY stuff
188  // ===================================================================================================
189 
190 #ifdef DQMGLOBAL
191 
192  private:
193  bool bCSCEventCounted;
194 
195  public:
196  void processEvent(const edm::Event& e, const edm::InputTag& inputTag);
197 
198 #endif
199  };
200 
201 } // namespace cscdqm
202 
203 #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.
tuple result
Definition: mps_fire.py:311
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.
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:79
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.