LOGO

RestFrames  v1.0.1
RestFrames HEP Event Analysis Software Library
example_H_to_hh_to_4Wlnu.C
Go to the documentation of this file.
1 // RestFrames: particle physics event analysis library
3 // --------------------------------------------------------------------
4 // Copyright (c) 2014-2019, Christopher Rogan
13 //
14 // This file is part of RestFrames.
15 //
16 // RestFrames is free software; you can redistribute it and/or modify
17 // it under the terms of the GNU General Public License as published by
18 // the Free Software Foundation; either version 2 of the License, or
19 // (at your option) any later version.
20 //
21 // RestFrames is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with RestFrames. If not, see <http://www.gnu.org/licenses/>.
29 
30 #if (!defined(__CINT__) && !defined(__CLING__))
31 #define COMPILER
32 #endif
33 #if defined(__MAKECINT__) || defined(__ROOTCLING__) || defined(COMPILER)
34 #include "RestFrames/RestFrames.hh"
35 #else
36 RestFrames::RFKey ensure_autoload(1);
37 #endif
38 
39 using namespace RestFrames;
40 
41 void example_H_to_hh_to_4Wlnu(const std::string& output_name = "output_example_04.root"){
42 
43  double mH = 750.;
44  double mh = 125.;
45  double mW = 80.385; // GeV, PDG 2016
46  double wW = 2.085;
47  double mL = 0.106;
48  double mN = 0.;
49 
50  // number of events to generate
51  int Ngen = 100000;
52 
54  g_Log << LogInfo << "Initializing generator frames and tree..." << LogEnd;
56  ppLabGenFrame LAB_G("LAB_G","LAB");
57  DecayGenFrame H_G("H_G","H");
58  DecayGenFrame ha_G("ha_G","h_{a}");
59  DecayGenFrame hb_G("hb_G","h_{b}");
60  ResonanceGenFrame Waa_G("Waa_G","W_{a,a}");
61  ResonanceGenFrame Wab_G("Wab_G","W_{a,b}");
62  ResonanceGenFrame Wba_G("Wba_G","W_{b,a}");
63  ResonanceGenFrame Wbb_G("Wbb_G","W_{b,b}");
64  VisibleGenFrame Laa_G("Laa_G","#it{l}_{a,a}");
65  InvisibleGenFrame Naa_G("Naa_G","#nu_{a,a}");
66  VisibleGenFrame Lab_G("Lab_G","#it{l}_{a,b}");
67  InvisibleGenFrame Nab_G("Nab_G","#nu_{a,b}");
68  VisibleGenFrame Lba_G("Lba_G","#it{l}_{b,a}");
69  InvisibleGenFrame Nba_G("Nba_G","#nu_{b,a}");
70  VisibleGenFrame Lbb_G("Lbb_G","#it{l}_{b,b}");
71  InvisibleGenFrame Nbb_G("Nbb_G","#nu_{b,b}");
72 
73  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
74 
75  LAB_G.SetChildFrame(H_G);
76  H_G.AddChildFrame(ha_G);
77  H_G.AddChildFrame(hb_G);
78  ha_G.AddChildFrame(Waa_G);
79  ha_G.AddChildFrame(Wab_G);
80  hb_G.AddChildFrame(Wba_G);
81  hb_G.AddChildFrame(Wbb_G);
82  Waa_G.AddChildFrame(Laa_G);
83  Waa_G.AddChildFrame(Naa_G);
84  Wab_G.AddChildFrame(Lab_G);
85  Wab_G.AddChildFrame(Nab_G);
86  Wba_G.AddChildFrame(Lba_G);
87  Wba_G.AddChildFrame(Nba_G);
88  Wbb_G.AddChildFrame(Lbb_G);
89  Wbb_G.AddChildFrame(Nbb_G);
90 
91  if(LAB_G.InitializeTree())
92  g_Log << LogInfo << "...Successfully initialized generator tree" << LogEnd;
93  else
94  g_Log << LogError << "...Failed initializing generator tree" << LogEnd;
95  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
96 
97  // set Higgs masses
98  H_G.SetMass(mH);
99  ha_G.SetMass(mh);
100  hb_G.SetMass(mh);
101  // set W masses and widths
102  Waa_G.SetMass(mW);
103  Wab_G.SetMass(mW);
104  Wba_G.SetMass(mW);
105  Wbb_G.SetMass(mW);
106  Waa_G.SetWidth(wW);
107  Wab_G.SetWidth(wW);
108  Wba_G.SetWidth(wW);
109  Wbb_G.SetWidth(wW);
110  // set lepton masses
111  Laa_G.SetMass(mL);
112  Lab_G.SetMass(mL);
113  Lba_G.SetMass(mL);
114  Lbb_G.SetMass(mL);
115  Laa_G.SetPtCut(8.);
116  Lab_G.SetPtCut(8.);
117  Lba_G.SetPtCut(8.);
118  Lbb_G.SetPtCut(8.);
119  Laa_G.SetEtaCut(2.5);
120  Lab_G.SetEtaCut(2.5);
121  Lba_G.SetEtaCut(2.5);
122  Lbb_G.SetEtaCut(2.5);
123 
124  if(LAB_G.InitializeAnalysis())
125  g_Log << LogInfo << "...Successfully initialized generator analysis" << LogEnd;
126  else
127  g_Log << LogError << "...Failed initializing generator analysis" << LogEnd;
130 
132  g_Log << LogInfo << "Initializing reconstruction frames and trees..." << LogEnd;
134 
135  LabRecoFrame LAB_R("LAB_R","LAB");
136  DecayRecoFrame H_R("H_R","H");
137  DecayRecoFrame ha_R("ha_R","h_{a}");
138  DecayRecoFrame hb_R("hb_R","h_{b}");
139  DecayRecoFrame Waa_R("Waa_R","W_{a,a}");
140  DecayRecoFrame Wab_R("Wab_R","W_{a,b}");
141  DecayRecoFrame Wba_R("Wba_R","W_{b,a}");
142  DecayRecoFrame Wbb_R("Wbb_R","W_{b,b}");
143  VisibleRecoFrame Laa_R("Laa_R","#it{l}_{a,a}");
144  InvisibleRecoFrame Naa_R("Naa_R","#nu_{a,a}");
145  VisibleRecoFrame Lab_R("Lab_R","#it{l}_{a,b}");
146  InvisibleRecoFrame Nab_R("Nab_R","#nu_{a,b}");
147  VisibleRecoFrame Lba_R("Lba_R","#it{l}_{b,a}");
148  InvisibleRecoFrame Nba_R("Nba_R","#nu_{b,a}");
149  VisibleRecoFrame Lbb_R("Lbb_R","#it{l}_{b,b}");
150  InvisibleRecoFrame Nbb_R("Nbb_R","#nu_{b,b}");
151 
152  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
153 
154  LAB_R.SetChildFrame(H_R);
155  H_R.AddChildFrame(ha_R);
156  H_R.AddChildFrame(hb_R);
157  ha_R.AddChildFrame(Waa_R);
158  ha_R.AddChildFrame(Wab_R);
159  hb_R.AddChildFrame(Wba_R);
160  hb_R.AddChildFrame(Wbb_R);
161  Waa_R.AddChildFrame(Laa_R);
162  Waa_R.AddChildFrame(Naa_R);
163  Wab_R.AddChildFrame(Lab_R);
164  Wab_R.AddChildFrame(Nab_R);
165  Wba_R.AddChildFrame(Lba_R);
166  Wba_R.AddChildFrame(Nba_R);
167  Wbb_R.AddChildFrame(Lbb_R);
168  Wbb_R.AddChildFrame(Nbb_R);
169 
170  if(LAB_R.InitializeTree())
171  g_Log << LogInfo << "...Successfully initialized reconstruction tree" << LogEnd;
172  else
173  g_Log << LogError << "...Failed initializing reconstruction tree" << LogEnd;
174 
175  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
176 
177  // Invisible Group
178  InvisibleGroup INV_R("INV_R","4 #nu Jigsaws");
179  INV_R.AddFrame(Naa_R);
180  INV_R.AddFrame(Nab_R);
181  INV_R.AddFrame(Nba_R);
182  INV_R.AddFrame(Nbb_R);
183 
184  // Set nu nu mass equal to l l mass
185  SetMassInvJigsaw NuNuM_R("NuNuM_R", "M_{4 #nu} = f(4 #it{l})");
186  INV_R.AddJigsaw(NuNuM_R);
187 
188  SetRapidityInvJigsaw NuNuR_R("NuNuR_R", "#eta_{4 #nu} = #eta_{4 #it{l}}");
189  INV_R.AddJigsaw(NuNuR_R);
190  NuNuR_R.AddVisibleFrames(LAB_R.GetListVisibleFrames());
191 
192  ContraBoostInvJigsaw MinMh_R("MinMh_R","min M_{h}, M_{h}^{ a}= M_{h}^{ b}");
193  INV_R.AddJigsaw(MinMh_R);
194  MinMh_R.AddInvisibleFrames(ha_R.GetListInvisibleFrames(), 0);
195  MinMh_R.AddInvisibleFrames(hb_R.GetListInvisibleFrames(), 1);
196  MinMh_R.AddVisibleFrames(ha_R.GetListVisibleFrames(), 0);
197  MinMh_R.AddVisibleFrames(hb_R.GetListVisibleFrames(), 1);
198 
199  //ContraBoostInvJigsaw MinMWa_R("MinMWa_R","min M_{W}, M_{W}^{a,a}= M_{W}^{a,b}");
200  MinMassesSqInvJigsaw MinMWa_R("MinMWa_R","min M_{W}, M_{W}^{a,a}= M_{W}^{a,b}", 2);
201  //MinMassDiffInvJigsaw MinMWa_R("MinMWa_R","min M_{W}, M_{W}^{a,a}= M_{W}^{a,b}", 2);
202  //MaxProbBreitWignerInvJigsaw MinMWa_R("MinMWa_R","min M_{W}, M_{W}^{a,a}= M_{W}^{a,b}", 2);
203  INV_R.AddJigsaw(MinMWa_R);
204  MinMWa_R.AddInvisibleFrames(Waa_R.GetListInvisibleFrames(), 0);
205  MinMWa_R.AddInvisibleFrames(Wab_R.GetListInvisibleFrames(), 1);
206  MinMWa_R.AddVisibleFrames(Waa_R.GetListVisibleFrames(), 0);
207  MinMWa_R.AddVisibleFrames(Wab_R.GetListVisibleFrames(), 1);
208  MinMWa_R.AddMassFrame(Lba_R, 0);
209  MinMWa_R.AddMassFrame(Lbb_R, 1);
210 
211  //ContraBoostInvJigsaw MinMWb_R("MinMWa_R","min M_{W}, M_{W}^{b,a}= M_{W}^{b,b}");
212  MinMassesSqInvJigsaw MinMWb_R("MinMWa_R","min M_{W}, M_{W}^{b,a}= M_{W}^{b,b}", 2);
213  //MinMassDiffInvJigsaw MinMWb_R("MinMWb_R","min M_{W}, M_{W}^{b,a}= M_{W}^{b,b}", 2);
214  //MaxProbBreitWignerInvJigsaw MinMWb_R("MinMWa_R","min M_{W}, M_{W}^{b,a}= M_{W}^{b,b}", 2);
215  INV_R.AddJigsaw(MinMWb_R);
216  MinMWb_R.AddInvisibleFrames(Wba_R.GetListInvisibleFrames(), 0);
217  MinMWb_R.AddInvisibleFrames(Wbb_R.GetListInvisibleFrames(), 1);
218  MinMWb_R.AddVisibleFrames(Wba_R.GetListVisibleFrames(), 0);
219  MinMWb_R.AddVisibleFrames(Wbb_R.GetListVisibleFrames(), 1);
220  MinMWb_R.AddMassFrame(Laa_R, 0);
221  MinMWb_R.AddMassFrame(Lab_R, 1);
222 
223  // Combinatoric Group for leptons
224  CombinatoricGroup VIS_R("VIS_R","Lepton Combinatoric Jigsaws");
225  VIS_R.AddFrame(Lab_R);
226  VIS_R.AddFrame(Lbb_R);
227  // lepton frames must have at least one element
228  VIS_R.SetNElementsForFrame(Lab_R, 1);
229  VIS_R.SetNElementsForFrame(Lbb_R, 1);
230 
231  MinMassesCombJigsaw MinMll_R("MinMll_R", "min Mll");
232  VIS_R.AddJigsaw(MinMll_R);
233  MinMll_R.AddCombFrame(Lab_R, 0);
234  MinMll_R.AddCombFrame(Lbb_R, 1);
235  MinMll_R.AddObjectFrames(Laa_R+Lab_R, 0);
236  MinMll_R.AddObjectFrames(Lba_R+Lbb_R, 1);
237 
238  if(LAB_R.InitializeAnalysis())
239  g_Log << LogInfo << "...Successfully initialized analysis" << LogEnd;
240  else
241  g_Log << LogError << "...Failed initializing analysis" << LogEnd;
242 
245 
246  TreePlot* treePlot = new TreePlot("TreePlot","TreePlot");
247 
248  treePlot->SetTree(LAB_G);
249  treePlot->Draw("GenTree", "Generator Tree", true);
250 
251  treePlot->SetTree(LAB_R);
252  treePlot->Draw("RecoTree", "Reconstruction Tree");
253 
254  treePlot->SetTree(INV_R);
255  treePlot->Draw("InvTree", "Invisible Jigsaws", true);
256 
257  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
258 
259  std::string plot_title = "H #rightarrow hh #rightarrow W(#it{l} #nu) W(#it{l} #nu)";
260  HistPlot* histPlot = new HistPlot("HistPlot", plot_title);
261 
262  const HistPlotCategory& cat_Gen = histPlot->GetNewCategory("Gen", "Generator");
263  const HistPlotCategory& cat_Reco = histPlot->GetNewCategory("Reco", "Reconstruction");
264 
265  const HistPlotVar& MH = histPlot->GetNewVar("MH", "M_{H} / m_{H}^{ true}", 0., 2.);
266  const HistPlotVar& Mha = histPlot->GetNewVar("Mha", "M_{h}^{ a} / m_{h}^{ true}", 0., 2.);
267  const HistPlotVar& Mhb = histPlot->GetNewVar("Mhb", "M_{h}^{ b} / m_{h}^{ true}", 0., 2.);
268  const HistPlotVar& MWaa = histPlot->GetNewVar("MWaa", "M_{W a,a}", 0., mW*2., "[GeV]");
269  const HistPlotVar& MWab = histPlot->GetNewVar("MWab", "M_{W a,b}", 0., mW*2., "[GeV]");
270  const HistPlotVar& MWba = histPlot->GetNewVar("MWba", "M_{W b,a}", 0., mW*2., "[GeV]");
271  const HistPlotVar& MWbb = histPlot->GetNewVar("MWbb", "M_{W b,b}", 0., mW*2., "[GeV]");
272  const HistPlotVar& cosH = histPlot->GetNewVar("cosH","cos #theta_{H}", -1., 1.);
273  const HistPlotVar& cosha = histPlot->GetNewVar("cosha","cos #theta_{h a}", -1., 1.);
274  const HistPlotVar& coshb = histPlot->GetNewVar("coshb","cos #theta_{h b}", -1., 1.);
275  const HistPlotVar& cosWaa = histPlot->GetNewVar("cosWaa","cos #theta_{Waa}", -1., 1.);
276  const HistPlotVar& cosWab = histPlot->GetNewVar("cosWab","cos #theta_{Wab}", -1., 1.);
277  const HistPlotVar& DcosH = histPlot->GetNewVar("DcosH","#theta_{H} - #theta_{H}^{ true}",
278  -acos(-1.)/2., acos(-1.)/2.);
279  const HistPlotVar& Dcosha = histPlot->GetNewVar("Dcosha","#theta_{h a} - #theta_{h a}^{ true}",
280  -acos(-1.)/2., acos(-1.)/2.);
281  const HistPlotVar& Dcoshb = histPlot->GetNewVar("Dcoshb","#theta_{h b} - #theta_{h b}^{ true}",
282  -acos(-1.)/2., acos(-1.)/2.);
283  const HistPlotVar& DcosWaa = histPlot->GetNewVar("DcosWaa","#theta_{W a,a} - #theta_{W a,a}^{ true}",
284  -acos(-1.)/2., acos(-1.)/2.);
285  const HistPlotVar& DcosWab = histPlot->GetNewVar("DcosWab","#theta_{W a,b} - #theta_{W a,b}^{ true}",
286  -acos(-1.)/2., acos(-1.)/2.);
287 
288  histPlot->AddPlot(MWaa, MWab, cat_Gen+cat_Reco);
289  histPlot->AddPlot(MWaa, MWbb, cat_Gen+cat_Reco);
290  histPlot->AddPlot(MH, cat_Reco);
291  histPlot->AddPlot(Mha, cat_Reco);
292  histPlot->AddPlot(Mha, Mhb, cat_Reco);
293  histPlot->AddPlot(MH, Mha, cat_Reco);
294  histPlot->AddPlot(Mha, MWaa, cat_Reco);
295  histPlot->AddPlot(DcosH, cat_Reco);
296  histPlot->AddPlot(Dcosha, cat_Reco);
297  histPlot->AddPlot(MH, DcosH, cat_Reco);
298  histPlot->AddPlot(Mha, Dcosha, cat_Reco);
299  histPlot->AddPlot(DcosH, Dcosha, cat_Reco);
300  histPlot->AddPlot(Dcosha, Dcoshb, cat_Reco);
301 
302  histPlot->AddPlot(DcosWaa, cat_Reco);
303  histPlot->AddPlot(DcosWaa, DcosWab, cat_Reco);
304 
307 
308  for(int igen = 0; igen < Ngen; igen++){
309  if(igen%((std::max(Ngen,10))/10) == 0)
310  g_Log << LogInfo << "Generating event " << igen << " of " << Ngen << LogEnd;
311 
312  // generate event
313  LAB_G.ClearEvent(); // clear the gen tree
314  double PTH = mH*gRandom->Rndm()*0.1;
315  LAB_G.SetTransverseMomentum(PTH); // give the Higgs some Pt
316  LAB_G.AnalyzeEvent(); // generate a new event
317 
318  // analyze event
319  LAB_R.ClearEvent(); // clear the reco tree
320 
321  TVector3 MET = LAB_G.GetInvisibleMomentum(); // Get the MET from gen tree
322  MET.SetZ(0.);
323  INV_R.SetLabFrameThreeVector(MET); // Set the MET in reco tree
324 
325  // std::vector<RFKey> L_ID; // ID for tracking leptons in tree
326  // L_ID.push_back(VIS_R.AddLabFrameFourVector(Laa_G.GetFourVector(),-1));
327  // L_ID.push_back(VIS_R.AddLabFrameFourVector(Lab_G.GetFourVector(), 1));
328  // L_ID.push_back(VIS_R.AddLabFrameFourVector(Lba_G.GetFourVector(),-1));
329  // L_ID.push_back(VIS_R.AddLabFrameFourVector(Lbb_G.GetFourVector(), 1));
330 
331  Laa_R.SetLabFrameFourVector(Laa_G.GetFourVector(),-1);
332  Lba_R.SetLabFrameFourVector(Lba_G.GetFourVector(),-1);
333  std::vector<RFKey> L_ID; // ID for tracking anti-leptons in tree
334  L_ID.push_back(VIS_R.AddLabFrameFourVector(Lab_G.GetFourVector(), 1));
335  L_ID.push_back(VIS_R.AddLabFrameFourVector(Lbb_G.GetFourVector(), 1));
336 
337  // Laa_R.SetLabFrameFourVector(Laa_G.GetFourVector(),-1);
338  // Lab_R.SetLabFrameFourVector(Lab_G.GetFourVector(), 1);
339  // Lba_R.SetLabFrameFourVector(Lba_G.GetFourVector(),-1);
340  // Lbb_R.SetLabFrameFourVector(Lbb_G.GetFourVector(), 1);
341 
342  LAB_R.AnalyzeEvent(); // analyze the event
343 
344  // Generator-level observables
345  double MHgen = H_G.GetMass();
346  double Mhagen = ha_G.GetMass();
347  double Mhbgen = hb_G.GetMass();
348  double cosHgen = H_G.GetCosDecayAngle();
349  double coshagen = ha_G.GetCosDecayAngle();
350  double coshbgen = hb_G.GetCosDecayAngle();
351  cosH = cosHgen;
352  cosha = coshagen;
353  coshb = coshbgen;
354  double cosWaagen = Waa_G.GetCosDecayAngle();
355  double cosWabgen = Wab_G.GetCosDecayAngle();
356  MWaa = Waa_G.GetMass();
357  MWab = Wab_G.GetMass();
358  MWba = Wba_G.GetMass();
359  MWbb = Wbb_G.GetMass();
360 
361  histPlot->Fill(cat_Gen);
362 
363  // Reconstruction-level observables
364  MH = H_R.GetMass()/MHgen;
365  //MH = 2.*H_R.GetListVisibleFrames().GetEnergy(H_R)/MHgen;
366  //Mha = ha_R.GetMass()/Mhagen;
367  //Mhb = hb_R.GetMass()/Mhbgen;
368  Mha = 2.*ha_R.GetListVisibleFrames().GetEnergy(ha_R)/Mhagen;
369  Mhb = 2.*hb_R.GetListVisibleFrames().GetEnergy(hb_R)/Mhbgen;
370  cosH = H_R.GetCosDecayAngle();
371  cosha = ha_R.GetCosDecayAngle();
372  coshb = hb_R.GetCosDecayAngle();
373  cosWaa = Waa_R.GetCosDecayAngle();
374  cosWab = Wab_R.GetCosDecayAngle();
375  MWaa = Waa_R.GetMass();
376  MWab = Wab_R.GetMass();
377  MWba = Wba_R.GetMass();
378  MWbb = Wbb_R.GetMass();
379  DcosH = asin(sqrt(1.-cosH*cosH)*cosHgen-sqrt(1.-cosHgen*cosHgen)*cosH);
380  Dcosha = asin(sqrt(1.-cosha*cosha)*coshagen-sqrt(1.-coshagen*coshagen)*cosha);
381  Dcoshb = asin(sqrt(1.-coshb*coshb)*coshbgen-sqrt(1.-coshbgen*coshbgen)*coshb);
382  DcosWaa = asin(sqrt(1.-cosWaa*cosWaa)*cosWaagen-sqrt(1.-cosWaagen*cosWaagen)*cosWaa);
383  DcosWab = asin(sqrt(1.-cosWab*cosWab)*cosWabgen-sqrt(1.-cosWabgen*cosWabgen)*cosWab);
384 
385 
386  histPlot->Fill(cat_Reco);
387 
388  // g_Log << LogInfo;
389 // g_Log << "Correct assignments " << std::endl;
390 // g_Log << (VIS_R.GetFrame(L_ID[0]).GetParentFrame().GetParentFrame() == VIS_R.GetFrame(L_ID[1]).GetParentFrame().GetParentFrame()) << " ";
391 // g_Log << (VIS_R.GetFrame(L_ID[2]).GetParentFrame().GetParentFrame() == VIS_R.GetFrame(L_ID[3]).GetParentFrame().GetParentFrame()) << " " << std::endl;
392 // g_Log << (VIS_R.GetFrame(L_ID[0]).GetParentFrame().GetParentFrame() == VIS_R.GetFrame(L_ID[3]).GetParentFrame().GetParentFrame()) << " ";
393 // g_Log << (VIS_R.GetFrame(L_ID[2]).GetParentFrame().GetParentFrame() == VIS_R.GetFrame(L_ID[1]).GetParentFrame().GetParentFrame()) << " " << std::endl;
394 // g_Log << (Laa_R+Lab_R).GetMass() << " " << (Lba_R+Lbb_R).GetMass() << std::endl;
395 // g_Log << (Laa_R+Lbb_R).GetMass() << " " << (Lba_R+Lab_R).GetMass() << std::endl;
396 // g_Log << LogEnd;
397  }
398 
399  histPlot->Draw();
400 
401  LAB_G.PrintGeneratorEfficiency();
402 
403  g_Log << LogInfo << "Finished" << LogEnd;
404 }
405 
406 # ifndef __CINT__ // main function for stand-alone compilation
407 int main(){
408  example_H_to_hh_to_4Wlnu();
409  return 0;
410 }
411 #endif
RestFrames::TreePlot::Draw
void Draw(const std::string &name="", const std::string &title="", bool invert_colors=false, bool invert_node_colors=false)
Draws tree.
Definition: TreePlot.cc:127
RestFrames::InvisibleGenFrame
Definition: InvisibleGenFrame.hh:41
RestFrames::HistPlotCategory
Definition: HistPlotCategory.hh:40
RestFrames::ppLabGenFrame
Definition: ppLabGenFrame.hh:40
RestFrames::DecayGenFrame
Definition: DecayGenFrame.hh:43
RestFrames::HistPlot
Definition: HistPlot.hh:44
RestFrames.hh
RestFrames::ContraBoostInvJigsaw
Definition: ContraBoostInvJigsaw.hh:37
RestFrames::VisibleRecoFrame
Definition: VisibleRecoFrame.hh:41
RestFrames::ResonanceGenFrame
Definition: ResonanceGenFrame.hh:40
RestFrames::MinMassesSqInvJigsaw
Definition: MinMassesSqInvJigsaw.hh:37
RestFrames::SetRapidityInvJigsaw
Definition: SetRapidityInvJigsaw.hh:37
RestFrames::LabRecoFrame
Definition: LabRecoFrame.hh:44
RestFrames::MinMassesCombJigsaw
Definition: MinMassesCombJigsaw.hh:40
RestFrames::HistPlot::GetNewVar
HistPlotVar const & GetNewVar(const std::string &name, const std::string &title, double minval, double maxval, const std::string &unit="")
Defines new variable for plotting.
Definition: HistPlot.cc:86
RestFrames::VisibleGenFrame
Definition: VisibleGenFrame.hh:41
RestFrames::HistPlot::Fill
void Fill(double weight=1.)
Fills an entry for every variable according to each's current value.
Definition: HistPlot.cc:242
RestFrames::InvisibleRecoFrame
Definition: InvisibleRecoFrame.hh:41
RestFrames::TreePlot::SetTree
void SetTree(const RestFrame &frame)
Sets a tree that shows how a frame is connected to its children.
Definition: TreePlot.cc:160
RestFrames::TreePlot
Definition: TreePlot.hh:50
RestFrames::HistPlot::GetNewCategory
HistPlotCategory const & GetNewCategory(const std::string &name, const std::string &title)
Defines new category for plotting.
Definition: HistPlot.cc:95
RestFrames::SetMassInvJigsaw
Definition: SetMassInvJigsaw.hh:37
RestFrames::CombinatoricGroup
Definition: CombinatoricGroup.hh:39
RestFrames::RFKey
Definition: RFKey.hh:38
RestFrames::HistPlotVar
Definition: HistPlotVar.hh:40
RestFrames::HistPlot::AddPlot
void AddPlot(const HistPlotVar &var, RestFrames::RFList< const HistPlotCategory > cats=RestFrames::RFList< const RestFrames::HistPlotCategory >(), bool invert_colors=false)
Defines new 1D plot.
Definition: HistPlot.cc:102
RestFrames::DecayRecoFrame
Definition: DecayRecoFrame.hh:43
RestFrames::InvisibleGroup
Definition: InvisibleGroup.hh:40
RestFrames::HistPlot::Draw
void Draw(bool invert_colors=false)
Draws every histogram.
Definition: HistPlot.cc:272