Examples#

For those who learn by example


Random Histogram#

This example plots a random histogram using the ATLAS Style.

import ROOT as root
import atlasplots as aplt


def main():
    # Set the ATLAS Style
    aplt.set_atlas_style()

    # Create a figure and axes
    fig, ax = aplt.subplots(1, 1, name="fig1", figsize=(800, 600))

    # Define a distribution
    sqroot = root.TF1("sqroot", "x*gaus(0) + [3]*abs(sin(x)/x)", 0, 10)
    sqroot.SetParameters(10, 4, 1, 20)

    # Randomly fill the histogram according to the above distribution
    hist = root.TH1F("hist", "Random Histogram", 50, 0, 10)
    hist.FillRandom("sqroot", 20000)

    # Fit the histogram with the original distribution; store graphics func but do not draw
    hist.Fit("sqroot", "0")

    # Draw the histogram on these axes
    ax.plot(hist, label="Random Hist", labelfmt="F")

    # Draw the fit function
    sqroot.SetNpx(1000)
    ax.plot(sqroot, label="Fit", labelfmt="L", linecolor=root.kRed+1)

    # Add extra space at top of plot to make room for labels
    ax.add_margins(top=0.16)

    # Set axis titles
    ax.set_xlabel("X [GeV]")
    ax.set_ylabel("Events / 0.2 GeV")

    # Add the ATLAS Label
    aplt.atlas_label(text="Internal", loc="upper left")
    ax.text(0.2, 0.86, "#sqrt{s} = 13 TeV, 139 fb^{-1}", size=22, align=13)

    # Add legend
    ax.legend(loc=(0.65, 0.8, 0.95, 0.92))

    # Save the plot as a PDF
    fig.savefig("hist.pdf")


if __name__ == '__main__':
    root.gROOT.SetBatch()
    main()
hist.svg

Ratio Plot#

This example plots two random histograms and their ratios using the ATLAS Style.

import ROOT as root
import atlasplots as aplt


def main():
    # Set the ATLAS Style
    aplt.set_atlas_style()

    # Create a figure and axes
    fig, (ax1, ax2) = aplt.ratio_plot(name="fig1", figsize=(800, 800), hspace=0.05)

    # Define a distribution
    sqroot = root.TF1("sqroot", "x*gaus(0) + [3]*abs(sin(x)/x)", 0, 10)
    sqroot.SetParameters(10, 4, 1, 20)

    # Randomly fill two histograms according to the above distribution
    hist1 = root.TH1F("hist1", "Random Histogram 1", 50, 0, 10)
    hist1.FillRandom("sqroot", 20000)

    sqroot.SetParameters(10, 4, 1.1, 20)
    hist2 = root.TH1F("hist2", "Random Histogram 2", 50, 0, 10)
    hist2.FillRandom("sqroot", 20000)

    # Draw the histograms on these axes
    ax1.plot(hist1, linecolor=root.kRed+1, label="Red", labelfmt="L")
    ax1.plot(hist2, linecolor=root.kBlue+1, label="Blue", labelfmt="L")

    # Draw line at y=1 in ratio panel
    line = root.TLine(ax1.get_xlim()[0], 1, ax1.get_xlim()[1], 1)
    ax2.plot(line)

    # Calculate and draw the ratio
    ratio_hist = hist1.Clone("ratio_hist")
    ratio_hist.Divide(hist2)
    ax2.plot(ratio_hist, linecolor=root.kBlack)

    # Add extra space at top of plot to make room for labels
    ax1.add_margins(top=0.16)

    # Set axis titles
    ax2.set_xlabel("X [GeV]")
    ax1.set_ylabel("Events / 0.2 GeV")
    ax2.set_ylabel("Red / Blue", loc="centre")

    # Add extra space at top and bottom of ratio panel
    ax2.add_margins(top=0.1, bottom=0.1)

    # Go back to top axes to add labels
    ax1.cd()

    # Add the ATLAS Label
    aplt.atlas_label(text="Internal", loc="upper left")
    ax1.text(0.2, 0.84, "#sqrt{s} = 13 TeV, 139 fb^{-1}", size=22, align=13)

    # Add legend
    ax1.legend(loc=(0.78, 0.78, 1, 0.90))

    # Save the plot as a PDF
    fig.savefig("ratio.pdf")


if __name__ == '__main__':
    root.gROOT.SetBatch()
    main()
ratio.svg

Fit Histogram and Show Residuals#

This example generates and plots a random histogram, fits the original function to the generated distribution, and displays the fit residuals in the lower panel.

import ROOT as root
import atlasplots as aplt


def main():
    # Set the ATLAS Style
    aplt.set_atlas_style()

    # Create a figure and axes
    fig, (ax1, ax2) = aplt.ratio_plot(name="fig1", figsize=(800, 800), hspace=0.05)

    # Define a distribution
    sqroot = root.TF1("sqroot", "x*gaus(0) + [3]*abs(sin(x)/x)", 0, 10)
    sqroot.SetParameters(10, 4, 1, 20)

    # Randomly fill the histogram according to the above distribution
    hist = root.TH1F("hist", "Random Histogram", 50, 0, 10)
    hist.FillRandom("sqroot", 20000)

    # Fit the histogram with the original distribution; store graphics func but do not draw
    hist.Fit("sqroot", "0")

    # Draw the histogram on these axes
    ax1.plot(hist, "EP X0", linewidth=1, label="Random Hist", labelfmt="EP")

    # Draw the fit function
    sqroot.SetNpx(1000)
    ax1.plot(sqroot, linecolor=root.kRed+1, label="Fit", labelfmt="L")

    # Draw line at y=0 in residuals panel
    line = root.TLine(ax1.get_xlim()[0], 0, ax1.get_xlim()[1], 0)
    ax2.plot(line)

    # Calculate and draw the fit residuals
    resids = hist.Clone("resids")
    for i in range(1, resids.GetNbinsX()):
        resids.SetBinContent(i, hist.GetBinContent(i) - sqroot.Eval(hist.GetBinCenter(i)))
        resids.SetBinError(i, hist.GetBinError(i))

    ax2.plot(resids, "EP X0", linewidth=1)

    # Add extra space at top of plot to make room for labels
    ax1.add_margins(top=0.16)

    # Add extra space at top and bottom of residuals panel
    ax2.add_margins(top=0.1, bottom=0.1)

    # Set axis titles
    ax2.set_xlabel("X [GeV]")
    ax1.set_ylabel("Events / 0.2 GeV")
    ax2.set_ylabel("Data - Fit", loc="center")

    # Go back to top axes to add labels
    ax1.cd()

    # Add the ATLAS Label
    aplt.atlas_label(text="Internal", loc="upper left")
    ax1.text(0.2, 0.84, "#sqrt{s} = 13 TeV, 139 fb^{-1}", size=22, align=13)

    # Add legend
    ax1.legend(loc=(0.65, 0.78, 0.92, 0.90))

    # Save the plot as a PDF
    fig.savefig("fit_and_resids.pdf")


if __name__ == '__main__':
    root.gROOT.SetBatch()
    main()
fit_and_resids.svg

Data vs. MC Plot#

This example plots a stacked histograms representing the “signal” and “backgrounds”, and a graph showing the “data” superimposed. The Data/MC ratio is also displayed.

import ROOT as root
import atlasplots as aplt


def main():
    # Set the ATLAS Style
    aplt.set_atlas_style()

    # Create a figure and axes
    fig, (ax1, ax2) = aplt.ratio_plot(name="fig1", figsize=(800, 800), hspace=0.05)

    # Define "background", "signal" and "data" distributions
    bkg_func = root.TF1("bkg_func", "[0]*exp([1]*x)", 0, 10)
    bkg_func.SetParameters(20, -0.5)

    sig_func = root.TF1("sig_func", "x*gaus(0)", 0, 10)
    sig_func.SetParameters(10, 4, 1)

    # Randomly fill histograms according to the above distributions
    bkg_hist = root.TH1F("bkg_hist", "Background", 50, 0, 10)
    bkg_hist.FillRandom("bkg_func", 10000)
    bkg_hist.SetFillColor(root.kRed+1)
    bkg_hist.SetLineWidth(0)
    bkg_hist.Sumw2()

    sig_hist = root.TH1F("sig_hist", "Signal", 50, 0, 10)
    sig_hist.FillRandom("sig_func", 10000)
    sig_hist.SetFillColor(root.kAzure+1)
    sig_hist.SetLineWidth(0)
    sig_hist.Sumw2()

    data_hist = root.TH1F("data_hist", "Data", 50, 0, 10)
    data_hist.FillRandom("bkg_func", 10000)
    data_hist.FillRandom("sig_func", 10000)
    data_hist.Sumw2()
    data_hist.SetBinErrorOption(root.TH1.EBinErrorOpt.kPoisson)  # Use 68% Poisson errors

    # Stack the background and signal histograms
    bkg_and_sig = root.THStack("bkg_and_sig", "")
    bkg_and_sig.Add(bkg_hist)
    bkg_and_sig.Add(sig_hist)

    # Draw the stacked histogram on these axes
    ax1.plot(bkg_and_sig)

    # Plot the MC stat error as a hatched band
    err_band = aplt.root_helpers.hist_to_graph(
        bkg_and_sig.GetStack().Last(),
        show_bin_width=True
    )
    ax1.plot(err_band, "2", fillcolor=1, fillstyle=3254, linewidth=0)

    # ax1.set_yscale("log")  # uncomment to use log scale for y axis

    # Plot the data as a graph
    data_graph = aplt.root_helpers.hist_to_graph(data_hist)
    ax1.plot(data_graph, "P")

    # Use same x-range in lower axes as upper axes
    ax2.set_xlim(ax1.get_xlim())

    # Draw line at y=1 in ratio panel
    line = root.TLine(ax1.get_xlim()[0], 1, ax1.get_xlim()[1], 1)
    ax2.plot(line)

    # Plot the relative error on the ratio axes
    err_band_ratio = aplt.root_helpers.hist_to_graph(
        bkg_and_sig.GetStack().Last(),
        show_bin_width=True,
        norm=True
    )
    ax2.plot(err_band_ratio, "2", fillcolor=1, fillstyle=3254)

    # Calculate and draw the ratio
    ratio_hist = data_hist.Clone("ratio_hist")
    ratio_hist.Divide(bkg_and_sig.GetStack().Last())
    ratio_graph = aplt.root_helpers.hist_to_graph(ratio_hist)
    ax2.plot(ratio_graph, "P0")

    # Add extra space at top of plot to make room for labels
    ax1.add_margins(top=0.16)

    # Set axis titles
    ax2.set_xlabel("X [GeV]")
    ax1.set_ylabel("Events / 0.2 GeV")
    ax2.set_ylabel("Data / Pred.", loc="centre")

    ax2.set_ylim(0.75, 1.25)
    ax2.draw_arrows_outside_range(ratio_graph)

    # Go back to top axes to add labels
    ax1.cd()

    # Add the ATLAS Label
    aplt.atlas_label(text="Internal", loc="upper left")
    ax1.text(0.2, 0.84, "#sqrt{s} = 13 TeV, 139 fb^{-1}", size=22, align=13)

    # Add legend
    legend = ax1.legend(
        loc=(0.68, 0.65, 1 - root.gPad.GetRightMargin() - 0.03, 1 - root.gPad.GetTopMargin() - 0.03),
        textsize=22
    )
    legend.AddEntry(data_graph, "Data", "EP")
    legend.AddEntry(sig_hist, "Signal", "F")
    legend.AddEntry(bkg_hist, "Background", "F")
    legend.AddEntry(err_band, "MC Stat. Unc.", "F")

    # Save the plot as a PDF
    fig.savefig("data_vs_mc.pdf")


if __name__ == '__main__':
    root.gROOT.SetBatch()
    main()
data_vs_mc.svg

Graph from Numpy Array#

This example plots a graph directly from numpy arrays using the ATLAS Style.

import numpy as np

import ROOT as root
import atlasplots as aplt


def main():
    # Set the ATLAS Style
    aplt.set_atlas_style()

    # Create a figure and axes
    fig, ax = aplt.subplots(1, 1, name="fig1", figsize=(800, 600))

    # Plot a line with Gaussian noise
    x = np.arange(20)
    y = 2 * x + np.random.normal(size=20)
    yerr = np.random.normal(loc=np.linspace(1, 2, num=20), scale=0.1, size=20)

    graph = ax.graph(x, y, yerr=yerr, label="Data", labelfmt="EP")

    # Fit the graph; store graphics func but do not draw
    graph.Fit("pol1", "0")

    # Add extra space at top of plot to make room for labels
    ax.add_margins(top=0.18, left=0.05, right=0.05, bottom=0.05)

    # Plot fit function and extend its range to fill plot area
    func = graph.GetFunction("pol1")
    func.SetRange(*ax.get_xlim())
    func.SetNpx(1000)
    ax.plot(func, linecolor=root.kRed+1, expand=False, label="Fit", labelfmt="L")

    # Set axis titles
    ax.set_xlabel("x")
    ax.set_ylabel("y")

    # Add the ATLAS Label
    aplt.atlas_label(text="Simulation Internal", loc="upper left")

    # Add legend
    ax.legend(loc=(0.65, 0.8, 0.95, 0.92))

    # Save the plot as a PDF
    fig.savefig("numpy_graph.pdf")


if __name__ == '__main__':
    root.gROOT.SetBatch()
    main()
numpy_graph.svg

Random 2D Histogram#

This example plots a random, 2D histogram using the ATLAS Style.

import ROOT as root
import atlasplots as aplt


root.gInterpreter.Declare("""
void fillHist2D(TH2F* hist, int N) {
    float px, py;
    for (Int_t i = 0; i < N; i++) {
        gRandom->Rannor(px, py);
        hist->Fill(px + 4, 25 * py + 100);
        hist->Fill(0.5 * px + 8, 12 * py + 50, 0.1);
    }
}
""")


def main():
    # Set the ATLAS Style
    aplt.set_atlas_style()

    # Create a figure and axes
    fig, ax = aplt.subplots(1, 1, name="fig1", figsize=(800, 600))

    # Randomly fill a 2D histogram
    hist = root.TH2F("hist", "Random Histogram", 50, 0, 10, 40, 0, 200)
    root.fillHist2D(hist, 50000)

    # Draw the histogram on these axes
    ax.plot2d(hist, "COLZ")

    # Change pad margins to allow space for z-axis colour bar and for ATLAS label
    ax.set_pad_margins(right=0.20, top=0.08)

    # Set axis titles
    ax.set_xlabel("X")
    ax.set_ylabel("Y")
    ax.set_zlabel("Events / (0.2 #times 5)", titleoffset=1.2)

    # Add the ATLAS Label
    aplt.atlas_label(ax.pad.GetLeftMargin(), 0.97, text="Internal", align=13)
    ax.text(1 - ax.pad.GetRightMargin(), 0.97, "#sqrt{s} = 13 TeV, 139 fb^{-1}", size=22, align=33)

    # Save the plot as a PDF
    fig.savefig("hist2D.pdf")


if __name__ == '__main__':
    root.gROOT.SetBatch()
    main()
hist2D.svg