CMS 3D CMS Logo

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 
22 class DQMGenericTnPClient : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
23 public:
26 
28  ~DQMGenericTnPClient() override;
29  void analyze(const edm::Event& event, const edm::EventSetup& eventSetup) override{};
30  void beginRun(edm::Run const&, edm::EventSetup const&) override {}
31  void endRun(const edm::Run& run, const edm::EventSetup& setup) override;
32  void calculateEfficiency(const std::string& dirName, const ParameterSet& pset);
33  void findAllSubdirectories(const std::string& dir, std::set<std::string>* myList, TString pattern);
34 
35 private:
37  TFile* plots;
40  bool verbose;
44 };
45 
47  : subDirs(pset.getUntrackedParameter<vstring>("subDirs", vstring())),
48  myDQMrootFolder(pset.getUntrackedParameter<std::string>("MyDQMrootFolder", "")),
49  verbose(pset.getUntrackedParameter<bool>("Verbose", false)),
50  efficiencies(pset.getUntrackedParameter<VParameterSet>("Efficiencies")) {
51  usesResource("DQMStore");
52  TString savePlotsInRootFileName = pset.getUntrackedParameter<string>("SavePlotsInRootFileName", "");
53  plots = savePlotsInRootFileName != "" ? new TFile(savePlotsInRootFileName, "recreate") : nullptr;
56 }
57 
59  TPRegexp metacharacters("[\\^\\$\\.\\*\\+\\?\\|\\(\\)\\{\\}\\[\\]]");
60 
62  if (!dqmStore) {
63  LogError("DQMGenericTnPClient") << "Could not find DQMStore service\n";
64  return;
65  }
67 
68  set<std::string> subDirSet;
69 
70  if (!myDQMrootFolder.empty())
71  subDirSet.insert(myDQMrootFolder);
72  else {
73  for (auto iSubDir = subDirs.begin(); iSubDir != subDirs.end(); ++iSubDir) {
74  std::string subDir = *iSubDir;
75  if (subDir[subDir.size() - 1] == '/')
76  subDir.erase(subDir.size() - 1);
77  if (TString(subDir).Contains(metacharacters)) {
78  const string::size_type shiftPos = subDir.rfind('/');
79  const string searchPath = subDir.substr(0, shiftPos);
80  const string pattern = subDir.substr(shiftPos + 1, subDir.length());
81  findAllSubdirectories(searchPath, &subDirSet, pattern);
82  } else {
83  subDirSet.insert(subDir);
84  }
85  }
86  }
87 
88  for (auto const& dirName : subDirSet) {
89  for (auto const& efficiencie : efficiencies) {
90  calculateEfficiency(dirName, efficiencie);
91  }
92  }
93 }
94 
96  //get hold of numerator and denominator histograms
97  string allMEname = dirName + "/" + pset.getUntrackedParameter<string>("DenominatorMEname");
98  string passMEname = dirName + "/" + pset.getUntrackedParameter<string>("NumeratorMEname");
99  MonitorElement* allME = dqmStore->get(allMEname);
100  MonitorElement* passME = dqmStore->get(passMEname);
101  if (allME == nullptr || passME == nullptr) {
102  LogDebug("DQMGenericTnPClient") << "Could not find MEs: " << allMEname << " or " << passMEname << endl;
103  return;
104  }
105  TH1* all = allME->getTH1();
106  TH1* pass = passME->getTH1();
107  //setup the fitter
108  string fitFunction = pset.getUntrackedParameter<string>("FitFunction");
109  AbstractFitter* fitter = nullptr;
110  if (fitFunction == "GaussianPlusLinear") {
111  GPLfitter->setup(pset.getUntrackedParameter<double>("ExpectedMean"),
112  pset.getUntrackedParameter<double>("FitRangeLow"),
113  pset.getUntrackedParameter<double>("FitRangeHigh"),
114  pset.getUntrackedParameter<double>("ExpectedSigma"));
115  fitter = GPLfitter;
116  } else if (fitFunction == "VoigtianPlusExponential") {
117  VPEfitter->setup(pset.getUntrackedParameter<double>("ExpectedMean"),
118  pset.getUntrackedParameter<double>("FitRangeLow"),
119  pset.getUntrackedParameter<double>("FitRangeHigh"),
120  pset.getUntrackedParameter<double>("ExpectedSigma"),
121  pset.getUntrackedParameter<double>("FixedWidth"));
122  fitter = VPEfitter;
123  } else {
124  LogError("DQMGenericTnPClient") << "Fit function: " << fitFunction << " does not exist" << endl;
125  return;
126  }
127  //check dimensions
128  int dimensions = all->GetDimension();
129  int massDimension = pset.getUntrackedParameter<int>("MassDimension");
130  if (massDimension > dimensions) {
131  LogError("DQMGenericTnPClient") << "Monitoring Element " << allMEname << " has smaller dimension than "
132  << massDimension << endl;
133  return;
134  }
135  //figure out directory and efficiency names
136  string effName = pset.getUntrackedParameter<string>("EfficiencyMEname");
137  string effDir = dirName;
138  string::size_type slashPos = effName.rfind('/');
139  if (string::npos != slashPos) {
140  effDir += "/" + effName.substr(0, slashPos);
141  effName.erase(0, slashPos + 1);
142  }
143  dqmStore->setCurrentFolder(effDir);
144  TString prefix(effDir.c_str());
145  prefix.ReplaceAll('/', '_');
146  //calculate and book efficiency
147  if (dimensions == 2) {
148  TProfile* eff = nullptr;
149  TProfile* effChi2 = nullptr;
150  TString error = fitter->calculateEfficiency(
151  (TH2*)pass, (TH2*)all, massDimension, eff, effChi2, plots ? prefix + effName.c_str() : "");
152  if (error != "") {
153  LogError("DQMGenericTnPClient") << error << endl;
154  return;
155  }
156  dqmStore->bookProfile(effName, eff);
157  dqmStore->bookProfile(effName + "Chi2", effChi2);
158  delete eff;
159  delete effChi2;
160  } else if (dimensions == 3) {
161  TProfile2D* eff = nullptr;
162  TProfile2D* effChi2 = nullptr;
163  TString error = fitter->calculateEfficiency(
164  (TH3*)pass, (TH3*)all, massDimension, eff, effChi2, plots ? prefix + effName.c_str() : "");
165  if (error != "") {
166  LogError("DQMGenericTnPClient") << error << endl;
167  return;
168  }
169  dqmStore->bookProfile2D(effName, eff);
170  dqmStore->bookProfile2D(effName + "Chi2", effChi2);
171  delete eff;
172  delete effChi2;
173  }
174 }
175 
177  delete GPLfitter;
178  if (plots) {
179  plots->Close();
180  }
181 }
182 
184  std::set<std::string>* myList,
185  TString pattern = "") {
186  if (!dqmStore->dirExists(dir)) {
187  LogError("DQMGenericTnPClient") << " DQMGenericTnPClient::findAllSubdirectories ==> Missing folder " << dir
188  << " !!!";
189  return;
190  }
191  TPRegexp nonPerlWildcard("\\w\\*|^\\*");
192  if (pattern != "") {
193  if (pattern.Contains(nonPerlWildcard))
194  pattern.ReplaceAll("*", ".*");
195  TPRegexp regexp(pattern);
196  dqmStore->cd(dir);
197  vector<string> foundDirs = dqmStore->getSubdirs();
198  for (auto iDir = foundDirs.begin(); iDir != foundDirs.end(); ++iDir) {
199  TString dirName = iDir->substr(iDir->rfind('/') + 1, iDir->length());
200  if (dirName.Contains(regexp))
202  }
203  } else if (dqmStore->dirExists(dir)) {
204  myList->insert(dir);
205  dqmStore->cd(dir);
207  } else {
208  LogInfo("DQMGenericClient") << "Trying to find sub-directories of " << dir << " failed because " << dir
209  << " does not exist";
210  }
211  return;
212 }
213 
214 //define this as a plug-in
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:485
dqm::legacy::DQMStore DQMStore
GaussianPlusLinearFitter * GPLfitter
def all(container)
workaround iterator generators for ROOT classes
Definition: cmstools.py:25
vector< string > vstring
Definition: ExoticaDQM.cc:7
bool verbose
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:35
double fitFunction(double *x, double *par)
void setCurrentFolder(std::string const &fullpath) override
Definition: DQMStore.h:656
void endRun(const edm::Run &run, const edm::EventSetup &setup) override
virtual bool dirExists(std::string const &path) const
Definition: DQMStore.cc:769
void setup(double expectedMean_, double massLow, double massHigh, double expectedSigma_)
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
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:408
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
Log< level::Info, false > LogInfo
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:712
virtual TH1 * getTH1() const
HLT enums.
void beginRun(edm::Run const &, edm::EventSetup const &) override
const VParameterSet efficiencies
Definition: event.py:1
Definition: Run.h:45
#define LogDebug(id)
dqm::legacy::MonitorElement MonitorElement
virtual DQM_DEPRECATED std::vector< std::string > getSubdirs() const
Definition: DQMStore.cc:739