ROOT Analysis for the Inclusive Detachted Dilepton Trigger Lines
You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
 
 

248 lines
10 KiB

#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 "TLorentzVector.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>
#include <filesystem>
const int nBins = 70;
const Double_t MASS_HIST_MIN = 5100.;
const Double_t MASS_HIST_MAX = 6000.;
const Double_t K_MASS = 493.677;
const int PID_KAON = 321;
const int PID_PION = 211;
const int PID_MUON = 13;
const char *NAME = "BuToHpMuMu_Fullstream_NoKTrueID";
const char *X_AXIS = "m(#pi^{+}_{#rightarrow K^{+}}#mu^{+}#mu^{-})";
const char *B_NAME = "B";
TString quickformat_name(const char *format)
{
return TString::Format(format, NAME);
};
void CreateRooFitAndDraw(TH1D *hist, double total_entries, double selected_entries);
int simulation_fullstream_Bu2KpMuMu()
{
TChain *data_chain = new TChain("BuToHpMuMu_noPID/DecayTree");
data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_fullstream_v0r0p6671378_BuToKpMuMu_12143001_magdown.root");
/*
rd_btoxll_simulation_fullstream_v0r0p6671378_B0ToKpPimMuMu_11144002_magdown.root
rd_btoxll_simulation_fullstream_v0r0p6671378_BuToKpMuMu_12143001_magdown.root
rd_btoxll_simulation_turbo_v0r0p6657752_B0ToKpPimMuMu_11144002_magdown.root
rd_btoxll_simulation_turbo_v0r0p6657752_BuToKpMuMu_12143001_magdown.root
*/
Int_t B_BKGCAT;
data_chain->SetBranchAddress(TString::Format("%s_BKGCAT", B_NAME).Data(), &B_BKGCAT);
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;
Double_t B_M;
Int_t L1_TRUEID, L2_TRUEID, Hp_TRUEID;
data_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
data_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
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("Hp_TRUEID", &Hp_TRUEID);
TH1D *h1_B_M = new TH1D("h1_B_M", "B Mass", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
TH1D *h1_B_M_trueid = new TH1D("h1_B_M_trueid", "B Mass, TrueID matched", nBins, MASS_HIST_MIN, MASS_HIST_MAX);
TH1I *h1_B_BKGCAT = new TH1I("h1_B_BKGCAT", "B Background Category", nBins, 0, 120);
TH2D *h2_B_M_vs_B_BKGCAT = new TH2D("h2_B_M_vs_B_BKGCAT", "B Mass vs Background Category", nBins, MASS_HIST_MIN, MASS_HIST_MAX, nBins, 0, 120);
h1_B_M->GetXaxis()->SetTitle(X_AXIS);
h1_B_M_trueid->GetXaxis()->SetTitle(X_AXIS);
unsigned int entries = data_chain->GetEntries();
unsigned int selected_entries = 0;
for (size_t i = 0; i < entries; i++)
{
data_chain->GetEntry(i);
Double_t reconstructed_B_Mass = 0;
if (B_BKGCAT == 30 && TMath::Abs(L1_TRUEID) == PID_MUON && L2_TRUEID == -L1_TRUEID && TMath::Abs(Hp_TRUEID) == PID_KAON)
{
TVector3 K_momentum(Hp_PX, Hp_PY, Hp_PZ);
double K_energy = TMath::Sqrt(TMath::Sq(K_MASS) + K_momentum.Mag2());
TLorentzVector K_4v(K_momentum, K_energy);
TLorentzVector l1_4v(L1_PX, L1_PY, L1_PZ, L1_ENERGY);
TLorentzVector l2_4v(L2_PX, L2_PY, L2_PZ, L2_ENERGY);
reconstructed_B_Mass = (K_4v + l1_4v + l2_4v).M();
selected_entries++;
h1_B_M_trueid->Fill(reconstructed_B_Mass);
h1_B_M->Fill(reconstructed_B_Mass);
h1_B_BKGCAT->Fill(B_BKGCAT);
h2_B_M_vs_B_BKGCAT->Fill(reconstructed_B_Mass, B_BKGCAT);
}
if (i % 10000 == 0)
{
std::cout << "["
<< NAME
<< "] Processed event: " << i << " (" << TString::Format("%.2f", ((double)i / (double)entries) * 100.) << "%), "
<< "Current Efficiency: " << TString::Format("%.2f%%", ((double)selected_entries / (double)entries) * 100.)
<< std::endl;
}
}
std::cout << "["
<< NAME
<< "] Processed event: " << entries << " (" << TString::Format("%.2f", 100.) << "%)" << std::endl;
h1_B_M->SetLineColor(kBlue);
h1_B_M_trueid->SetLineColor(kMagenta + 2);
h1_B_M_trueid->SetFillColor(kMagenta + 2);
h1_B_M_trueid->SetFillStyle(3244);
h1_B_BKGCAT->SetLineColor(kBlue);
h1_B_BKGCAT->SetFillColor(kBlue);
h1_B_BKGCAT->SetFillStyle(3351);
std::filesystem::create_directory(TString::Format("images/sim/%s", NAME).Data());
TCanvas *c1 = new TCanvas("c1", "c1", 0, 0, 800, 600);
h1_B_M->SetStats(0);
h1_B_M_trueid->SetStats(0);
h1_B_M->Draw();
h1_B_M_trueid->Draw("SAME");
c1->BuildLegend(0.58, 0.65, 0.85, 0.87);
c1->Draw();
c1->SaveAs(quickformat_name("images/sim/%s/h1_B_M.pdf"));
TCanvas *c2 = new TCanvas("c2", "c2", 0, 0, 800, 600);
h1_B_BKGCAT->SetStats(0);
h1_B_BKGCAT->Draw();
c2->BuildLegend(0.58, 0.77, 0.85, 0.87);
c2->Draw();
c2->SaveAs(quickformat_name("images/sim/%s/h1_B_BKGCAT.pdf"));
TCanvas *c3 = new TCanvas("c3", "c3", 0, 0, 800, 600);
h2_B_M_vs_B_BKGCAT->SetStats(0);
h2_B_M_vs_B_BKGCAT->Draw("COLZ");
c3->Draw();
c3->SaveAs(quickformat_name("images/sim/%s/h2_B_M_vs_B_BKGCAT.pdf"));
CreateRooFitAndDraw(h1_B_M_trueid, entries, selected_entries);
return 0;
}
void CreateRooFitAndDraw(TH1D *hist, double total_entries, double selected_entries)
{
// std::cout << "### Entries: " << hist->GetEntries() << std::endl;
RooRealVar roo_var_mass("var_mass", "B Mass Variable", MASS_HIST_MIN, MASS_HIST_MAX);
RooDataHist roohist_B_M("roohist_B_M", "B Mass Histogram", roo_var_mass, RooFit::Import(*hist));
RooPlot *roo_frame_mass = roo_var_mass.frame(RooFit::Title(quickformat_name("Simulation of %s")));
roohist_B_M.plotOn(roo_frame_mass, RooFit::Binning(nBins), RooFit::Name("B Mass Distribution"));
roo_frame_mass->GetXaxis()->SetTitle(X_AXIS);
// 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., 25.);
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.);
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", "Background C", -0.00231049, -0.003, -0.002);
// RooExponential bkg_exp("bkg_exp", "Exp Background", var_mass, var_mass_bkg_c);
// RooRealVar var_mass_nsig(TString::Format("var_%s_nsig", var_id.c_str()), "Mass N Signal", 0., hist->GetEntries());
// RooRealVar var_mass_nbkg(TString::Format("var_%s_nbkg", var_id.c_str()), "Mass N Background", 0., hist->GetEntries());
// TString pdf_name = TString::Format("%s_sigplusbkg", var_id.c_str());
// RooAddPdf sigplusbkg(pdf_name, "Sig and Bkg PDF", RooArgList(sig_cb, bkg_exp), RooArgList(var_mass_nsig, var_mass_nbkg));
RooFitResult *fitres = sig_cb.fitTo(roohist_B_M, RooFit::Save(), RooFit::PrintLevel(1), RooFit::Range("fitting_range"));
sig_cb.plotOn(roo_frame_mass, RooFit::VisualizeError(*fitres, 1), RooFit::FillColor(kOrange + 1), RooFit::FillStyle(3144));
sig_cb.plotOn(roo_frame_mass, RooFit::LineColor(kRed), RooFit::LineStyle(kSolid), RooFit::Range("fitting_range"), RooFit::Name("sig_cb"));
// sigplusbkg.plotOn(roo_frame_mass, RooFit::Components(RooArgSet(bkg_exp)), RooFit::LineColor(kBlue), RooFit::LineStyle(kDashed), RooFit::Range(5170., 5500.));
TCanvas *c = new TCanvas("roofit_c", "roofit_c", 0, 0, 800, 600);
roo_frame_mass->Draw();
TLegend *leg1 = new TLegend(0.58, 0.50, 0.87, 0.87);
// leg1->SetFillColor(kWhite);
leg1->SetLineColor(kWhite);
// leg1->AddEntry(roo_frame_mass->findObject(pdf_name.c_str()), "Signal + Background", "LP");
leg1->AddEntry(roo_frame_mass->findObject("sig_cb"), "Signal", "LP");
// leg1->AddEntry(roo_frame_mass->findObject(name_fit_func_bkg.c_str()), "Background", "LP");
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 = %.8f #pm %.8f", roo_sig_tail.getTitle().Data(), roo_sig_tail.getVal(), roo_sig_tail.getError()).Data(), "");
leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_alphaL.getTitle().Data(), var_mass_alphaL.getVal(), var_mass_alphaL.getError()).Data(), "");
leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_nL.getTitle().Data(), var_mass_nL.getVal(), var_mass_nL.getError()).Data(), "");
leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_alphaR.getTitle().Data(), var_mass_alphaR.getVal(), var_mass_alphaR.getError()).Data(), "");
leg1->AddEntry((TObject *)0, TString::Format("%s = %.3f #pm %.3f", var_mass_nR.getTitle().Data(), var_mass_nR.getVal(), var_mass_nR.getError()).Data(), "");
leg1->AddEntry((TObject *)0, TString::Format("#epsilon = %.2f%%", (selected_entries / total_entries) * 100).Data(), "");
leg1->Draw();
c->Draw();
c->SaveAs(quickformat_name("images/sim/%s/roof_B_M.pdf"));
}