LOGO

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