{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "
\n", " \n", "
" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "

Searching for the Higgs boson in the H→γγ channel

\n", "\n", "\n", "SM Higgs signal Feynman diagram:\n", "
" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Introduction**\n", "Let's take a current ATLAS Open Data sample and create a histogram:" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Welcome to JupyROOT 6.18/04\n" ] } ], "source": [ "import ROOT\n", "from ROOT import TMath\n", "import time" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [], "source": [ "%jsroot on" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [], "source": [ "start = time.time()" ] }, { "cell_type": "code", "execution_count": 4, "metadata": { "scrolled": true }, "outputs": [], "source": [ "f = ROOT.TFile.Open(\"https://atlas-opendata.web.cern.ch/atlas-opendata/samples/2020/GamGam/Data/data_A.GamGam.root\")" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [], "source": [ "canvas = ROOT.TCanvas(\"Canvas\",\"cz\",800,600)" ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [], "source": [ "tree = f.Get(\"mini\")" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "430344" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "tree.GetEntries()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we're going to extract the photons variables" ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [], "source": [ "#Invariant mass histograms definition\n", "hist = ROOT.TH1F(\"h_M_Hyy\",\"Diphoton invariant-mass ; Invariant Mass m_{yy} [GeV] ; events\",30,105,160)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Here we're filling the variables defined above with the content of those inside the input ntuples." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We're creating a histogram for this example. The plan is to fill them with events." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We are selecting below a simple look for them." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The Higgs boson analysis implemented here considers Higgs boson decays into a photon-photon pair. The event selection criteria are:" ] }, { "cell_type": "code", "execution_count": 9, "metadata": { "scrolled": true }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "10000\n", "20000\n", "30000\n", "40000\n", "50000\n", "60000\n", "70000\n", "80000\n", "90000\n", "100000\n", "110000\n", "120000\n", "130000\n", "140000\n", "150000\n", "160000\n", "170000\n", "180000\n", "190000\n", "200000\n", "210000\n", "220000\n", "230000\n", "240000\n", "250000\n", "260000\n", "270000\n", "280000\n", "290000\n", "300000\n", "310000\n", "320000\n", "330000\n", "340000\n", "350000\n", "360000\n", "370000\n", "380000\n", "390000\n", "400000\n", "410000\n", "420000\n", "430000\n" ] } ], "source": [ "Photon_1 = ROOT.TLorentzVector()\n", "Photon_2 = ROOT.TLorentzVector()\n", "n = 0\n", "for event in tree:\n", " n += 1\n", " ## printing the evolution in number of events\n", " if(n%10000==0):\n", " print(n)\n", " ## checking the trigger \n", " if(tree.trigP):\n", " goodphoton_index = [0]*5\n", " goodphoton_n = 0\n", " photon_index = 0\n", " ## \n", " j=0\n", " ## looping the photons per event\n", " for j in range(tree.photon_n):\n", " ##\n", " if(tree.photon_isTightID[j]):\n", " ##\n", " if(tree.photon_pt[j] > 25000 and (TMath.Abs(tree.photon_eta[j]) < 2.37)\\\n", " and (TMath.Abs(tree.photon_eta[j]) < 1.37 or TMath.Abs(tree.photon_eta[j]) > 1.52)):\n", " ##\n", " goodphoton_n += 1 #count\n", " goodphoton_index[photon_index]=j\n", " photon_index += 1\n", " ## end Pt and eta pre-selection\n", " ## end on request of quality of the photon\n", " ## end looping photons in the current event\n", " \n", " ## Using the two selected photons\n", " if(goodphoton_n==2):\n", " ##\n", " goodphoton1_index = goodphoton_index[0]\n", " goodphoton2_index = goodphoton_index[1]\n", " ## Getting couple of photons with good isolation \n", " if((tree.photon_ptcone30[goodphoton1_index]/tree.photon_pt[goodphoton1_index] < 0.065)\\\n", " and (tree.photon_etcone20[goodphoton1_index] / tree.photon_pt[goodphoton1_index] < 0.065)):\n", " ##\n", " if((tree.photon_ptcone30[goodphoton2_index]/tree.photon_pt[goodphoton2_index] < 0.065)\\\n", " and (tree.photon_etcone20[goodphoton2_index] / tree.photon_pt[goodphoton2_index] < 0.065)):\n", " ##\n", " Photon_1.SetPtEtaPhiE(tree.photon_pt[goodphoton1_index]/1000., tree.photon_eta[goodphoton1_index],\\\n", " tree.photon_phi[goodphoton1_index],tree.photon_E[goodphoton1_index]/1000.)\n", " Photon_2.SetPtEtaPhiE(tree.photon_pt[goodphoton2_index]/1000., tree.photon_eta[goodphoton2_index],\\\n", " tree.photon_phi[goodphoton2_index],tree.photon_E[goodphoton2_index]/1000.)\n", " ## Adding the two TLorentz vectors\n", " Photon_12 = Photon_1 + Photon_2\n", " ## Filling with the mass of the gamma-gamma system\n", " hist.Fill(Photon_12.M())\n", " ## end isolation photon #2\n", " ## end isolation photon #1\n", " ## end 2-good photons\n", " ## end of trigger request\n", "## End loop in the events" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Final plot" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "
\n", "
\n", "\n", "\n" ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "hist.Draw(\"E\")\n", "canvas.Draw()" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "#### Log Scale" ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/html": [ "\n", "
\n", "
\n", "\n", "\n" ], "text/plain": [ "" ] }, "metadata": {}, "output_type": "display_data" } ], "source": [ "hist.Draw(\"E\")\n", "hist.SetMinimum(10)\n", "canvas.SetLogy()\n", "canvas.Draw()" ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Finished in 0 min 41 s\n" ] } ], "source": [ "end = time.time()\n", "duration = end-start\n", "print(\"Finished in {} min {} s\".format(int(duration//60),int(duration%60))) # Python3" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "**Done!**" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" }, "language_info": { "codemirror_mode": { "name": "ipython", "version": 3 }, "file_extension": ".py", "mimetype": "text/x-python", "name": "python", "nbconvert_exporter": "python", "pygments_lexer": "ipython3", "version": "3.6.9" } }, "nbformat": 4, "nbformat_minor": 1 }