|
|
#include "TH1D.h"
#include "TH2D.h"
#include "THStack.h"
#include "TGraph.h"
#include "TTree.h"
#include "TChain.h"
#include "TFile.h"
#include "TCanvas.h"
#include "TROOT.h"
#include "TStyle.h"
#include "TColor.h"
#include "TLorentzVector.h"
#include "TRandom3.h"
#include "TLegend.h"
#include "TLatex.h"
#include "RooDataHist.h"
#include "RooRealVar.h"
#include "RooPlot.h"
#include "RooGaussian.h"
#include "RooExponential.h"
#include "RooRealConstant.h"
#include "RooAddPdf.h"
#include "RooFitResult.h"
#include "RooProduct.h"
#include "RooCrystalBall.h"
#include "RooBreitWigner.h"
#include "RooArgSet.h"
#include "RooFFTConvPdf.h"
#include "RooNovosibirsk.h"
#include <string>
#include <iostream>
#include <cmath>
const int nBins = 70;
const Double_t MASS_HIST_MIN = 5100.; const Double_t MASS_HIST_MAX = 5700.; const Double_t MASS_HIST_FIT_MIN = 5100.; const Double_t MASS_HIST_FIT_MAX = 5700.;
const char *TITLE = "SpruceRD_B0ToHpHmMuMu (#pi^{-} #rightarrow K^{-})"; const char *FILE_NAME = "SpruceRD_B0ToHpHmMuMu_Pip2Kp"; const char *MASS_LITERAL = "m(#pi^{+}#pi^{-}_{(#rightarrow K^{-})}#mu^{+}#mu^{-})";
const Double_t K_MASS = 493.677; const Double_t PSI2S_MASS = 3686.093; const Double_t JPSI_MASS = 3096.9; const Double_t KSTAR_MASS = 891.76;
void CreateRooFitAndDraw(TH1D *hist);
int analysis_fullstream_b02hphmmumu() { TChain *data_chain = new TChain("B0ToHpHmMuMu/DecayTree"); data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0r1_tuple_90000000_2023_v0r0p6288631.root"); data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/spruce_magdown_2023_v0_tuple_90000000_v0r0p6288631.root");
Float_t B0_BPVFDCHI2, B0_BPVIPCHI2, L1_BPVIPCHI2, L2_BPVIPCHI2, L1_PT, L2_PT, Jpsi_BPVFDCHI2, Hp_PT, Hp_BPVIPCHI2, Hp_P, Hm_PT, Hm_BPVIPCHI2, Hm_P, Res_PT, Res_BPVIPCHI2, Res_P;
Double_t L1_PID_MU, L2_PID_MU, B0_CHI2VXNDOF, Jpsi_MAXDOCACHI2, Jpsi_CHI2DOF, Res_CHI2DOF, Hp_PID_K, Hm_PID_K, Jpsi_M, B0_M, Res_MAXDOCACHI2, Res_M;
Bool_t L1_ISMUON, L2_ISMUON, Hlt2RD_B0ToKpPimMuMuDecision, Hlt2_InclDetDiMuon_4BodyDecision, Hlt2_InclDetDiMuon_3BodyDecision, Hlt2_InclDetDiMuonDecision, Hlt1TrackMVADecision, Hlt1TwoTrackMVADecision;
data_chain->SetBranchAddress("B0_M", &B0_M);
data_chain->SetBranchAddress("B0_BPVFDCHI2", &B0_BPVFDCHI2);
data_chain->SetBranchAddress("B0_BPVIPCHI2", &B0_BPVIPCHI2);
data_chain->SetBranchAddress("L1_BPVIPCHI2", &L1_BPVIPCHI2); data_chain->SetBranchAddress("L2_BPVIPCHI2", &L2_BPVIPCHI2); data_chain->SetBranchAddress("L1_PID_MU", &L1_PID_MU); data_chain->SetBranchAddress("L2_PID_MU", &L2_PID_MU); data_chain->SetBranchAddress("L1_ISMUON", &L1_ISMUON); data_chain->SetBranchAddress("L2_ISMUON", &L2_ISMUON); data_chain->SetBranchAddress("L1_PT", &L1_PT); data_chain->SetBranchAddress("L2_PT", &L2_PT);
data_chain->SetBranchAddress("B0_CHI2VXNDOF", &B0_CHI2VXNDOF); data_chain->SetBranchAddress("Jpsi_MAXDOCACHI2", &Jpsi_MAXDOCACHI2); data_chain->SetBranchAddress("Jpsi_CHI2DOF", &Jpsi_CHI2DOF); data_chain->SetBranchAddress("Res_CHI2DOF", &Res_CHI2DOF);
data_chain->SetBranchAddress("Hp_PT", &Hp_PT); data_chain->SetBranchAddress("Hp_BPVIPCHI2", &Hp_BPVIPCHI2); data_chain->SetBranchAddress("Hp_P", &Hp_P); data_chain->SetBranchAddress("Hp_PID_K", &Hp_PID_K); data_chain->SetBranchAddress("Hm_PT", &Hm_PT); data_chain->SetBranchAddress("Hm_BPVIPCHI2", &Hm_BPVIPCHI2); data_chain->SetBranchAddress("Hm_P", &Hm_P); data_chain->SetBranchAddress("Hm_PID_K", &Hm_PID_K); data_chain->SetBranchAddress("Res_PT", &Res_PT); data_chain->SetBranchAddress("Res_BPVIPCHI2", &Res_BPVIPCHI2); data_chain->SetBranchAddress("Res_P", &Res_P); data_chain->SetBranchAddress("Res_MAXDOCACHI2", &Res_MAXDOCACHI2); data_chain->SetBranchAddress("Res_M", &Res_M);
data_chain->SetBranchAddress("Jpsi_BPVFDCHI2", &Jpsi_BPVFDCHI2); data_chain->SetBranchAddress("Jpsi_M", &Jpsi_M);
data_chain->SetBranchAddress("Hlt2RD_B0ToKpPimMuMuDecision", &Hlt2RD_B0ToKpPimMuMuDecision); data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_4BodyDecision", &Hlt2_InclDetDiMuon_4BodyDecision); data_chain->SetBranchAddress("Hlt2_InclDetDiMuon_3BodyDecision", &Hlt2_InclDetDiMuon_3BodyDecision); data_chain->SetBranchAddress("Hlt2_InclDetDiMuonDecision", &Hlt2_InclDetDiMuonDecision); data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision); data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
Float_t L1_PX, L1_PY, L1_PZ, L1_ENERGY, L2_PX, L2_PY, L2_PZ, L2_ENERGY, Hp_PX, Hp_PY, Hp_PZ, Hp_ENERGY, Hm_PX, Hm_PY, Hm_PZ, Hm_ENERGY;
data_chain->SetBranchAddress("L1_PX", &L1_PX); data_chain->SetBranchAddress("L1_PY", &L1_PY); data_chain->SetBranchAddress("L1_PZ", &L1_PZ); data_chain->SetBranchAddress("L1_ENERGY", &L1_ENERGY); data_chain->SetBranchAddress("L2_PX", &L2_PX); data_chain->SetBranchAddress("L2_PY", &L2_PY); data_chain->SetBranchAddress("L2_PZ", &L2_PZ); data_chain->SetBranchAddress("L2_ENERGY", &L2_ENERGY); data_chain->SetBranchAddress("Hp_PX", &Hp_PX); data_chain->SetBranchAddress("Hp_PY", &Hp_PY); data_chain->SetBranchAddress("Hp_PZ", &Hp_PZ); data_chain->SetBranchAddress("Hp_ENERGY", &Hp_ENERGY); data_chain->SetBranchAddress("Hm_PX", &Hm_PX); data_chain->SetBranchAddress("Hm_PY", &Hm_PY); data_chain->SetBranchAddress("Hm_PZ", &Hm_PZ); data_chain->SetBranchAddress("Hm_ENERGY", &Hm_ENERGY);
TH1D *h1_B0_M = new TH1D("h1_B0_M", TITLE, nBins, MASS_HIST_MIN, MASS_HIST_MAX); TH1D *h1_Res_M = new TH1D("h1_Res_M", "K^{*0} Mass", nBins, 0., 3500.); TH1D *h1_Hm_PID_K = new TH1D("h1_Hm_PID_K", "H^{-} PID K", nBins, -10., 10.); TH1D *h1_Hp_PID_K = new TH1D("h1_Hp_PID_K", "H^{+} PID K", nBins, -10., 10.);
unsigned int entries = data_chain->GetEntries(); unsigned int fitting_entries = 0; for (unsigned int i = 0; i < entries; i++) { data_chain->GetEntry(i);
TLorentzVector l1_4v(L1_PX, L1_PY, L1_PZ, L1_ENERGY); TLorentzVector l2_4v(L2_PX, L2_PY, L2_PZ, L2_ENERGY);
// Pi+ -> K+
TVector3 Kp_momentum(Hp_PX, Hp_PY, Hp_PZ); double Kp_energy = TMath::Sqrt(TMath::Sq(K_MASS) + Kp_momentum.Mag2()); TLorentzVector Pim_4v(Hm_PX, Hm_PY, Hm_PZ, Hm_ENERGY); TLorentzVector Kp_4v(Kp_momentum, Kp_energy); Double_t reconstructed_Kp_Pim_Mass = (Kp_4v + Pim_4v).M(); Double_t reconstructed_B0_Kp_Pim_Mass = (Kp_4v + Pim_4v + l1_4v + l2_4v).M();
// Pi- -> K-
TVector3 Km_momentum(Hm_PX, Hm_PY, Hm_PZ); double Km_energy = TMath::Sqrt(TMath::Sq(K_MASS) + Km_momentum.Mag2()); TLorentzVector Pip_4v(Hp_PX, Hp_PY, Hp_PZ, Hp_ENERGY); TLorentzVector Km_4v(Km_momentum, Km_energy); Double_t reconstructed_Km_Pip_Mass = (Km_4v + Pip_4v).M(); Double_t reconstructed_B0_Km_Pip_Mass = (Km_4v + Pip_4v + l1_4v + l2_4v).M();
h1_Res_M->Fill(reconstructed_Kp_Pim_Mass);
h1_Hm_PID_K->Fill(Hm_PID_K); h1_Hp_PID_K->Fill(Hp_PID_K);
if (( (B0_BPVFDCHI2 > 64) & (B0_CHI2VXNDOF < 9) & (B0_BPVIPCHI2 < 25) & //
(Jpsi_MAXDOCACHI2 < 36) & //
(L1_BPVIPCHI2 > 9) & (L2_BPVIPCHI2 > 9) & (L1_PID_MU > -3) & (L2_PID_MU > -3) & (L1_ISMUON) & (L2_ISMUON) & (L1_PT > 350) & (L2_PT > 350) & //
(Jpsi_CHI2DOF < 9) & (Jpsi_BPVFDCHI2 > 16) & (Jpsi_M < 5500) & //
(Res_PT > 400) & (Res_MAXDOCACHI2 < 36) & (reconstructed_Kp_Pim_Mass < 2600) & (reconstructed_Km_Pip_Mass < 2600) & (Res_CHI2DOF < 25) & (Res_BPVIPCHI2 > 4) & //
(Hm_BPVIPCHI2 > 6) & (Hm_PT > 250) & (Hm_P > 2000) & (Hm_PID_K > 0) & //
(Hp_BPVIPCHI2 > 6) & (Hp_PT > 250) & (Hp_P > 2000) & (Hp_PID_K < 0) ) & //
(Hlt2RD_B0ToKpPimMuMuDecision) & (((Hlt2_InclDetDiMuon_4BodyDecision) | (Hlt2_InclDetDiMuon_3BodyDecision) | (Hlt2_InclDetDiMuonDecision)))) { if ((TMath::Abs(Jpsi_M - JPSI_MASS) < 100) & ((Hlt1TrackMVADecision) | (Hlt1TwoTrackMVADecision))) { h1_B0_M->Fill(reconstructed_B0_Km_Pip_Mass); } }
if ((i + 1) % 10000 == 0 || i + 1 == entries) { std::cout << "[" << FILE_NAME << "] Processed event: " << i + 1 << " (" << TString::Format("%.2f", ((double)(i + 1) / (double)entries) * 100.) << "%)" << std::endl; } }
h1_B0_M->GetXaxis()->SetTitle(MASS_LITERAL);
TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600); h1_B0_M->Draw(); c1->Draw(); c1->SaveAs(TString::Format("images/root_hist_%s_bmass.pdf", FILE_NAME).Data());
TCanvas *c2 = new TCanvas("c2", "c2", 0, 0, 800, 600); h1_Res_M->Draw(); c2->Draw();
TCanvas *c3 = new TCanvas("c3", "c3", 0, 0, 800, 600);
h1_Hm_PID_K->SetStats(0); h1_Hm_PID_K->SetLineColor(kRed);
h1_Hp_PID_K->SetStats(0);
h1_Hm_PID_K->Draw(); h1_Hp_PID_K->Draw("SAME");
c3->BuildLegend(); c3->Draw();
CreateRooFitAndDraw(h1_B0_M);
return 0; }
void CreateRooFitAndDraw(TH1D *hist) { RooRealVar roo_var_mass("var_mass", "B0 Mass Variable", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX); roo_var_mass.setRange("fitting_range", MASS_HIST_FIT_MIN, MASS_HIST_FIT_MAX);
RooDataHist roohist_B0_M("roohist_B0_M", "B0 Mass Histogram", roo_var_mass, RooFit::Import(*hist));
RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(TITLE)); roohist_B0_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name("B Mass Distribution")); roo_frame_mass->GetXaxis()->SetTitle(MASS_LITERAL);
// Crystal Ball for Signal
RooRealVar var_mass_x0("var_mass_x0", "#mu", 5278., 5170., 5500.); RooRealVar var_mass_sigmaLR("var_mass_sigmaLR", "#sigma_{LR}", 16., 5., 40.);
// Same Variables for Left and Right Tail
// RooRealVar var_mass_alphaL("var_mass_alphaL", "#alpha_{L}", 2., 0., 4.);
// RooRealVar var_mass_nL("var_mass_nL", "n_{L}", 5., 0., 15.);
// RooRealVar var_mass_alphaR("var_mass_alphaR", "#alpha_{R}", 2., 0., 4.);
// RooRealVar var_mass_nR("var_mass_nR", "n_{R}", 5., 0., 15.);
auto var_mass_alphaL = RooRealConstant::value(1.915); auto var_mass_nL = RooRealConstant::value(0.573); auto var_mass_alphaR = RooRealConstant::value(2.257); auto var_mass_nR = RooRealConstant::value(0.418);
RooCrystalBall sig_cb("sig_cb", "Signal Crystal Ball", roo_var_mass, var_mass_x0, var_mass_sigmaLR, var_mass_alphaL, var_mass_nL, var_mass_alphaR, var_mass_nR);
// Exponential for Background
RooRealVar var_mass_bkg_c("var_mass_bkg_c", "#lambda", -0.0014, -0.004, -0.000); RooExponential bkg_exp("bkg_exp", "Exp Background", roo_var_mass, var_mass_bkg_c);
RooRealVar var_mass_nsig("nsig", "Mass N Signal", 0., hist->GetEntries()); RooRealVar var_mass_nbkg("nbkg", "Mass N Background", 0., hist->GetEntries());
// TString pdf_name = TString::Format("%s_sigplusbkg", var_id.c_str());
RooAddPdf sigplusbkg("sigplusbkg", "Sig and Bkg PDF", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
RooFitResult *fitres = sigplusbkg.fitTo(roohist_B0_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range("fitting_range"));
sigplusbkg.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144)); sigplusbkg.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range("fitting_range"), RooFit::Name("sigplusbkg")); sigplusbkg.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(bkg_exp)), RooFit::LineColor(kBlue - 7), RooFit::LineStyle(kDashed), RooFit::Range("fitting_range"), RooFit::Name("bkg_exp")); sigplusbkg.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(sig_cb)), RooFit::LineColor(kRed - 7), RooFit::LineStyle(kDashed), RooFit::Range("fitting_range"), RooFit::Name("sig_cb"));
TCanvas *c = new TCanvas("roofit_c", "roofit_c", 0, 0, 800, 600);
roo_frame_mass->Draw();
TLegend *leg1 = new TLegend(0.50, 0.58, 0.87, 0.89); leg1->SetFillColor(kWhite); leg1->SetLineColor(kBlack);
leg1->AddEntry(roo_frame_mass->findObject("sigplusbkg"), "Signal + Background", "LP"); leg1->AddEntry(roo_frame_mass->findObject("sig_cb"), "Signal", "LP"); leg1->AddEntry(roo_frame_mass->findObject("bkg_exp"), "Background", "LP"); leg1->AddEntry((TObject *)0, "", ""); leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_x0.getTitle().Data(), var_mass_x0.getVal(), var_mass_x0.getError()).Data(), ""); leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_sigmaLR.getTitle().Data(), var_mass_sigmaLR.getVal(), var_mass_sigmaLR.getError()).Data(), ""); leg1->AddEntry((TObject *)0, TString::Format("%s = %.6f #pm %.6f", var_mass_bkg_c.getTitle().Data(), var_mass_bkg_c.getVal(), var_mass_bkg_c.getError()).Data(), ""); leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "S", var_mass_nsig.getVal(), var_mass_nsig.getError()).Data(), ""); leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", "B", var_mass_nbkg.getVal(), var_mass_nbkg.getError()).Data(), "");
leg1->AddEntry((TObject *)0, TString::Format("Entries: %.1f", var_mass_nsig.getVal() + var_mass_nbkg.getVal()), "");
leg1->Draw();
c->Draw(); c->SaveAs(TString::Format("images/data/roofit_hist_%s_fitted_bmass.pdf", FILE_NAME).Data()); }
|