File indexing completed on 2024-04-06 11:57:25
0001 #include <TH2.h>
0002 #include <TStyle.h>
0003 #include <TCanvas.h>
0004 #include <TMath.h>
0005 #include <cmath>
0006 #include <TString.h>
0007
0008 void Get_Plot(TH1D, TString);
0009
0010 void CosmicRateTool_Kinematics(const char *fileName) {
0011 TString InputFile = Form("%s", fileName);
0012 TFile *file = new TFile(InputFile);
0013
0014 bool IsFileExist;
0015 IsFileExist = file->IsZombie();
0016 if (IsFileExist) {
0017 cout << endl
0018 << "====================================================================================================="
0019 << endl;
0020 cout << fileName << " is not found. Check the file!" << endl;
0021 cout << "====================================================================================================="
0022 << endl
0023 << endl;
0024 exit(EXIT_FAILURE);
0025 }
0026
0027 TTree *tree = (TTree *)file->Get("cosmicRateAnalyzer/Event");
0028
0029 vector<double> *pt;
0030 vector<double> *charge;
0031 vector<double> *chi2;
0032 vector<double> *chi2_ndof;
0033 vector<double> *eta;
0034 vector<double> *theta;
0035 vector<double> *phi;
0036 vector<double> *p;
0037 vector<double> *d0;
0038 vector<double> *dz;
0039
0040 pt = 0;
0041 charge = 0;
0042 chi2 = 0;
0043 chi2_ndof = 0;
0044 eta = 0;
0045 theta = 0;
0046 phi = 0;
0047 p = 0;
0048 d0 = 0;
0049 dz = 0;
0050
0051 tree->SetBranchAddress("pt", &pt);
0052 tree->SetBranchAddress("charge", &charge);
0053 tree->SetBranchAddress("chi2", &chi2);
0054 tree->SetBranchAddress("chi2_ndof", &chi2_ndof);
0055 tree->SetBranchAddress("eta", &eta);
0056 tree->SetBranchAddress("theta", &theta);
0057 tree->SetBranchAddress("phi", &phi);
0058 tree->SetBranchAddress("p", &p);
0059 tree->SetBranchAddress("d0", &d0);
0060 tree->SetBranchAddress("dz", &dz);
0061
0062
0063
0064
0065
0066 TH1D h_pt("h_pt", "h_pt", 100, 0, 100);
0067 TH1D h_charge("h_charge", "h_charge", 10, -5, 5);
0068 TH1D h_chi2("h_chi2", "h_chi2", 200, 0, 100);
0069 TH1D h_chi2_ndof("h_chi2_ndof", "h_chi2_ndof", 100, 0, 10);
0070 TH1D h_eta("h_eta", "h_eta", 500, -3, 3);
0071 TH1D h_theta("h_theta", "h_theta", 500, -3, 3);
0072 TH1D h_phi("h_phi", "h_phi", 400, -3.5, 3.5);
0073 TH1D h_d0("h_d0", "h_d0", 1000, -85, 85);
0074 TH1D h_dz("h_dz", "h_dz", 1500, -350, 350);
0075
0076
0077
0078 int nTotalEvents = 0, nTotalTracks = 0;
0079 Long64_t n = tree->GetEntriesFast();
0080 for (Long64_t jentry = 0; jentry < n; jentry++)
0081 {
0082 tree->GetEntry(jentry);
0083
0084 for (int k = 0; k < pt->size(); k++)
0085 {
0086 h_pt.Fill(pt->at(k));
0087 h_charge.Fill(charge->at(k));
0088 h_chi2.Fill(chi2->at(k));
0089 h_chi2_ndof.Fill(chi2_ndof->at(k));
0090 h_eta.Fill(eta->at(k));
0091 h_theta.Fill(theta->at(k));
0092 h_phi.Fill(phi->at(k));
0093 h_d0.Fill(d0->at(k));
0094 h_dz.Fill(dz->at(k));
0095
0096 nTotalTracks++;
0097 }
0098
0099 nTotalEvents++;
0100 }
0101
0102 std::cout << "Total Events: " << nTotalEvents << std::endl;
0103 std::cout << "Total Tracks: " << nTotalTracks << std::endl;
0104
0105
0106 gSystem->Exec("mkdir -p Kinematical_Plots");
0107
0108
0109 Get_Plot(h_pt, "pt");
0110 Get_Plot(h_eta, "eta");
0111 Get_Plot(h_phi, "phi");
0112 Get_Plot(h_theta, "theta");
0113 Get_Plot(h_d0, "d0");
0114 Get_Plot(h_dz, "dz");
0115 Get_Plot(h_chi2, "chi2");
0116 Get_Plot(h_chi2_ndof, "chi2_ndof");
0117 Get_Plot(h_charge, "charge");
0118 }
0119
0120 void Get_Plot(TH1D h1, TString variable) {
0121 TCanvas c("c", "c", 556, 214, 661, 641);
0122 gStyle->SetOptStat(0);
0123 gStyle->SetOptTitle(0);
0124 c.Range(-7.156863, -810349, 5.764706, 4951034);
0125 c.SetFillColor(0);
0126 c.SetBorderMode(0);
0127 c.SetBorderSize(3);
0128 c.SetGridx();
0129 c.SetGridy();
0130 c.SetTickx(1);
0131 c.SetTicky(1);
0132 c.SetLeftMargin(0.1669196);
0133 c.SetRightMargin(0.05918058);
0134 c.SetTopMargin(0.08233276);
0135 c.SetBottomMargin(0.1406518);
0136 c.SetFrameLineWidth(3);
0137 c.SetFrameBorderMode(0);
0138 c.SetFrameLineWidth(3);
0139 c.SetFrameBorderMode(0);
0140
0141 TGaxis::SetMaxDigits(3);
0142
0143 h1.SetLineColor(kRed);
0144 h1.SetLineWidth(3);
0145
0146
0147 TString TempEta = "eta", TempChi2 = "chi2";
0148 if (variable.Contains("pt")) {
0149 h1.SetXTitle("Track p_{T} (GeV)");
0150 } else if (variable.Contains("charge")) {
0151 h1.SetXTitle("Track charge (e)");
0152 } else if (variable.Data() == TempEta) {
0153 h1.SetXTitle("Track #eta");
0154 } else if (variable.Contains("phi")) {
0155 h1.SetXTitle("Track #phi (rad)");
0156 } else if (variable.Contains("theta")) {
0157 h1.SetXTitle("Track #theta (rad)");
0158 } else if (variable.Contains("d0")) {
0159 h1.SetXTitle("Track d_{0} (cm)");
0160 } else if (variable.Contains("dz")) {
0161 h1.SetXTitle("Track d_{z} (cm)");
0162 } else if (variable.Data() == TempChi2) {
0163 h1.SetXTitle("Track #chi^{2}");
0164 } else if (variable.Contains("ndof")) {
0165 h1.SetXTitle("Track #chi^{2} per NDF");
0166 } else {
0167 std::cout << "Title does not match anything in the categories defined!" << std::endl;
0168 }
0169
0170 h1.SetYTitle("Tracks (#)");
0171 h1.SetLabelSize(0.05);
0172 h1.GetXaxis()->SetLabelSize(0.05);
0173 h1.GetXaxis()->SetTitleSize(0.05);
0174 h1.GetYaxis()->SetLabelSize(0.05);
0175 h1.GetYaxis()->SetTitleSize(0.06);
0176 h1.GetXaxis()->SetTitleOffset(1.12);
0177 h1.Draw();
0178
0179
0180 TLatex Title = TLatex();
0181 Title.SetTextFont(42);
0182 Title.SetTextSize(0.039);
0183 Title.DrawLatexNDC(0.76, 0.94, "cosmic rays");
0184
0185
0186 TString PlotFormat[] = {"png", "pdf", "C"};
0187 for (int k = 0; k < 3; k++) {
0188 TString Format = variable + "." + PlotFormat[k];
0189 c.SaveAs(Format.Data());
0190 TString mv_folder_string = "mv " + Format + " Kinematical_Plots";
0191 gSystem->Exec(mv_folder_string.Data());
0192 }
0193 c.Close();
0194 }