10 #include "RooAddPdf.h" 11 #include "RooDataHist.h" 12 #include "RooExponential.h" 13 #include "RooGaussian.h" 15 #include "RooRealVar.h" 16 #include "RooVoigtian.h" 17 #include "RooCBShape.h" 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")),
29 edm::LogInfo(
"DiMuonMassBiasClient") <<
"DiMuonMassBiasClient::Constructing DiMuonMassBiasClient ";
51 edm::LogInfo(
"DiMuonMassBiasClient") <<
"DiMuonMassBiasClient::Deleting DiMuonMassBiasClient ";
58 edm::LogInfo(
"DiMuonMassBiasClient") <<
"DiMuonMassBiasClient::beginJob done";
65 edm::LogInfo(
"DiMuonMassBiasClient") <<
"DiMuonMassBiasClient:: Begining of Run";
75 edm::LogError(
"DiMuonMassBiasClient") <<
"could not find MonitorElement for key: " <<
key << std::endl;
79 const auto&
title =
ME->getTitle();
80 const auto&
xtitle =
ME->getAxisTitle(1);
81 const auto&
ytitle =
ME->getAxisTitle(2);
83 const auto& nxbins =
ME->getNbinsX();
84 const auto&
xmin =
ME->getAxisMin(1);
85 const auto&
xmax =
ME->getAxisMax(1);
107 if (toHarvest ==
nullptr) {
108 edm::LogError(
"DiMuonMassBiasClient") <<
"could not find input MonitorElement: " << inFolder + hname << std::endl;
120 edm::LogInfo(
"DiMuonMassBiasClient") <<
"DiMuonMassBiasClient::endLuminosityBlock";
127 edm::LogPrint(
"DiMuonMassBiasClient") <<
"dealing with key: " <<
key << std::endl;
128 TH2F* bareHisto =
ME->getTH2F();
130 const auto&
xaxis = bareHisto->GetXaxis();
131 const auto& low_edge =
xaxis->GetBinLowEdge(
bin);
132 const auto& high_edge =
xaxis->GetBinUpEdge(
bin);
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));
143 edm::LogWarning(
"DiMuonMassBiasClient") <<
"the current bin has invalid data" << std::endl;
165 edm::LogWarning(
"DiMuonMassBiasClient") <<
" Input histogram:" <<
hist->GetName() <<
" has not enough entries (" 166 <<
hist->GetEntries() <<
") for a meaningful Voigtian fit!\n" 172 TCanvas*
c1 =
new TCanvas();
175 c1->SetLeftMargin(0.15);
176 c1->SetRightMargin(0.10);
182 Double_t
xmin =
hist->GetXaxis()->GetXmin();
183 Double_t
xmax =
hist->GetXaxis()->GetXmax();
186 edm::LogPrint(
"DiMuonMassBiasClient") <<
"fitting range: (" <<
xmin <<
"-" <<
xmax <<
")" << std::endl;
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());
195 RooRealVar
mean(
"#mu",
"mean", meanConfig_[0], meanConfig_[1], meanConfig_[2]);
196 RooRealVar
width(
"width",
"width", widthConfig_[0], widthConfig_[1], widthConfig_[2]);
197 RooRealVar sigma(
"#sigma",
"sigma", sigmaConfig_[0], sigmaConfig_[1], sigmaConfig_[2]);
198 RooVoigtian voigt(
"voigt",
"voigt", InvMass,
mean,
width, sigma);
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);
208 RooRealVar lambda(
"#lambda",
"slope", 0., -50., 50.);
209 RooExponential expo(
"expo",
"expo", InvMass, lambda);
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);
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());
222 if (fitBackground_) {
224 if (useRooCBShape_) {
225 if (useRooCMSShape_) {
228 listPdf.add(exp_pdf);
235 if (useRooCMSShape_) {
238 listPdf.add(exp_pdf);
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));
252 fullModel.plotOn(
frame.get(), RooFit::Components(expo), RooFit::LineStyle(kDashed));
254 fullModel.paramOn(
frame.get(), RooFit::Layout(0.65, 0.90, 0.90));
256 if (useRooCBShape_) {
258 crystalball.fitTo(datahist, RooFit::PrintLevel(-1));
261 RooFit::Layout(0.65, 0.90, 0.90));
264 voigt.fitTo(datahist, RooFit::PrintLevel(-1));
265 voigt.plotOn(
frame.get(), RooFit::LineColor(kRed));
266 voigt.paramOn(
frame.get(), RooFit::Layout(0.65, 0.90, 0.90));
271 datahist.plotOn(
frame.get());
272 frame->GetYaxis()->SetTitle(
"n. of events");
273 TString histName =
hist->GetName();
274 frame->SetName(
"frame" + histName);
279 c1->Print(
"fit_debug" + histName +
".pdf");
283 float mass_mean = useRooCBShape_ ? peakCB.getVal() :
mean.getVal();
284 float mass_sigma = useRooCBShape_ ? sigmaCB.getVal() : sigma.getVal();
286 float mass_mean_err = useRooCBShape_ ? peakCB.getError() :
mean.getError();
287 float mass_sigma_err = useRooCBShape_ ? sigmaCB.getError() : sigma.getError();
301 desc.add<
bool>(
"fitBackground",
false);
302 desc.add<
bool>(
"useRooCMSShape",
false);
303 desc.add<
bool>(
"useRooCBShape",
false);
304 desc.add<
bool>(
"debugMode",
false);
307 fit_par.add<std::vector<double>>(
"mean_par",
312 fit_par.add<std::vector<double>>(
"width_par",
317 fit_par.add<std::vector<double>>(
"sigma_par",
324 desc.add<std::vector<std::string>>(
"MEtoHarvest",
325 {
"DiMuMassVsMuMuPhi",
327 "DiMuMassVsMuPlusPhi",
328 "DiMuMassVsMuPlusEta",
329 "DiMuMassVsMuMinusPhi",
330 "DiMuMassVsMuMinusEta",
331 "DiMuMassVsMuMuDeltaEta",
332 "DiMuMassVsCosThetaCS"});
void getMEsToHarvest(DQMStore::IGetter &igetter)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
virtual void setCurrentFolder(std::string const &fullpath)
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_
#define DEFINE_FWK_MODULE(type)
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)
static constexpr int minimumHits
virtual MonitorElement * get(std::string const &fullpath) const
DiMuonMassBiasClient(const edm::ParameterSet &ps)
Constructor.
void bookMEs(DQMStore::IBooker &ibooker)
book MEs
double crystalball(double *x, double *par)
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())
void fillArrayF(float *x, const edm::ParameterSet &cfg, const char *name)