CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DQMGenericTnPClient.cc
Go to the documentation of this file.
10 
12 
13 #include <TString.h>
14 #include <TPRegexp.h>
15 
16 using namespace edm;
17 using namespace dqmTnP;
18 using namespace std;
19 
20 using vstring = std::vector<std::string>;
21 
23 public:
26 
28  ~DQMGenericTnPClient() override;
29  void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override{};
30  void endRun(const edm::Run& run, const edm::EventSetup& setup) override;
31  void calculateEfficiency(const std::string& dirName, const ParameterSet& pset);
32  void findAllSubdirectories(const std::string& dir, std::set<std::string>* myList, TString pattern);
33 
34 private:
36  TFile* plots;
39  bool verbose;
43 };
44 
46  : subDirs(pset.getUntrackedParameter<vstring>("subDirs", vstring())),
47  myDQMrootFolder(pset.getUntrackedParameter<std::string>("MyDQMrootFolder", "")),
48  verbose(pset.getUntrackedParameter<bool>("Verbose", false)),
49  efficiencies(pset.getUntrackedParameter<VParameterSet>("Efficiencies")) {
50  TString savePlotsInRootFileName = pset.getUntrackedParameter<string>("SavePlotsInRootFileName", "");
51  plots = savePlotsInRootFileName != "" ? new TFile(savePlotsInRootFileName, "recreate") : nullptr;
54 }
55 
57  TPRegexp metacharacters("[\\^\\$\\.\\*\\+\\?\\|\\(\\)\\{\\}\\[\\]]");
58 
60  if (!dqmStore) {
61  LogError("DQMGenericTnPClient") << "Could not find DQMStore service\n";
62  return;
63  }
65 
66  set<std::string> subDirSet;
67 
68  if (!myDQMrootFolder.empty())
69  subDirSet.insert(myDQMrootFolder);
70  else {
71  for (auto iSubDir = subDirs.begin(); iSubDir != subDirs.end(); ++iSubDir) {
72  std::string subDir = *iSubDir;
73  if (subDir[subDir.size() - 1] == '/')
74  subDir.erase(subDir.size() - 1);
75  if (TString(subDir).Contains(metacharacters)) {
76  const string::size_type shiftPos = subDir.rfind('/');
77  const string searchPath = subDir.substr(0, shiftPos);
78  const string pattern = subDir.substr(shiftPos + 1, subDir.length());
79  findAllSubdirectories(searchPath, &subDirSet, pattern);
80  } else {
81  subDirSet.insert(subDir);
82  }
83  }
84  }
85 
86  for (auto const& dirName : subDirSet) {
87  for (auto const& efficiencie : efficiencies) {
88  calculateEfficiency(dirName, efficiencie);
89  }
90  }
91 }
92 
94  //get hold of numerator and denominator histograms
95  string allMEname = dirName + "/" + pset.getUntrackedParameter<string>("DenominatorMEname");
96  string passMEname = dirName + "/" + pset.getUntrackedParameter<string>("NumeratorMEname");
97  MonitorElement* allME = dqmStore->get(allMEname);
98  MonitorElement* passME = dqmStore->get(passMEname);
99  if (allME == nullptr || passME == nullptr) {
100  LogDebug("DQMGenericTnPClient") << "Could not find MEs: " << allMEname << " or " << passMEname << endl;
101  return;
102  }
103  TH1* all = allME->getTH1();
104  TH1* pass = passME->getTH1();
105  //setup the fitter
106  string fitFunction = pset.getUntrackedParameter<string>("FitFunction");
107  AbstractFitter* fitter = nullptr;
108  if (fitFunction == "GaussianPlusLinear") {
109  GPLfitter->setup(pset.getUntrackedParameter<double>("ExpectedMean"),
110  pset.getUntrackedParameter<double>("FitRangeLow"),
111  pset.getUntrackedParameter<double>("FitRangeHigh"),
112  pset.getUntrackedParameter<double>("ExpectedSigma"));
113  fitter = GPLfitter;
114  } else if (fitFunction == "VoigtianPlusExponential") {
115  VPEfitter->setup(pset.getUntrackedParameter<double>("ExpectedMean"),
116  pset.getUntrackedParameter<double>("FitRangeLow"),
117  pset.getUntrackedParameter<double>("FitRangeHigh"),
118  pset.getUntrackedParameter<double>("ExpectedSigma"),
119  pset.getUntrackedParameter<double>("FixedWidth"));
120  fitter = VPEfitter;
121  } else {
122  LogError("DQMGenericTnPClient") << "Fit function: " << fitFunction << " does not exist" << endl;
123  return;
124  }
125  //check dimensions
126  int dimensions = all->GetDimension();
127  int massDimension = pset.getUntrackedParameter<int>("MassDimension");
128  if (massDimension > dimensions) {
129  LogError("DQMGenericTnPClient") << "Monitoring Element " << allMEname << " has smaller dimension than "
130  << massDimension << endl;
131  return;
132  }
133  //figure out directory and efficiency names
134  string effName = pset.getUntrackedParameter<string>("EfficiencyMEname");
135  string effDir = dirName;
136  string::size_type slashPos = effName.rfind('/');
137  if (string::npos != slashPos) {
138  effDir += "/" + effName.substr(0, slashPos);
139  effName.erase(0, slashPos + 1);
140  }
141  dqmStore->setCurrentFolder(effDir);
142  TString prefix(effDir.c_str());
143  prefix.ReplaceAll('/', '_');
144  //calculate and book efficiency
145  if (dimensions == 2) {
146  TProfile* eff = nullptr;
147  TProfile* effChi2 = nullptr;
148  TString error = fitter->calculateEfficiency(
149  (TH2*)pass, (TH2*)all, massDimension, eff, effChi2, plots ? prefix + effName.c_str() : "");
150  if (error != "") {
151  LogError("DQMGenericTnPClient") << error << endl;
152  return;
153  }
154  dqmStore->bookProfile(effName, eff);
155  dqmStore->bookProfile(effName + "Chi2", effChi2);
156  delete eff;
157  delete effChi2;
158  } else if (dimensions == 3) {
159  TProfile2D* eff = nullptr;
160  TProfile2D* effChi2 = nullptr;
161  TString error = fitter->calculateEfficiency(
162  (TH3*)pass, (TH3*)all, massDimension, eff, effChi2, plots ? prefix + effName.c_str() : "");
163  if (error != "") {
164  LogError("DQMGenericTnPClient") << error << endl;
165  return;
166  }
167  dqmStore->bookProfile2D(effName, eff);
168  dqmStore->bookProfile2D(effName + "Chi2", effChi2);
169  delete eff;
170  delete effChi2;
171  }
172 }
173 
175  delete GPLfitter;
176  if (plots) {
177  plots->Close();
178  }
179 }
180 
182  std::set<std::string>* myList,
183  TString pattern = "") {
184  if (!dqmStore->dirExists(dir)) {
185  LogError("DQMGenericTnPClient") << " DQMGenericTnPClient::findAllSubdirectories ==> Missing folder " << dir
186  << " !!!";
187  return;
188  }
189  TPRegexp nonPerlWildcard("\\w\\*|^\\*");
190  if (pattern != "") {
191  if (pattern.Contains(nonPerlWildcard))
192  pattern.ReplaceAll("*", ".*");
193  TPRegexp regexp(pattern);
194  dqmStore->cd(dir);
195  vector<string> foundDirs = dqmStore->getSubdirs();
196  for (auto iDir = foundDirs.begin(); iDir != foundDirs.end(); ++iDir) {
197  TString dirName = iDir->substr(iDir->rfind('/') + 1, iDir->length());
198  if (dirName.Contains(regexp))
199  findAllSubdirectories(*iDir, myList);
200  }
201  } else if (dqmStore->dirExists(dir)) {
202  myList->insert(dir);
203  dqmStore->cd(dir);
204  findAllSubdirectories(dir, myList, "*");
205  } else {
206  LogInfo("DQMGenericClient") << "Trying to find sub-directories of " << dir << " failed because " << dir
207  << " does not exist";
208  }
209  return;
210 }
211 
212 //define this as a plug-in
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * bookProfile2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, double lowZ, double highZ, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:399
dqm::legacy::DQMStore DQMStore
GaussianPlusLinearFitter * GPLfitter
virtual DQM_DEPRECATED std::vector< std::string > getSubdirs() const
Definition: DQMStore.cc:700
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:34
void setCurrentFolder(std::string const &fullpath) override
Definition: DQMStore.h:569
void endRun(const edm::Run &run, const edm::EventSetup &setup) override
void setup(double expectedMean_, double massLow, double massHigh, double expectedSigma_)
virtual bool dirExists(std::string const &path) const
Definition: DQMStore.cc:730
void calculateEfficiency(const std::string &dirName, const ParameterSet &pset)
Log< level::Error, false > LogError
uint16_t size_type
void setup(double expectedMean_, double massLow, double massHigh, double expectedSigma_, double width_)
DQMGenericTnPClient(const edm::ParameterSet &pset)
void findAllSubdirectories(const std::string &dir, std::set< std::string > *myList, TString pattern)
VoigtianPlusExponentialFitter * VPEfitter
static constexpr int verbose
MonitorElement * bookProfile(TString const &name, TString const &title, int nchX, double lowX, double highX, int, double lowY, double highY, char const *option="s", FUNC onbooking=NOOP())
Definition: DQMStore.h:322
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:673
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
def all
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
Log< level::Info, false > LogInfo
const VParameterSet efficiencies
std::vector< std::string > vstring
Definition: Schedule.cc:667
virtual TH1 * getTH1() const
Definition: Run.h:45
#define LogDebug(id)
dqm::legacy::MonitorElement MonitorElement