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.
 
 

179 lines
7.4 KiB

#include <string>
#include <iostream>
#include <cmath>
#include <algorithm>
#include <filesystem>
#include <string_view>
#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 "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 "constants.h"
#include "basic_analysis.h"
#include "hlt1_decision_analysis.h"
#include "bdt_classification.h"
// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Exclusive B2CC Bu To Kp Mu Mu %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
int new_analysis_bu2jpsikplus()
{
const char *analysis_name = "BuToJpsiKplus";
const char *data_tree_name = "Hlt2B2CC_BuToJpsiKplus_JpsiToMuMu_Detached";
const char *sim_tree_name = "BuToKpMuMu_noPID";
const char *end_state_mass_literal = "m(K^{+}#mu^{+}#mu^{-})";
const bool retrain_bdt = false;
TChain *data_chain = new TChain(TString::Format("%s/DecayTree", data_tree_name));
data_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/data_samples/b2cc_turbo_superliaison_collision23_butojpsikplus_turbopass_magdown_v0r0p6205916.root");
FourVect *l14v_data = FourVect::Init(data_chain, "muplus");
FourVect *l24v_data = FourVect::Init(data_chain, "muminus");
FourVect *hp4v_data = FourVect::Init(data_chain, "Kplus");
Bool_t muplus_ISMUON, muminus_ISMUON, Hlt1TrackMVADecision, Hlt1TwoTrackMVADecision;
data_chain->SetBranchAddress("muplus_ISMUON", &muplus_ISMUON);
data_chain->SetBranchAddress("muminus_ISMUON", &muminus_ISMUON);
data_chain->SetBranchAddress("Hlt1TrackMVADecision", &Hlt1TrackMVADecision);
data_chain->SetBranchAddress("Hlt1TwoTrackMVADecision", &Hlt1TwoTrackMVADecision);
TChain *sim_chain = new TChain(TString::Format("%s/DecayTree", sim_tree_name));
sim_chain->Add("/auto/data/pfeiffer/inclusive_detached_dilepton/MC/rd_btoxll_simulation_turbo_v0r0p6657752_BuToKpMuMu_12143001_magdown.root");
FourVect *l14v_sim = FourVect::Init(sim_chain, "L1");
FourVect *l24v_sim = FourVect::Init(sim_chain, "L2");
FourVect *hp4v_sim = FourVect::Init(sim_chain, "K");
Int_t B_BKGCAT, L1_TRUEID, L2_TRUEID, Hp_TRUEID;
sim_chain->SetBranchAddress("L1_TRUEID", &L1_TRUEID);
sim_chain->SetBranchAddress("L2_TRUEID", &L2_TRUEID);
sim_chain->SetBranchAddress("K_TRUEID", &Hp_TRUEID);
sim_chain->SetBranchAddress("B_BKGCAT", &B_BKGCAT);
Double_t B_Mass_jpsi_var, B_Mass_psi2s_var, B_Mass_sim_var;
TString B_Mass_jpsi_var_name = "B_Mass_jpsi_var";
TString B_Mass_psi2s_var_name = "B_Mass_psi2s_var";
TString B_Mass_sim_var_name = "B_Mass_sim_var";
TTree *tree_B_Mass_jpsi = new TTree("tree_B_Mass_jpsi", TString::Format("B Mass, J/#psi Mode (%s)", analysis_name));
TTree *tree_B_Mass_psi2s = new TTree("tree_B_Mass_psi2s", TString::Format("B Mass, #psi(2S) Mode (%s)", analysis_name));
TTree *tree_B_Mass_sim = new TTree("tree_B_Mass_sim", TString::Format("B Mass, Simualted (%s_noPID)", analysis_name));
tree_B_Mass_jpsi->Branch(B_Mass_jpsi_var_name, &B_Mass_jpsi_var);
tree_B_Mass_psi2s->Branch(B_Mass_psi2s_var_name, &B_Mass_psi2s_var);
tree_B_Mass_sim->Branch(B_Mass_sim_var_name, &B_Mass_sim_var);
TH1D *h1_B_Mass_unf = new TH1D("h1_B_Mass_unf", TString::Format("B Mass (%s), Unfiltered", data_tree_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
TH1D *h1_B_Mass_sim_unf = new TH1D("h1_B_Mass_sim_unf", TString::Format("B Mass, Simualted (%s), Unfiltered", sim_tree_name), B_MASS_HIST_BINS, B_MASS_VAR_MIN, B_MASS_VAR_MAX);
TH2D *h2_Hlt1_flags_B_Mass = new TH2D("h2_Hlt1_flags_B_Mass", "Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
TH2D *h2_Hlt1_flags_excl_B_Mass = new TH2D("h2_Hlt1_flags_excl_B_Mass", "Excl Hlt1 Decision vs B Mass", 50, 5100, 5400, 13, 1., 14.);
h1_B_Mass_unf->GetXaxis()->SetTitle(end_state_mass_literal);
h1_B_Mass_sim_unf->GetXaxis()->SetTitle(end_state_mass_literal);
h2_Hlt1_flags_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
h2_Hlt1_flags_excl_B_Mass->GetXaxis()->SetTitle(end_state_mass_literal);
// ConnectHlt1Decisions(data_chain, h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass);
// auto hlt1_decision_histos = CreateHlt1DecisionHistos(analysis_name);
// std::map<std::string, int> exclusive_hits{};
unsigned int data_entries = data_chain->GetEntries();
unsigned int sim_entries = sim_chain->GetEntries();
for (unsigned int i = 0; i < sim_entries; i++)
{
sim_chain->GetEntry(i);
Double_t reconstructed_B_Mass = (hp4v_sim->LorentzVector() + l14v_sim->LorentzVector() + l24v_sim->LorentzVector()).M();
h1_B_Mass_sim_unf->Fill(reconstructed_B_Mass);
if (B_BKGCAT == 0 && TMath::Abs(L1_TRUEID) == PID_MUON && L2_TRUEID == -L1_TRUEID && TMath::Abs(Hp_TRUEID) == PID_KAON)
{
B_Mass_sim_var = reconstructed_B_Mass;
tree_B_Mass_sim->Fill();
}
PrintProgress(TString::Format("%s SIG Collection", analysis_name), sim_entries, 10000, i);
}
std::cout << "# Starting evaluation of data." << std::endl;
for (unsigned int i = 0; i < data_entries; i++)
{
data_chain->GetEntry(i);
TLorentzVector dimuon = l14v_data->LorentzVector() + l24v_data->LorentzVector();
Double_t reconstructed_B_Mass = (hp4v_data->LorentzVector() + dimuon).M();
// if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
// {
// CheckHlt1Decisioins(h2_Hlt1_flags_B_Mass, h2_Hlt1_flags_excl_B_Mass, exclusive_hits, reconstructed_B_Mass);
// FillHlt1DecisionHistos(hlt1_decision_histos, reconstructed_B_Mass);
// }
h1_B_Mass_unf->Fill(reconstructed_B_Mass);
if (muplus_ISMUON && muminus_ISMUON && (Hlt1TrackMVADecision || Hlt1TwoTrackMVADecision))
{
if (TMath::Abs(dimuon.M() - JPSI_MASS) < 100.)
{
B_Mass_jpsi_var = reconstructed_B_Mass;
tree_B_Mass_jpsi->Fill();
}
else if (TMath::Abs(dimuon.M() - PSI2S_MASS) < 100.)
{
B_Mass_psi2s_var = reconstructed_B_Mass;
tree_B_Mass_psi2s->Fill();
}
}
PrintProgress(TString::Format("%s BDT Evaluation", analysis_name), data_entries, 10000, i);
}
// std::cout << "# Exclusive Hits" << std::endl;
// for (const auto &[line, hits] : exclusive_hits)
// {
// std::cout << line << ": " << hits << std::endl;
// }
DrawInDefaultCanvas(h1_B_Mass_unf, analysis_name, 0.1);
DrawInDefaultCanvas(h1_B_Mass_sim_unf, analysis_name, 0.1);
auto roofit_hist_sim = CreateRooDataSetAndFitCB(tree_B_Mass_sim, B_Mass_sim_var_name, end_state_mass_literal, false, false, ShapeParamters{});
auto roofit_hist_jpsi_fitsum = CreateRooDataSetAndFitCB(tree_B_Mass_jpsi, B_Mass_jpsi_var_name, end_state_mass_literal, true, true, roofit_hist_sim.shape_parameters, 5140., 5460.);
auto roofit_hist_psi2s_fitsum = CreateRooDataSetAndFitCB(tree_B_Mass_psi2s, B_Mass_psi2s_var_name, end_state_mass_literal, true, true, roofit_hist_sim.shape_parameters, 5140., 5460.);
DrawInDefaultCanvas(roofit_hist_jpsi_fitsum, analysis_name);
DrawInDefaultCanvas(roofit_hist_psi2s_fitsum, analysis_name);
DrawInDefaultCanvas(roofit_hist_sim, analysis_name);
// DrawHlt1DecisionHistos(analysis_name, hlt1_decision_histos);
return 0;
}