28parser = argparse.ArgumentParser()
29parser.add_argument(
"--lumi-scale", type=float, default=0.05,
30 help=
"Run only on a fraction of the total available 10 fb^-1 (only usable together with --full-dataset)")
31parser.add_argument(
"--full-dataset", action=
"store_true", default=
False,
32 help=
"Use the full dataset (use --lumi-scale to run only on a fraction of it)")
33parser.add_argument(
"-b", action=
"store_true", default=
False, help=
"Use ROOT batch mode")
34parser.add_argument(
"-t", action=
"store_true", default=
False, help=
"Use implicit multi threading (for the full dataset only possible with --lumi-scale 1.0)")
35args = parser.parse_args()
37if args.b: ROOT.gROOT.SetBatch(
True)
40if not args.full_dataset: lumi_scale = 0.05
41else: lumi_scale = args.lumi_scale
43print(
'Run on data corresponding to {:.1f} fb^-1 ...'.format(lumi * lumi_scale / 1000.0))
45if args.full_dataset: dataset_path =
"root://eospublic.cern.ch//eos/opendata/atlas/OutreachDatasets/2020-01-22"
46else: dataset_path =
"root://eospublic.cern.ch//eos/root-eos/reduced_atlas_opendata/singletop"
50files = json.load(open(os.path.join(os.path.dirname(os.path.abspath(__file__)),
"df107_SingleTopAnalysis.json")))
51processes = files.keys()
64 samples.append(sample)
65 df[sample] =
ROOT.RDataFrame(
"mini",
"{}/1lep/{}/{}.1lep.root".format(dataset_path, folder, sample))
68 if args.full_dataset
and lumi_scale < 1.0:
69 df[sample] = df[sample].
Range(
int(num_events * lumi_scale))
74ROOT.gInterpreter.Declare(
"""
75using cRVecF = const ROOT::RVecF &;
76using cRVecI = const ROOT::RVecI &;
77int FindGoodLepton(cRVecI goodlep, cRVecI type, cRVecF lep_pt, cRVecF lep_eta, cRVecF lep_phi, cRVecF lep_e, cRVecF trackd0pv, cRVecF tracksigd0pv, cRVecF z0)
79 int idx = -1; // Return -1 if no good lepton is found.
80 for(auto i = 0; i < type.size(); i++) {
81 if(!goodlep[i]) continue;
82 if (type[i] == 11 && abs(lep_eta[i]) < 2.47 && (abs(lep_eta[i]) < 1.37 || abs(lep_eta[i]) > 1.52) && abs(trackd0pv[i] / tracksigd0pv[i]) < 5) {
83 const ROOT::Math::PtEtaPhiEVector p(lep_pt[i], lep_eta[i], lep_phi[i], lep_e[i]);
84 if (abs(z0[i] * sin(p.Theta())) < 0.5) {
85 if (idx == -1) idx = i;
86 else return -1; // Accept only events with exactly one good lepton
89 if (type[i] == 13 && abs(lep_eta[i]) < 2.5 && abs(trackd0pv[i] / tracksigd0pv[i]) < 3) {
90 const ROOT::Math::PtEtaPhiEVector p(lep_pt[i], lep_eta[i], lep_phi[i], lep_e[i]);
91 if (abs(z0[i] * sin(p.Theta())) < 0.5) {
92 if (idx == -1) idx = i;
93 else return -1; // Accept only events with exactly one good lepton
103 df[s] = df[s].Filter(
"trigE || trigM")\
104 .Filter(
"met_et > 30000")
107 df[s] = df[s].Define(
"goodlep",
"lep_isTightID && lep_pt > 35000 && lep_ptcone30 / lep_pt < 0.1 && lep_etcone20 / lep_pt < 0.1")\
108 .Filter(
"ROOT::VecOps::Sum(goodlep) > 0")
111 df[s] = df[s].Define(
"idx_lep",
"FindGoodLepton(goodlep, lep_type, lep_pt, lep_eta, lep_phi, lep_E, lep_trackd0pvunbiased, lep_tracksigd0pvunbiased, lep_z0)")\
112 .Filter(
"idx_lep != -1")
116 df[s] = df[s].Define(
"mtw",
"sqrt(2 * lep_pt[idx_lep] * met_et * (1 - cos(lep_phi[idx_lep] - met_phi)))")\
117 .Filter(
"mtw > 60000")
120 df[s] = df[s].Filter(
"ROOT::VecOps::Sum(jet_pt > 30000 && abs(jet_eta) < 2.5) > 0")
123 df[s] = df[s].Define(
"goodjet",
"jet_pt > 60000 || abs(jet_eta) > 2.4 || jet_jvt > 0.59")\
124 .Filter(
"ROOT::VecOps::Sum(goodjet) == 2")\
125 .Define(
"goodbjet",
"goodjet && jet_MV2c10 > 0.8244273")\
126 .Filter(
"ROOT::VecOps::Sum(goodbjet) == 1")\
127 .Define(
"idx_tagged",
"ROOT::VecOps::ArgMax(goodjet && goodbjet)")\
128 .Define(
"idx_untagged",
"ROOT::VecOps::ArgMax(goodjet && !goodbjet)")
132 df[s] = df[s].Filter(
"abs(jet_eta[idx_untagged]) > 1.5 && abs(jet_eta[idx_tagged] - jet_eta[idx_untagged]) > 1.5")\
133 .Filter(
"lep_pt[idx_lep] + jet_pt[idx_tagged] + jet_pt[idx_untagged] + met_et > 195000")
138 df[s] = df[s].Define(
"weight",
"1.0")
141 if "single" in s: stop_norm =
"mcWeight / abs(mcWeight)"
142 else: stop_norm =
"mcWeight"
143 df[s] = df[s].Define(
"weight",
"scaleFactor_ELE * scaleFactor_MUON * scaleFactor_LepTRIGGER * scaleFactor_PILEUP * scaleFactor_BTAG * {} * {} / {} * {}".format(stop_norm, xsecs[s], sumws[s], lumi))
148ROOT.gInterpreter.Declare(
"""
149float ComputeTopMass(float lep_pt, float lep_eta, float lep_phi, float lep_e, float jet_pt, float jet_eta, float jet_phi, float jet_e, float met_et, float met_phi)
151 const ROOT::Math::PtEtaPhiEVector lep(lep_pt / 1000.0, lep_eta, lep_phi, lep_e / 1000.0);
152 const ROOT::Math::PtEtaPhiEVector met(met_et / 1000.0, 0, met_phi, met_et / 1000.0);
153 const ROOT::Math::PtEtaPhiEVector bjet(jet_pt / 1000.0, jet_eta, jet_phi, jet_e / 1000.0);
154 // Please note that we treat here the missing transverse energy as the neutrino, even though the z component is missing!
155 return (lep + met + bjet).M();
161 df[s] = df[s].Define(
"top_mass",
"ComputeTopMass(lep_pt[idx_lep], lep_eta[idx_lep], lep_phi[idx_lep], lep_E[idx_lep], jet_pt[idx_tagged], jet_eta[idx_tagged], jet_phi[idx_tagged], jet_E[idx_tagged], met_et, met_phi)")
162 histos[s] = df[s].Histo1D(
ROOT.RDF.TH1DModel(
"top_mass",
"", 10, 100, 400),
"top_mass",
"weight")
172def merge_histos(label):
174 for i, d
in enumerate(files[label]):
175 t = histos[d[1]].GetValue()
176 if i == 0: h = t.Clone()
178 h.SetNameTitle(label, label)
181data = merge_histos(
"data")
182twtb = merge_histos(
"twtb")
183singletop = merge_histos(
"singletop")
184wjets = merge_histos(
"wjets")
189ROOT.gROOT.SetStyle(
"ATLAS")
192c = ROOT.TCanvas(
"c",
"", 600, 600)
193pad = ROOT.TPad(
"upper_pad",
"", 0, 0, 1, 1)
200stack = ROOT.THStack()
203 [wjets, twtb, singletop],
204 [(222, 90, 106), (155, 152, 204), (208, 240, 193)]):
207 h.SetFillColor(ROOT.TColor.GetColor(*color))
210stack.GetXaxis().SetTitle(
"m_{W(l#nu)+b} [GeV]")
211stack.GetYaxis().SetTitle(
"Events")
212stack.GetYaxis().SetLabelSize(0.04)
213stack.GetYaxis().SetTitleSize(0.045)
214stack.GetXaxis().SetLabelSize(0.04)
215stack.GetXaxis().SetTitleSize(0.045)
217stack.SetMaximum(5000 * lumi_scale)
218stack.GetYaxis().ChangeLabel(1, -1, 0)
221data.SetMarkerStyle(20)
222data.SetMarkerSize(1.2)
224data.SetLineColor(ROOT.kBlack)
228legend = ROOT.TLegend(0.60, 0.65, 0.92, 0.92)
229legend.SetTextFont(42)
230legend.SetFillStyle(0)
231legend.SetBorderSize(0)
232legend.SetTextSize(0.035)
233legend.SetTextAlign(32)
234legend.AddEntry(data,
"Data" ,
"lep")
235legend.AddEntry(singletop,
"Single top + jet",
"f")
236legend.AddEntry(twtb,
"t#bar{t},Wt,t#bar{b}",
"f")
237legend.AddEntry(wjets,
"W+jets",
"f")
244text.SetTextSize(0.045)
245text.DrawLatex(0.21, 0.86,
"ATLAS")
247text.DrawLatex(0.21 + 0.16, 0.86,
"Open Data")
248text.SetTextSize(0.04)
249text.DrawLatex(0.21, 0.80,
"#sqrt{{s}} = 13 TeV, {:.1f} fb^{{-1}}".format(lumi * lumi_scale / 1000.0))
252c.SaveAs(
"df107_SingleTopAnalysis.png")
253print(
"Saved figure to df107_SingleTopAnalysis.png")
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTree,...
void RunGraphs(std::vector< RResultHandle > handles)
Trigger the event loop of multiple RDataFrames concurrently.
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
A struct which stores the parameters of a TH1D.