CMS 3D CMS Logo

SusyPostProcessor.cc
Go to the documentation of this file.
1 #include <iostream>
2 
7 
8 using namespace std;
9 
10 const char* SusyPostProcessor::messageLoggerCatregory = "SusyDQMPostProcessor";
11 
13 {
14  iConfig = pSet;
15 
16  SUSYFolder = iConfig.getParameter<string>("folderName");
17  _quantile = iConfig.getParameter<double>("quantile");
18 
19 }
20 
22 }
23 
24 
25 
27 {
28  if(ME->getTH1()->GetEntries()>0.)
29  {
30  Quantile q(static_cast<const TH1*>(ME->getTH1()));
31  Float_t mean=q[q_value].first;
32  Float_t RMS=q[q_value].second;
33 
34  Float_t xLow=-5.5;
35  Float_t xUp=9.5;
36  Int_t NBin=15;
37 
38  if(mean>0.)
39  {
40  Float_t DBin=RMS*TMath::Sqrt(12.)/2.;
41  xLow=mean-int(mean/DBin+2)*DBin;
42  xUp=int(0.2*mean/DBin)*DBin+mean+5*DBin;
43  NBin=(xUp-xLow)/DBin;
44  }
45 
46  ibooker_.setCurrentFolder(ME->getPathname());
47  TString name=ME->getTH1()->GetName();
48  name+="_quant";
49  ME=ibooker_.book1D(name,"",NBin, xLow, xUp);
50  ME->Fill(mean-RMS);
51  ME->Fill(mean+RMS);
52  }
53 }
54 
55 
56 
58 {
59  // MET
60  //----------------------------------------------------------------------------
61  iget_.setCurrentFolder("JetMET/MET");
62 
63  Dirs = iget_.getSubdirs();
64 
65  std::vector<std::string> metFolders;
66 
67  metFolders.push_back("Uncleaned/");
68  metFolders.push_back("Cleaned/");
69 
70  //Need our own copy for thread safety
71  TF1 mygaus("mygaus","gaus");
72 
73  for (int i=0; i<int(Dirs.size()); i++) {
74 
75  std::string prefix = "dummy";
76 
77  if (size_t(Dirs[i].find("met")) != string::npos) prefix = "met";
78  if (size_t(Dirs[i].find("pfMet")) != string::npos) prefix = "pfMET";
79 
80  for (std::vector<std::string>::const_iterator ic=metFolders.begin();
81  ic!=metFolders.end(); ic++) {
82 
83  std::string dirName = Dirs[i] +"/" + *ic;
84 
85  MEx = iget_.get(dirName + "/"+"MEx");
86  MEy = iget_.get(dirName + "/"+"MEy");
87 
88  if (MEx && MEx->kind() == MonitorElement::DQM_KIND_TH1F) {
89  if (MEx->getTH1F()->GetEntries() > 50) MEx->getTH1F()->Fit(&mygaus, "q");
90  }
91 
92  if (MEy && MEy->kind() == MonitorElement::DQM_KIND_TH1F) {
93  if (MEy->getTH1F()->GetEntries() > 50) MEy->getTH1F()->Fit(&mygaus, "q");
94  }
95  }
96  }
97 
98 
99  // SUSY
100  //----------------------------------------------------------------------------
101  iget_.setCurrentFolder(SUSYFolder);
102  Dirs = iget_.getSubdirs();
103  for (int i=0; i<int(Dirs.size()); i++)
104  {
105  size_t found = Dirs[i].find("Alpha");
106  if (found!=string::npos) continue;
107  if(!iget_.dirExists(Dirs[i])){
108  edm::LogError(messageLoggerCatregory)<< "Directory "<<Dirs[i]<<" doesn't exist!!";
109  continue;
110  }
111  vector<MonitorElement*> histoVector = iget_.getContents(Dirs[i]);
112  for (int i=0; i<int(histoVector.size()); i++) {
113  QuantilePlots(histoVector[i],_quantile,ibook_);
114  }
115  }
116 }
T getParameter(std::string const &) const
const std::string & getPathname() const
get pathname of parent folder
TH1F * getTH1F() const
TH1 * getTH1() const
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:361
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
static const char * messageLoggerCatregory
void Fill(long long x)
Definition: ME.h:11
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
std::vector< MonitorElement * > getContents(Args &&...args)
Definition: DQMStore.h:192
MonitorElement * get(std::string const &path)
Definition: DQMStore.cc:303
~SusyPostProcessor() override
void QuantilePlots(MonitorElement *&, double, DQMStore::IBooker &)
bool dirExists(std::string const &path)
Definition: DQMStore.cc:343
std::vector< std::string > getSubdirs()
Definition: DQMStore.cc:325
void dqmEndJob(DQMStore::IBooker &, DQMStore::IGetter &) override
SusyPostProcessor(const edm::ParameterSet &pSet)