#include #include #include #include #include #include #include #include "TTree.h" #include "TFile.h" #include "TString.h" #include "TSystem.h" #include "TH1.h" #include "TH2.h" #include "TH3.h" #include "TCanvas.h" #include "TStyle.h" #include "TChain.h" #include "TEllipse.h" #include "TLine.h" #include "TMarker.h" #include "TVector3.h" #include "TRotation.h" #include "TMath.h" #include "TROOT.h" #include "TF1.h" #include "TGraph.h" #include "TGraphErrors.h" #include "TLatex.h" #include "TPolyLine.h" #define DEUTERON_MASS 1.875614 #define PROTON_MASS 0.93827 #define RAD2DEG 57.295779 #define DEG2RAD 0.017453293 using namespace std; void HallA_style() { gROOT->SetStyle("Plain"); gStyle->SetPaperSize(TStyle::kUSLetter); gStyle->SetPaperSize(18,22); gStyle->SetOptFit(1111); gStyle->SetPalette(1); gStyle->SetNdivisions(508); gStyle->SetCanvasColor(10); gStyle->SetPadTopMargin(.05); gStyle->SetPadLeftMargin(.15); gStyle->SetPadRightMargin(.1); gStyle->SetPadBottomMargin(.10); gStyle->SetTitleYOffset(1.7); gStyle->SetTitleXOffset(1.2); gStyle->SetLabelFont(42,"X"); gStyle->SetLabelFont(42,"Y"); gStyle->SetTitleX(0.2); gStyle->SetTitleW(0.7); gStyle->SetTitleFillColor(18); // prepare gStyle to be useful // 1 = solid // 2 = long dash (30 10) // 3 = dotted (4 8) // 4 = dash-dot (15 12 4 12) // 5 = short dash ( 15 15 ) // 6 = dash-dot-dot gStyle->SetLineStyleString(1,"[]"); gStyle->SetLineStyleString(2,"[30 10]"); gStyle->SetLineStyleString(3,"[4 8]"); gStyle->SetLineStyleString(4,"[15 12 4 12]"); gStyle->SetLineStyleString(5,"[15 15]"); gStyle->SetLineStyleString(6,"[15 12 4 12 4 12]"); gStyle->SetLabelSize(0.04,"X"); gStyle->SetLabelSize(0.04,"Y"); gStyle->SetNdivisions(508,"Y"); gStyle->SetOptDate(0); gStyle->SetDateY(.98); gStyle->SetStripDecimals(kFALSE); gStyle->SetOptStat(1111); } //void LiveTimeCalculation(int runnumber = 2300) //void LiveTimeCalculation(int runnumber = 2414) //void LiveTimeCalculation(int runnumber = 3374) //void LiveTimeCalculation(int runnumber = 2618) void LiveTimeCalculation(int runnumber = 2418) { HallA_style(); //--- Module variables Double_t DL_t1[100]; Double_t DL_t2[100]; Double_t DL_t3[100]; Double_t DL_t5[100]; Double_t DL_t8[100]; Double_t DL_bbretime[100]; Int_t Ndata_DL_t3; Int_t Ndata_DL_t1; Int_t Ndata_DL_t2; Int_t Ndata_DL_t5; Int_t Ndata_DL_t8; Int_t Ndata_DL_bbretime; Int_t Ndata_DL_edtmbb; Double_t DL_evtypebits; Double_t g0hel_L_helicity; Double_t N_clocktime; Double_t N_PS1, N_PS2, N_PS3, N_PS4, N_PS5, N_PS6, N_PS7, N_PS8; Int_t Count_bbretime = 0; Int_t Count_t3 = 0; Int_t Count_t3_Plus = 0; Int_t Count_t3_Minus = 0; Int_t Count_TDC_t3 = 0; Int_t Count_t5 = 0; Int_t Count_t5_Plus = 0; Int_t Count_t5_Minus = 0; Int_t Count_TDC_t5 = 0; double corrfac = 0.0; double t3ratio = 0.0; double t3ratio_Plus = 0.0; double t3ratio_Minus = 0.0; double t5ratio = 0.0; double t5ratio_Plus = 0.0; double t5ratio_Minus = 0.0; //--- Epics variables Double_t IBC0L02Current; Double_t hac_bcm_average; //--- Scaler variables Double_t clock1024rate; Double_t clock1024count; Double_t BBretimerate; Double_t BBretimecount; Double_t scalerT3count; Double_t scalerT3rate; Double_t scalerT5count; Double_t scalerT5rate; Double_t scaler_G1_T3rate; Double_t scaler_G1_T3count; Double_t scaler_G1_T5rate; Double_t scaler_G1_T5count; Double_t scaler_G2_T3rate; Double_t scaler_G2_T3count; Double_t scaler_G2_T5rate; Double_t scaler_G2_T5count; TChain *t1 =new TChain("T"); TString filename = Form("./eed_BigBiteScalar_%d.root", runnumber); t1->Add(filename.Data()); assert(t1); //############################ VARIABLES ############################ t1->SetBranchAddress( "DL.evtypebits", &DL_evtypebits); t1->SetBranchAddress( "DL.t1", &DL_t1); t1->SetBranchAddress( "DL.t2", &DL_t2); t1->SetBranchAddress( "DL.t3", &DL_t3); t1->SetBranchAddress( "DL.t5", &DL_t5); t1->SetBranchAddress( "Ndata.DL.t1", &Ndata_DL_t1); t1->SetBranchAddress( "Ndata.DL.t2", &Ndata_DL_t2); t1->SetBranchAddress( "Ndata.DL.t3", &Ndata_DL_t3); t1->SetBranchAddress( "Ndata.DL.t5", &Ndata_DL_t5); t1->SetBranchAddress( "Ndata.DL.bbretime", &Ndata_DL_bbretime); t1->SetBranchAddress( "Ndata.DL.edtmbb", &Ndata_DL_edtmbb); t1->SetBranchAddress( "adchel.L.helicity", &g0hel_L_helicity); t1->SetBranchAddress( "N.clocktime", &N_clocktime); t1->SetBranchAddress( "N.PS1", &N_PS1); t1->SetBranchAddress( "N.PS2", &N_PS2); t1->SetBranchAddress( "N.PS3", &N_PS3); t1->SetBranchAddress( "N.PS4", &N_PS4); t1->SetBranchAddress( "N.PS5", &N_PS5); t1->SetBranchAddress( "N.PS6", &N_PS6); t1->SetBranchAddress( "N.PS7", &N_PS7); t1->SetBranchAddress( "N.PS8", &N_PS8); t1->SetBranchAddress( "evbbite_clock", &clock1024rate); t1->SetBranchAddress( "evbbite_clockcnt", &clock1024count); t1->SetBranchAddress( "evbbite_BBretime", &BBretimerate); t1->SetBranchAddress( "evbbite_BBretimecnt", &BBretimecount); t1->SetBranchAddress( "evbbite_T3", &scalerT3rate); t1->SetBranchAddress( "evbbite_T3cnt", &scalerT3count); t1->SetBranchAddress( "evbbite_T5", &scalerT5rate); t1->SetBranchAddress( "evbbite_T5cnt", &scalerT5count); //--- This thing works only for the case, then Target spin // signal is negative. One scaler module is missing // due to address mixup. t1->SetBranchAddress( "evbbite_G1_T3", &scaler_G1_T3rate); t1->SetBranchAddress( "evbbite_G1_T3cnt", &scaler_G1_T3count); t1->SetBranchAddress( "evbbite_G1_T5", &scaler_G1_T5rate); t1->SetBranchAddress( "evbbite_G1_T5cnt", &scaler_G1_T5count); t1->SetBranchAddress( "evbbite_G2_T3", &scaler_G2_T3rate); t1->SetBranchAddress( "evbbite_G2_T3cnt", &scaler_G2_T3count); t1->SetBranchAddress( "evbbite_G2_T5", &scaler_G2_T5rate); t1->SetBranchAddress( "evbbite_G2_T5cnt", &scaler_G2_T5count); t1->SetBranchAddress( "IBC0L02Current", &IBC0L02Current); t1->SetBranchAddress( "hac_bcm_average", &hac_bcm_average); //############################ Histograms ############################ TGraph *gr_bbretime_ratio_vs_time = new TGraph(); TGraph *gr_t3_ratio_vs_time = new TGraph(); TGraph *gr_TDCt3_ratio_vs_time = new TGraph(); TGraph *gr_t3_ratio_vs_time_Plus = new TGraph(); TGraph *gr_t3_ratio_vs_time_Minus = new TGraph(); TGraph *gr_t3_ratio_vs_time_PMRatio = new TGraph(); TGraph *gr_t5_ratio_vs_time = new TGraph(); TGraph *gr_TDCt5_ratio_vs_time = new TGraph(); TGraph *gr_t5_ratio_vs_time_Plus = new TGraph(); TGraph *gr_t5_ratio_vs_time_Minus = new TGraph(); TGraph *gr_t5_ratio_vs_time_PMRatio = new TGraph(); //############################ Analysis ############################ int num = t1->GetEntries(); //num = 30000; cout<<"Number of entries : "<GetEntry(i); if ((i%1000==0) && i>0) cout<1) cerr<<"We have more than one CODA event per Event !!! \a"<0) corrfac = (1.0*i)/(1.0*BBretimecount); gr_bbretime_ratio_vs_time->SetPoint(i, time, corrfac); //########################### T3 calculation ################## //--- Number of events accepted as T3 events if (((Int_t)DL_evtypebits&(1<<3))==(1<<3)) Count_t3++; if (scalerT3count >0) t3ratio = Count_t3/scalerT3count; gr_t3_ratio_vs_time->SetPoint(i, time, t3ratio); //--- Number of T3s in TDCs Count_TDC_t3 += Ndata_DL_t3; double TDCt3ratio = 0.0; if (scalerT3count >0) TDCt3ratio = Count_TDC_t3/scalerT3count; gr_TDCt3_ratio_vs_time->SetPoint(i, time, TDCt3ratio); // First lets determine beam helicity; if (g0hel_L_helicity > 0 ){ if (((Int_t)DL_evtypebits&(1<<3))==(1<<3)) Count_t3_Plus++; } else if (g0hel_L_helicity < 0 ) { if (((Int_t)DL_evtypebits&(1<<3))==(1<<3)) Count_t3_Minus++; } if (scaler_G1_T3count >0) t3ratio_Plus = Count_t3_Plus/scaler_G1_T3count; if (scaler_G2_T3count >0) t3ratio_Minus = Count_t3_Minus/scaler_G2_T3count; //gr_t3_ratio_vs_time_Plus->SetPoint(i, time, 3*Count_t3_Minus); //gr_t3_ratio_vs_time_Minus->SetPoint(i, time, scaler_G2_T3count-20000); gr_t3_ratio_vs_time_Plus->SetPoint(i, time, t3ratio_Plus ); gr_t3_ratio_vs_time_Minus->SetPoint(i, time, t3ratio_Minus); if (t3ratio_Minus > 0)gr_t3_ratio_vs_time_PMRatio->SetPoint(i, time, t3ratio_Plus/t3ratio_Minus); //########################### T5 calculation ################## //--- Number of events accepted as T5 events if (((Int_t)DL_evtypebits&(1<<5))==(1<<5)) Count_t5++; if (scalerT5count >0) t5ratio = Count_t5/scalerT5count; gr_t5_ratio_vs_time->SetPoint(i, time, t5ratio); //--- Number of T5s in TDCs Count_TDC_t5 += Ndata_DL_t5; double TDCt5ratio = 0.0; if (scalerT5count >0) TDCt5ratio = Count_TDC_t5/scalerT5count; gr_TDCt5_ratio_vs_time->SetPoint(i, time, TDCt5ratio); // First lets determine beam helicity; if (g0hel_L_helicity > 0 ){ if (((Int_t)DL_evtypebits&(1<<5))==(1<<5)) Count_t5_Plus++; } else if (g0hel_L_helicity < 0 ) { if (((Int_t)DL_evtypebits&(1<<5))==(1<<5)) Count_t5_Minus++; } if (scaler_G1_T5count >0) t5ratio_Plus = Count_t5_Plus/scaler_G1_T5count; if (scaler_G2_T5count >0) t5ratio_Minus = Count_t5_Minus/scaler_G2_T5count; gr_t5_ratio_vs_time_Plus->SetPoint(i, time, t5ratio_Plus ); gr_t5_ratio_vs_time_Minus->SetPoint(i, time, t5ratio_Minus); if (t5ratio_Minus > 0)gr_t5_ratio_vs_time_PMRatio->SetPoint(i, time, t5ratio_Plus/t5ratio_Minus); } //--- ReadPrescale factors: cout<<"######################################################################\n"<Delete(); delete t1; t1=NULL; //--- Canvas C1 TCanvas *c1 = new TCanvas("c1","CorrFac",800,800); c1->Divide(1,1); c1->cd(1); gPad->SetGrid(); gr_bbretime_ratio_vs_time->SetMarkerStyle(0); gr_bbretime_ratio_vs_time->SetMarkerColor(kRed); gr_bbretime_ratio_vs_time->SetMarkerSize(1); gr_bbretime_ratio_vs_time->GetYaxis()->SetTitle(" bbretime_{CODA}/bbretime_{Scaler} [1]"); gr_bbretime_ratio_vs_time->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_bbretime_ratio_vs_time->SetTitle("Correction Factor"); gr_bbretime_ratio_vs_time->Draw("AP"); double vx,vy; gr_bbretime_ratio_vs_time->GetPoint(gr_bbretime_ratio_vs_time->GetN()-1, vx, vy); TString CorrFacString = Form("Correction Factor = %f ", vy); TLatex TextCorrFac; TextCorrFac.SetTextAlign(12); TextCorrFac.SetTextSize(0.04); TextCorrFac.SetTextColor(kRed); TextCorrFac.DrawLatex(10.0,1.1,CorrFacString); TString FigureFileName = Form("./Figures/figure_CorrFac_%d.png",runnumber); c1->Print(FigureFileName.Data()); //--- Canvas C2 TCanvas *c2 = new TCanvas("c2","T3 Plots",1200,800); c2->Divide(2,1); c2->cd(1); gPad->SetGrid(); gr_t3_ratio_vs_time->SetMarkerStyle(0); gr_t3_ratio_vs_time->SetMarkerColor(kRed); gr_t3_ratio_vs_time->SetMarkerSize(1); gr_t3_ratio_vs_time->GetYaxis()->SetTitle(" T3_{CODA}/T3_{Scaler} [1]"); gr_t3_ratio_vs_time->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t3_ratio_vs_time->SetTitle("T3 ratio"); gr_t3_ratio_vs_time->GetYaxis()->SetRangeUser(0.0, 1.0); gr_t3_ratio_vs_time->Draw("AP"); gr_TDCt3_ratio_vs_time->SetMarkerStyle(0); gr_TDCt3_ratio_vs_time->SetMarkerColor(kBlue); gr_TDCt3_ratio_vs_time->SetMarkerSize(1); gr_TDCt3_ratio_vs_time->Draw("P"); //--- Now lets determine T3 live time... need to correct for // non-sync of scalers and events and for prescalefactor c2->cd(2); gPad->SetGrid(); TGraph *gr_t3_livetime = new TGraph(); double livetime=0.0; for (int i = 0; iGetN(); i++) { double t1, t2,y, fac; gr_t3_ratio_vs_time->GetPoint(i, t1, y); gr_bbretime_ratio_vs_time->GetPoint(i, t2, fac); if (fabs(t1-t2)>1E-6) cerr<<"Joining two different points!!!\a"<SetPoint(i,t1, livetime); } cout<<"--> T3 Live time at the end is: "<SetMarkerStyle(0); gr_t3_livetime->SetMarkerColor(kGreen); gr_t3_livetime->SetMarkerSize(1); gr_t3_livetime->GetYaxis()->SetTitle(" T3 Live-Time [1]"); gr_t3_livetime->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t3_livetime->SetTitle("T3 Live-Time"); gr_t3_livetime->GetYaxis()->SetRangeUser(0.0, 1.0); gr_t3_livetime->Draw("AP"); TString T3LiveTimeString = Form("T3 LiveTime = %f ", livetime); TLatex TextT3LiveTime; TextT3LiveTime.SetTextAlign(12); TextT3LiveTime.SetTextSize(0.04); TextT3LiveTime.SetTextColor(kGreen); TextT3LiveTime.DrawLatex(10.0,0.5,T3LiveTimeString); FigureFileName = Form("./Figures/figure_T3Plots_%d.png",runnumber); c2->Print(FigureFileName.Data()); //--- Canvas C3 TCanvas *c3 = new TCanvas("c3","T5 Plots",1200,800); c3->Divide(2,1); c3->cd(1); gPad->SetGrid(); gr_t5_ratio_vs_time->SetMarkerStyle(0); gr_t5_ratio_vs_time->SetMarkerColor(kRed); gr_t5_ratio_vs_time->SetMarkerSize(1); gr_t5_ratio_vs_time->GetYaxis()->SetTitle(" T5_{CODA}/T5_{Scaler} [1]"); gr_t5_ratio_vs_time->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t5_ratio_vs_time->SetTitle("T5 ratio"); gr_t5_ratio_vs_time->GetYaxis()->SetRangeUser(0.0, 1.0); gr_t5_ratio_vs_time->Draw("AP"); gr_TDCt5_ratio_vs_time->SetMarkerStyle(0); gr_TDCt5_ratio_vs_time->SetMarkerColor(kBlue); gr_TDCt5_ratio_vs_time->SetMarkerSize(1); gr_TDCt5_ratio_vs_time->Draw("P"); //--- Now lets determine T5 live time... need to correct for // non-sync of scalers and events and for prescalefactor c3->cd(2); gPad->SetGrid(); TGraph *gr_t5_livetime = new TGraph(); for (int i = 0; iGetN(); i++) { double t1, t2,y, fac; gr_t5_ratio_vs_time->GetPoint(i, t1, y); gr_bbretime_ratio_vs_time->GetPoint(i, t2, fac); if (fabs(t1-t2)>1E-6) cerr<<"Joining two different points!!!\a"<SetPoint(i,t1, livetime); } cout<<"--> T5 Live time at the end is: "<SetMarkerStyle(0); gr_t5_livetime->SetMarkerColor(kBlue); gr_t5_livetime->SetMarkerSize(1); gr_t5_livetime->GetYaxis()->SetTitle(" T5 Live-Time [1]"); gr_t5_livetime->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t5_livetime->SetTitle("T5 Live-Time"); gr_t5_livetime->GetYaxis()->SetRangeUser(0.0, 1.0); gr_t5_livetime->Draw("AP"); TString T5LiveTimeString = Form("T5 LiveTime = %f ", livetime); TLatex TextT5LiveTime; TextT5LiveTime.SetTextAlign(12); TextT5LiveTime.SetTextSize(0.04); TextT5LiveTime.SetTextColor(kBlue); TextT5LiveTime.DrawLatex(10.0,0.5,T5LiveTimeString); FigureFileName = Form("./Figures/figure_T5Plots_%d.png",runnumber); c3->Print(FigureFileName.Data()); //--- Canvas C4 TCanvas *c4 = new TCanvas("c4","T3 Gated Plots",1200,800); c4->Divide(2,1); c4->cd(1); gPad->SetGrid(); gr_t3_ratio_vs_time_Plus->SetMarkerStyle(0); gr_t3_ratio_vs_time_Plus->SetMarkerColor(kRed); gr_t3_ratio_vs_time_Plus->SetLineColor(kRed); gr_t3_ratio_vs_time_Plus->SetMarkerSize(1); gr_t3_ratio_vs_time_Plus->GetYaxis()->SetTitle(" T3_{CODA}/T3_{Scaler} [1]"); gr_t3_ratio_vs_time_Plus->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t3_ratio_vs_time_Plus->SetTitle("Gated (+/-) T3 ratio"); gr_t3_ratio_vs_time_Plus->GetYaxis()->SetRangeUser(0.0, 1.0); gr_t3_ratio_vs_time_Plus->Draw("AP"); gr_t3_ratio_vs_time_Minus->SetMarkerStyle(0); gr_t3_ratio_vs_time_Minus->SetMarkerColor(kBlue); gr_t3_ratio_vs_time_Minus->SetLineColor(kBlue); gr_t3_ratio_vs_time_Minus->SetMarkerSize(1); gr_t3_ratio_vs_time_Minus->Draw("P"); c4->cd(2); gPad->SetGrid(); gr_t3_ratio_vs_time_PMRatio->SetMarkerStyle(0); gr_t3_ratio_vs_time_PMRatio->SetMarkerColor(kRed); gr_t3_ratio_vs_time_PMRatio->SetLineColor(kRed); gr_t3_ratio_vs_time_PMRatio->SetMarkerSize(1); gr_t3_ratio_vs_time_PMRatio->GetYaxis()->SetTitle(" LT(+)/LT(-) [1]"); gr_t3_ratio_vs_time_PMRatio->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t3_ratio_vs_time_PMRatio->SetTitle("T3 (+/-) Live Time ratio"); gr_t3_ratio_vs_time_PMRatio->GetYaxis()->SetRangeUser(0.5, 1.5); gr_t3_ratio_vs_time_PMRatio->Draw("AP"); double t, livetimeratio; gr_t3_ratio_vs_time_PMRatio->GetPoint(gr_t3_ratio_vs_time_PMRatio->GetN()-1, t, livetimeratio); TString T3LiveTimeRatioString = Form("T3 (+/-) Livetime ratio = %f ", livetimeratio); TLatex TextT3LiveTimeRatio; TextT3LiveTimeRatio.SetTextAlign(12); TextT3LiveTimeRatio.SetTextSize(0.04); TextT3LiveTimeRatio.SetTextColor(kRed); TextT3LiveTimeRatio.DrawLatex(10.0,1.05,T3LiveTimeRatioString); FigureFileName = Form("./Figures/figure_T3PlotsGated_%d.png",runnumber); c4->Print(FigureFileName.Data()); //--- Canvas C5 TCanvas *c5 = new TCanvas("c5","T5 Gated Plots",1200,800); c5->Divide(2,1); c5->cd(1); gPad->SetGrid(); gr_t5_ratio_vs_time_Plus->SetMarkerStyle(0); gr_t5_ratio_vs_time_Plus->SetMarkerColor(kRed); gr_t5_ratio_vs_time_Plus->SetLineColor(kRed); gr_t5_ratio_vs_time_Plus->SetMarkerSize(1); gr_t5_ratio_vs_time_Plus->GetYaxis()->SetTitle(" T5_{CODA}/T5_{Scaler} [1]"); gr_t5_ratio_vs_time_Plus->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t5_ratio_vs_time_Plus->SetTitle("Gated (+/-) T5 ratio"); gr_t5_ratio_vs_time_Plus->GetYaxis()->SetRangeUser(0.0, 1.0); gr_t5_ratio_vs_time_Plus->Draw("AP"); gr_t5_ratio_vs_time_Minus->SetMarkerStyle(0); gr_t5_ratio_vs_time_Minus->SetMarkerColor(kBlue); gr_t5_ratio_vs_time_Minus->SetLineColor(kBlue); gr_t5_ratio_vs_time_Minus->SetMarkerSize(1); gr_t5_ratio_vs_time_Minus->Draw("P"); c5->cd(2); gPad->SetGrid(); gr_t5_ratio_vs_time_PMRatio->SetMarkerStyle(0); gr_t5_ratio_vs_time_PMRatio->SetMarkerColor(kRed); gr_t5_ratio_vs_time_PMRatio->SetLineColor(kRed); gr_t5_ratio_vs_time_PMRatio->SetMarkerSize(1); gr_t5_ratio_vs_time_PMRatio->GetYaxis()->SetTitle(" LT(+)/LT(-) [1]"); gr_t5_ratio_vs_time_PMRatio->GetXaxis()->SetTitle(" Scaler-Time [s]"); gr_t5_ratio_vs_time_PMRatio->SetTitle("T5 (+/-) Live Time ratio"); gr_t5_ratio_vs_time_PMRatio->GetYaxis()->SetRangeUser(0.5, 1.5); gr_t5_ratio_vs_time_PMRatio->Draw("AP"); gr_t5_ratio_vs_time_PMRatio->GetPoint(gr_t5_ratio_vs_time_PMRatio->GetN()-1, t, livetimeratio); TString T5LiveTimeRatioString = Form("T5 (+/-) Livetime ratio = %f ", livetimeratio); TLatex TextT5LiveTimeRatio; TextT5LiveTimeRatio.SetTextAlign(12); TextT5LiveTimeRatio.SetTextSize(0.04); TextT5LiveTimeRatio.SetTextColor(kRed); TextT5LiveTimeRatio.DrawLatex(10.0,1.05,T5LiveTimeRatioString); FigureFileName = Form("./Figures/figure_T5PlotsGated_%d.png",runnumber); c5->Print(FigureFileName.Data()); }