CMS 3D CMS Logo

DiMuonMassBiasClient.cc
Go to the documentation of this file.
1 // user includes
7 
8 // RooFit includes
9 #include "TCanvas.h"
10 #include "RooAddPdf.h"
11 #include "RooDataHist.h"
12 #include "RooExponential.h"
13 #include "RooGaussian.h"
14 #include "RooPlot.h"
15 #include "RooRealVar.h"
16 #include "RooVoigtian.h"
17 #include "RooCBShape.h"
18 
19 //-----------------------------------------------------------------------------------
21  : TopFolder_(iConfig.getParameter<std::string>("FolderName")),
22  fitBackground_(iConfig.getParameter<bool>("fitBackground")),
23  useRooCBShape_(iConfig.getParameter<bool>("useRooCBShape")),
24  useRooCMSShape_(iConfig.getParameter<bool>("useRooCMSShape")),
25  debugMode_(iConfig.getParameter<bool>("debugMode")),
26  MEtoHarvest_(iConfig.getParameter<std::vector<std::string>>("MEtoHarvest"))
27 //-----------------------------------------------------------------------------------
28 {
29  edm::LogInfo("DiMuonMassBiasClient") << "DiMuonMassBiasClient::Constructing DiMuonMassBiasClient ";
30 
31  // fill the parameters for the fit
36 
37  if (debugMode_) {
38  edm::LogPrint("DiMuonMassBiasClient")
39  << "mean: " << meanConfig_[0] << " (" << meanConfig_[1] << "," << meanConfig_[2] << ") " << std::endl;
40  edm::LogPrint("DiMuonMassBiasClient")
41  << "width: " << widthConfig_[0] << " (" << widthConfig_[1] << "," << widthConfig_[2] << ")" << std::endl;
42  edm::LogPrint("DiMuonMassBiasClient")
43  << "sigma: " << sigmaConfig_[0] << " (" << sigmaConfig_[1] << "," << sigmaConfig_[2] << ")" << std::endl;
44  }
45 }
46 
47 //-----------------------------------------------------------------------------------
49 //-----------------------------------------------------------------------------------
50 {
51  edm::LogInfo("DiMuonMassBiasClient") << "DiMuonMassBiasClient::Deleting DiMuonMassBiasClient ";
52 }
53 
54 //-----------------------------------------------------------------------------------
56 //-----------------------------------------------------------------------------------
57 {
58  edm::LogInfo("DiMuonMassBiasClient") << "DiMuonMassBiasClient::beginJob done";
59 }
60 
61 //-----------------------------------------------------------------------------------
63 //-----------------------------------------------------------------------------------
64 {
65  edm::LogInfo("DiMuonMassBiasClient") << "DiMuonMassBiasClient:: Begining of Run";
66 }
67 
68 //-----------------------------------------------------------------------------------
70 //-----------------------------------------------------------------------------------
71 {
72  iBooker.setCurrentFolder(TopFolder_ + "/DiMuonMassBiasMonitor/MassBias/Profiles");
73  for (const auto& [key, ME] : harvestTargets_) {
74  if (ME == nullptr) {
75  edm::LogError("DiMuonMassBiasClient") << "could not find MonitorElement for key: " << key << std::endl;
76  continue;
77  }
78 
79  const auto& title = ME->getTitle();
80  const auto& xtitle = ME->getAxisTitle(1);
81  const auto& ytitle = ME->getAxisTitle(2);
82 
83  const auto& nxbins = ME->getNbinsX();
84  const auto& xmin = ME->getAxisMin(1);
85  const auto& xmax = ME->getAxisMax(1);
86 
87  MonitorElement* meanToBook =
88  iBooker.book1D(("Mean" + key), (title + ";" + xtitle + ";" + ytitle), nxbins, xmin, xmax);
89  meanProfiles_.insert({key, meanToBook});
90 
91  MonitorElement* sigmaToBook =
92  iBooker.book1D(("Sigma" + key), (title + ";" + xtitle + ";" + "#sigma of " + ytitle), nxbins, xmin, xmax);
93  widthProfiles_.insert({key, sigmaToBook});
94  }
95 }
96 
97 //-----------------------------------------------------------------------------------
99 //-----------------------------------------------------------------------------------
100 {
101  std::string inFolder = TopFolder_ + "/DiMuonMassBiasMonitor/MassBias/";
102 
103  //loop on the list of histograms to harvest
104  for (const auto& hname : MEtoHarvest_) {
105  MonitorElement* toHarvest = iGetter.get(inFolder + hname);
106 
107  if (toHarvest == nullptr) {
108  edm::LogError("DiMuonMassBiasClient") << "could not find input MonitorElement: " << inFolder + hname << std::endl;
109  continue;
110  }
111 
112  harvestTargets_.insert({hname, toHarvest});
113  }
114 }
115 
116 //-----------------------------------------------------------------------------------
118 //-----------------------------------------------------------------------------------
119 {
120  edm::LogInfo("DiMuonMassBiasClient") << "DiMuonMassBiasClient::endLuminosityBlock";
121 
122  getMEsToHarvest(igetter);
123  bookMEs(ibooker);
124 
125  for (const auto& [key, ME] : harvestTargets_) {
126  if (debugMode_)
127  edm::LogPrint("DiMuonMassBiasClient") << "dealing with key: " << key << std::endl;
128  TH2F* bareHisto = ME->getTH2F();
129  for (int bin = 1; bin <= ME->getNbinsX(); bin++) {
130  const auto& xaxis = bareHisto->GetXaxis();
131  const auto& low_edge = xaxis->GetBinLowEdge(bin);
132  const auto& high_edge = xaxis->GetBinUpEdge(bin);
133 
134  if (debugMode_)
135  edm::LogPrint("DiMuonMassBiasClient") << "dealing with bin: " << bin << " range: (" << std::setprecision(2)
136  << low_edge << "," << std::setprecision(2) << high_edge << ")";
137  TH1D* Proj = bareHisto->ProjectionY(Form("%s_proj_%i", key.c_str(), bin), bin, bin);
138  Proj->SetTitle(Form("%s #in (%.2f,%.2f), bin: %i", Proj->GetTitle(), low_edge, high_edge, bin));
139 
141 
142  if (results.isInvalid()) {
143  edm::LogWarning("DiMuonMassBiasClient") << "the current bin has invalid data" << std::endl;
144  continue;
145  }
146 
147  // fill the mean profiles
148  const Measurement1D& bias = results.getBias();
149  meanProfiles_[key]->setBinContent(bin, bias.value());
150  meanProfiles_[key]->setBinError(bin, bias.error());
151 
152  // fill the width profiles
153  const Measurement1D& width = results.getWidth();
154  widthProfiles_[key]->setBinContent(bin, width.value());
155  widthProfiles_[key]->setBinError(bin, width.error());
156  }
157  }
158 }
159 
160 //-----------------------------------------------------------------------------------
162 //-----------------------------------------------------------------------------------
163 {
164  if (hist->GetEntries() < diMuonMassBias::minimumHits) {
165  edm::LogWarning("DiMuonMassBiasClient") << " Input histogram:" << hist->GetName() << " has not enough entries ("
166  << hist->GetEntries() << ") for a meaningful Voigtian fit!\n"
167  << "Skipping!";
168 
169  return diMuonMassBias::fitOutputs(Measurement1D(0., 0.), Measurement1D(0., 0.));
170  }
171 
172  TCanvas* c1 = new TCanvas();
173  if (debugMode_) {
174  c1->Clear();
175  c1->SetLeftMargin(0.15);
176  c1->SetRightMargin(0.10);
177  }
178 
179  // silence messages
180  RooMsgService::instance().setGlobalKillBelow(RooFit::FATAL);
181 
182  Double_t xmin = hist->GetXaxis()->GetXmin();
183  Double_t xmax = hist->GetXaxis()->GetXmax();
184 
185  if (debugMode_) {
186  edm::LogPrint("DiMuonMassBiasClient") << "fitting range: (" << xmin << "-" << xmax << ")" << std::endl;
187  }
188 
189  RooRealVar InvMass("InvMass", "di-muon mass M(#mu^{+}#mu^{-}) [GeV]", xmin, xmax);
190  std::unique_ptr<RooPlot> frame{InvMass.frame()};
191  RooDataHist datahist("datahist", "datahist", InvMass, RooFit::Import(*hist));
192  datahist.plotOn(frame.get());
193 
194  // parameters of the Voigtian
195  RooRealVar mean("#mu", "mean", meanConfig_[0], meanConfig_[1], meanConfig_[2]); //90.0, 60.0, 120.0 (for Z)
196  RooRealVar width("width", "width", widthConfig_[0], widthConfig_[1], widthConfig_[2]); // 5.0, 0.0, 120.0 (for Z)
197  RooRealVar sigma("#sigma", "sigma", sigmaConfig_[0], sigmaConfig_[1], sigmaConfig_[2]); // 5.0, 0.0, 120.0 (for Z)
198  RooVoigtian voigt("voigt", "voigt", InvMass, mean, width, sigma);
199 
200  // parameters of the Crystal-ball
201  RooRealVar peakCB("peakCB", "peakCB", meanConfig_[0], meanConfig_[1], meanConfig_[2]);
202  RooRealVar sigmaCB("#sigma", "sigma", sigmaConfig_[0], sigmaConfig_[1], sigmaConfig_[2]);
203  RooRealVar alphaCB("#alpha", "alpha", 1., 0., 10.);
204  RooRealVar nCB("n", "n", 1., 0., 100.);
205  RooCBShape crystalball("crystalball", "crystalball", InvMass, peakCB, sigmaCB, alphaCB, nCB);
206 
207  // for the simple background fit
208  RooRealVar lambda("#lambda", "slope", 0., -50., 50.);
209  RooExponential expo("expo", "expo", InvMass, lambda);
210 
211  // for the more refined background fit
212  RooRealVar exp_alpha("#alpha", "alpha", 40.0, 20.0, 160.0);
213  RooRealVar exp_beta("#beta", "beta", 0.05, 0.0, 2.0);
214  RooRealVar exp_gamma("#gamma", "gamma", 0.02, 0.0, 0.1);
215  RooRealVar exp_peak("peak", "peak", meanConfig_[0]);
216  RooCMSShape exp_pdf("exp_pdf", "bkg shape", InvMass, exp_alpha, exp_beta, exp_gamma, exp_peak);
217 
218  // define the signal and background fractions
219  RooRealVar b("N_{b}", "Number of background events", 0, hist->GetEntries() / 10.);
220  RooRealVar s("N_{s}", "Number of signal events", 0, hist->GetEntries());
221 
222  if (fitBackground_) {
223  RooArgList listPdf;
224  if (useRooCBShape_) {
225  if (useRooCMSShape_) {
226  // crystal-ball + CMS-shape fit
227  listPdf.add(crystalball);
228  listPdf.add(exp_pdf);
229  } else {
230  // crystal-ball + exponential fit
231  listPdf.add(crystalball);
232  listPdf.add(expo);
233  }
234  } else {
235  if (useRooCMSShape_) {
236  // voigtian + CMS-shape fit
237  listPdf.add(voigt);
238  listPdf.add(exp_pdf);
239  } else {
240  // voigtian + exponential fit
241  listPdf.add(voigt);
242  listPdf.add(expo);
243  }
244  }
245 
246  RooAddPdf fullModel("fullModel", "Signal + Background Model", listPdf, RooArgList(s, b));
247  fullModel.fitTo(datahist, RooFit::PrintLevel(-1));
248  fullModel.plotOn(frame.get(), RooFit::LineColor(kRed));
249  if (useRooCMSShape_) {
250  fullModel.plotOn(frame.get(), RooFit::Components(exp_pdf), RooFit::LineStyle(kDashed)); //Other option
251  } else {
252  fullModel.plotOn(frame.get(), RooFit::Components(expo), RooFit::LineStyle(kDashed)); //Other option
253  }
254  fullModel.paramOn(frame.get(), RooFit::Layout(0.65, 0.90, 0.90));
255  } else {
256  if (useRooCBShape_) {
257  // use crystal-ball for a fit-only signal
258  crystalball.fitTo(datahist, RooFit::PrintLevel(-1));
259  crystalball.plotOn(frame.get(), RooFit::LineColor(kRed)); //this will show fit overlay on canvas
260  crystalball.paramOn(frame.get(),
261  RooFit::Layout(0.65, 0.90, 0.90)); //this will display the fit parameters on canvas
262  } else {
263  // use voigtian for a fit-only signal
264  voigt.fitTo(datahist, RooFit::PrintLevel(-1));
265  voigt.plotOn(frame.get(), RooFit::LineColor(kRed)); //this will show fit overlay on canvas
266  voigt.paramOn(frame.get(), RooFit::Layout(0.65, 0.90, 0.90)); //this will display the fit parameters on canvas
267  }
268  }
269 
270  // Redraw data on top and print / store everything
271  datahist.plotOn(frame.get());
272  frame->GetYaxis()->SetTitle("n. of events");
273  TString histName = hist->GetName();
274  frame->SetName("frame" + histName);
275  frame->SetTitle(hist->GetTitle());
276  frame->Draw();
277 
278  if (debugMode_) {
279  c1->Print("fit_debug" + histName + ".pdf");
280  }
281  delete c1;
282 
283  float mass_mean = useRooCBShape_ ? peakCB.getVal() : mean.getVal();
284  float mass_sigma = useRooCBShape_ ? sigmaCB.getVal() : sigma.getVal();
285 
286  float mass_mean_err = useRooCBShape_ ? peakCB.getError() : mean.getError();
287  float mass_sigma_err = useRooCBShape_ ? sigmaCB.getError() : sigma.getError();
288 
289  Measurement1D resultM(mass_mean, mass_mean_err);
290  Measurement1D resultW(mass_sigma, mass_sigma_err);
291 
292  return diMuonMassBias::fitOutputs(resultM, resultW);
293 }
294 
295 //-----------------------------------------------------------------------------------
297 //-----------------------------------------------------------------------------------
298 {
300  desc.add<std::string>("FolderName", "DiMuonMassBiasMonitor");
301  desc.add<bool>("fitBackground", false);
302  desc.add<bool>("useRooCMSShape", false);
303  desc.add<bool>("useRooCBShape", false);
304  desc.add<bool>("debugMode", false);
305 
307  fit_par.add<std::vector<double>>("mean_par",
310  std::numeric_limits<float>::max()}); // par = mean
311 
312  fit_par.add<std::vector<double>>("width_par",
315  std::numeric_limits<float>::max()}); // par = width
316 
317  fit_par.add<std::vector<double>>("sigma_par",
320  std::numeric_limits<float>::max()}); // par = sigma
321 
322  desc.add<edm::ParameterSetDescription>("fit_par", fit_par);
323 
324  desc.add<std::vector<std::string>>("MEtoHarvest",
325  {"DiMuMassVsMuMuPhi",
326  "DiMuMassVsMuMuEta",
327  "DiMuMassVsMuPlusPhi",
328  "DiMuMassVsMuPlusEta",
329  "DiMuMassVsMuMinusPhi",
330  "DiMuMassVsMuMinusEta",
331  "DiMuMassVsMuMuDeltaEta",
332  "DiMuMassVsCosThetaCS"});
333  descriptions.addWithDefaultLabel(desc);
334 }
335 
void getMEsToHarvest(DQMStore::IGetter &igetter)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
static PFTauRenderPlugin instance
diMuonMassBias::fitOutputs fitLineShape(TH1 *hist, const bool &fitBackground=false) const
~DiMuonMassBiasClient() override
Destructor.
const std::string TopFolder_
void beginRun(edm::Run const &run, edm::EventSetup const &eSetup) override
BeginRun.
Log< level::Error, false > LogError
std::map< std::string, MonitorElement * > harvestTargets_
std::map< std::string, MonitorElement * > meanProfiles_
std::map< std::string, MonitorElement * > widthProfiles_
Definition: ME.h:11
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< std::string > MEtoHarvest_
Log< level::Warning, true > LogPrint
void dqmEndJob(DQMStore::IBooker &ibooker_, DQMStore::IGetter &igetter_) override
EndJob.
void beginJob(void) override
BeginJob.
Log< level::Info, false > LogInfo
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double b
Definition: hdecay.h:118
static constexpr int minimumHits
virtual MonitorElement * get(std::string const &fullpath) const
Definition: DQMStore.cc:697
DiMuonMassBiasClient(const edm::ParameterSet &ps)
Constructor.
double value() const
Definition: Measurement1D.h:25
void bookMEs(DQMStore::IBooker &ibooker)
book MEs
results
Definition: mysort.py:8
double error() const
Definition: Measurement1D.h:27
double crystalball(double *x, double *par)
Definition: myFunctions.cc:98
Log< level::Warning, false > LogWarning
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
void fillArrayF(float *x, const edm::ParameterSet &cfg, const char *name)
Definition: Run.h:45