CMS 3D CMS Logo

GEMDQMEfficiencySourceBase.h
Go to the documentation of this file.
1 #ifndef DQM_GEM_GEMDQMEfficiencySourceBase_h
2 #define DQM_GEM_GEMDQMEfficiencySourceBase_h
3 
17 
19 public:
20  using MEMap = std::map<GEMDetId, MonitorElement*>;
21 
23 
24  // Re: region / St: station, La: layer, Ch: chamber parity, Et: eta partition
25  inline GEMDetId getReStKey(const int, const int);
26  inline GEMDetId getReStKey(const GEMDetId&);
27  inline GEMDetId getReStLaKey(const GEMDetId&);
28  inline GEMDetId getReStEtKey(const GEMDetId&);
29  inline GEMDetId getReStLaChKey(const GEMDetId&);
30  inline GEMDetId getKey(const GEMDetId&); // == getReStLaChEtKey
31 
33 
36 
37  std::tuple<bool, int, int> getChamberRange(const GEMStation*);
38  std::tuple<bool, int, int> getEtaPartitionRange(const GEMStation*);
39 
40  MonitorElement* bookChamber(DQMStore::IBooker&, const TString&, const TString&, const GEMStation*);
41  MonitorElement* bookChamberEtaPartition(DQMStore::IBooker&, const TString&, const TString&, const GEMStation*);
42 
43  bool skipGEMStation(const int);
44 
45  bool maskChamberWithError(const GEMDetId& chamber_id, const GEMOHStatusCollection*, const GEMVFATStatusCollection*);
46 
47  bool hasMEKey(const MEMap&, const GEMDetId&);
48  void fillME(MEMap&, const GEMDetId&, const double);
49  void fillME(MEMap&, const GEMDetId&, const double, const double);
50 
51  double clampWithAxis(const double, const TAxis* axis);
52  void fillMEWithinLimits(MonitorElement*, const double);
53  void fillMEWithinLimits(MonitorElement*, const double, const double);
54  void fillMEWithinLimits(MEMap&, const GEMDetId&, const double);
55  void fillMEWithinLimits(MEMap&, const GEMDetId&, const double, const double);
56 
57  template <typename T>
58  bool checkRefs(const std::vector<T*>&);
59 
62 
63  const bool kMonitorGE11_;
64  const bool kMonitorGE21_;
65  const bool kMonitorGE0_;
68 };
69 
70 template <typename T>
71 bool GEMDQMEfficiencySourceBase::checkRefs(const std::vector<T*>& refs) {
72  if (refs.empty())
73  return false;
74  for (T* each : refs) {
75  if (each == nullptr) {
76  return false;
77  }
78  }
79  return true;
80 }
81 
83  // region, ring, station, layer, chamber, ieta
84  return GEMDetId{region, 1, station, 0, 0, 0};
85 }
86 
88  return getReStKey(id.region(), id.station());
89 }
90 
92  return GEMDetId{id.region(), 1, id.station(), id.layer(), 0, 0};
93 }
94 
96  return GEMDetId{id.region(), 1, id.station(), 0, 0, id.ieta()};
97 }
98 
100  return GEMDetId{id.region(), 1, id.station(), id.layer(), id.chamber() % 2, 0};
101 }
102 
104  return GEMDetId{id.region(), 1, id.station(), id.layer(), id.chamber() % 2, id.ieta()};
105 }
106 
107 #endif // DQM_GEM_GEMDQMEfficiencySourceBase_h
constexpr int region() const
Definition: GEMDetId.h:171
MonitorElement * bookChamberEtaPartition(DQMStore::IBooker &, const TString &, const TString &, const GEMStation *)
MonitorElement * bookNumerator2D(DQMStore::IBooker &, MonitorElement *)
void fillMEWithinLimits(MonitorElement *, const double)
GEMDetId getReStLaChKey(const GEMDetId &)
std::tuple< bool, int, int > getEtaPartitionRange(const GEMStation *)
bool maskChamberWithError(const GEMDetId &chamber_id, const GEMOHStatusCollection *, const GEMVFATStatusCollection *)
const edm::EDGetTokenT< GEMVFATStatusCollection > kGEMVFATStatusCollectionToken_
std::map< GEMDetId, MonitorElement * > MEMap
void fillME(MEMap &, const GEMDetId &, const double)
MonitorElement * bookNumerator1D(DQMStore::IBooker &, MonitorElement *)
MonitorElement * bookChamber(DQMStore::IBooker &, const TString &, const TString &, const GEMStation *)
GEMDetId getReStKey(const int, const int)
std::string nameNumerator(const std::string &)
GEMDetId getKey(const GEMDetId &)
const edm::EDGetTokenT< GEMOHStatusCollection > kGEMOHStatusCollectionToken_
std::tuple< bool, int, int > getChamberRange(const GEMStation *)
bool checkRefs(const std::vector< T *> &)
bool hasMEKey(const MEMap &, const GEMDetId &)
double clampWithAxis(const double, const TAxis *axis)
GEMDetId getReStLaKey(const GEMDetId &)
GEMDetId getReStEtKey(const GEMDetId &)
long double T
GEMDQMEfficiencySourceBase(const edm::ParameterSet &)
A container for a generic type of digis indexed by some index, implemented with a map<IndexType...