CMS 3D CMS Logo

SmartSelectionMonitor.h
Go to the documentation of this file.
1 #ifndef ALIGNMENT_OFFLINEVALIDATION_JETHTANALYZER_SMARTSELECTIONMONITOR_H
2 #define ALIGNMENT_OFFLINEVALIDATION_JETHTANALYZER_SMARTSELECTIONMONITOR_H
3 
4 // system include files
5 #include <iostream>
6 #include <string>
7 #include <map>
8 #include <algorithm>
9 #include <vector>
10 #include <memory>
11 #include <unordered_map>
12 
13 // user include files
14 #include "TH1D.h"
15 #include "TH2D.h"
16 #include "TProfile.h"
17 #include "TString.h"
18 #include "TROOT.h"
19 
21 public:
23  ~SmartSelectionMonitor() = default;
24 
25  typedef std::unordered_map<std::string, std::map<std::string, TH1*>*> Monitor_t;
26 
27  //short getters
28  inline Monitor_t& getAllMonitors() { return allMonitors_; }
29 
30  //checks if base Histo Exist
32  if (allMonitors_.find(histo) == allMonitors_.end())
33  return false;
34  return true;
35  }
36 
37  //checks if tag Exist for a given histo
38  inline bool hasTag(std::map<std::string, TH1*>* map, std::string tag) {
39  if (map->find(tag) != map->end())
40  return true;
41  if (map->find("all") == map->end())
42  return false;
43 
44  TH1* base = (*map)["all"];
45  TString allName = base->GetName();
46  TString name = tag + "_" + allName.Data();
47  TH1* h = (TH1*)base->Clone(name);
48  h->SetName(name);
49  h->SetTitle(name);
50  h->Reset("ICE");
51  h->SetDirectory(gROOT);
52  (*map)[tag] = h;
53  return true;
54  }
55 
56  //checks if tag Exist for a given histo
58  if (!hasBaseHisto(histo))
59  return false;
60  std::map<std::string, TH1*>* map = allMonitors_[histo];
61  return hasTag(map, tag);
62  }
63 
64  //get histo
65  inline TH1* getHisto(std::string histo, std::string tag = "all") {
66  if (!hasBaseHisto(histo))
67  return nullptr;
68  std::map<std::string, TH1*>* map = allMonitors_[histo];
69  if (!hasTag(map, tag))
70  return nullptr;
71  return (*map)[tag];
72  }
73 
74  //write all histo
75  inline void Write() {
76  for (Monitor_t::iterator it = allMonitors_.begin(); it != allMonitors_.end(); it++) {
77  std::map<std::string, TH1*>* map = it->second;
78  bool neverFilled = true;
79 
80  for (std::map<std::string, TH1*>::iterator h = map->begin(); h != map->end(); h++) {
81  if (!(h->second)) {
82  printf("histo = %30s %15s IS NULL", it->first.c_str(), h->first.c_str());
83  continue;
84  }
85  if (h->second->GetEntries() > 0)
86  neverFilled = false;
87 
88  if (h->first == "all") {
89  h->second->SetName(Form("%s_%s", h->first.c_str(), h->second->GetName()));
90  }
91  h->second->Write();
92  }
93 
94  if (neverFilled) {
95  printf(
96  "SmartSelectionMonitor: histo = '%s' is empty for all categories, you may want to cleanup your project to "
97  "remove this histogram\n",
98  it->first.c_str());
99  }
100  }
101  }
102 
103  //scale all histo by w
104  inline void Scale(double w) {
105  for (Monitor_t::iterator it = allMonitors_.begin(); it != allMonitors_.end(); it++) {
106  std::map<std::string, TH1*>* map = it->second;
107  for (std::map<std::string, TH1*>::iterator h = map->begin(); h != map->end(); h++) {
108  if (!(h->second)) {
109  continue;
110  }
111  h->second->Scale(w);
112  }
113  }
114  }
115 
116  //takes care of filling an histogram
117  bool fillHisto(std::string name, std::string tag, double valx, double weight, bool useBinWidth = false);
118  bool fillHisto(std::string name, std::string tag, double valx, double valy, double weight, bool useBinWidth = false);
119  bool fillProfile(std::string name, std::string tag, double valx, double valy, double weight);
120 
121  bool fillHisto(std::string name, std::vector<std::string> tags, double valx, double weight, bool useBinWidth = false);
123  std::vector<std::string> tags,
124  double valx,
125  double valy,
126  double weight,
127  bool useBinWidth = false);
128  bool fillProfile(std::string name, std::vector<std::string> tags, double valx, double valy, double weight);
129 
131  std::vector<std::string> tags,
132  double valx,
133  std::vector<double> weights,
134  bool useBinWidth = false);
136  std::vector<std::string> tags,
137  double valx,
138  double valy,
139  std::vector<double> weights,
140  bool useBinWidth = false);
141  bool fillProfile(
142  std::string name, std::vector<std::string> tags, double valx, double valy, std::vector<double> weights);
143 
144  //short inits the monitor plots for a new step
146 
147  //short add new histogram
148  TH1* addHistogram(TH1* h, std::string tag);
149  TH1* addHistogram(TH1* h);
150 
151 private:
152  //all the selection step monitors
154 };
155 
156 #endif
TH1 * getHisto(std::string histo, std::string tag="all")
bool fillHisto(std::string name, std::string tag, double valx, double weight, bool useBinWidth=false)
bool hasTag(std::map< std::string, TH1 *> *map, std::string tag)
T w() const
std::unordered_map< std::string, std::map< std::string, TH1 * > * > Monitor_t
Definition: weight.py:1
void initMonitorForStep(std::string tag)
TH1 * addHistogram(TH1 *h, std::string tag)
~SmartSelectionMonitor()=default
bool fillProfile(std::string name, std::string tag, double valx, double valy, double weight)
bool hasBaseHisto(std::string histo)
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
bool hasTag(std::string histo, std::string tag)