LOGO

RestFrames  v1.0.1
RestFrames HEP Event Analysis Software Library
example_DiStop_to_bXp_bXm_to_blNblN.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_DiStop_to_bXp_bXm_to_blNblN(std::string output_name =
42  "output_DiStop_to_bXp_bXm_to_blNblN.root"){
43 
44  double mT = 750.; // stop mass
45  double mN = 100.; // sneutrino mass
46  double mB = 4.18;
47  double mL = 0.106;
48 
49  int NMX = 8; // number of different chargino masses to evaluate
50 
51  bool fix_MXb = false; // fix the "b" chargino mass while varying "a"?
52  double mXb = 600.;
53 
54  // number of events to generate (per chargino mass)
55  int Ngen = 10000;
56 
58  g_Log << LogInfo << "Initializing generator frames and tree..." << LogEnd;
60  ppLabGenFrame LAB_Gen("LAB_Gen","LAB");
61  DecayGenFrame CM_Gen("CM_Gen","CM");
62  DecayGenFrame Ta_Gen("Ta_Gen","#tilde{t}_{a}");
63  DecayGenFrame Tb_Gen("Tb_Gen","#tilde{t}_{b}");
64  DecayGenFrame Xa_Gen("Xa_Gen","#tilde{#chi}_{a}^{ #pm}");
65  DecayGenFrame Xb_Gen("Xb_Gen","#tilde{#chi}_{b}^{ #mp}");
66  VisibleGenFrame Ba_Gen("Ba_Gen","b_{a}");
67  VisibleGenFrame La_Gen("La_Gen","#it{l}_{a}");
68  InvisibleGenFrame Na_Gen("Na_Gen","#tilde{#nu}_{a}");
69  VisibleGenFrame Bb_Gen("Bb_Gen","b_{b}");
70  VisibleGenFrame Lb_Gen("Lb_Gen","#it{l}_{b}");
71  InvisibleGenFrame Nb_Gen("Nb_Gen","#tilde{#nu}_{b}");
72 
73  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
74 
75  LAB_Gen.SetChildFrame(CM_Gen);
76  CM_Gen.AddChildFrame(Ta_Gen);
77  CM_Gen.AddChildFrame(Tb_Gen);
78  Ta_Gen.AddChildFrame(Ba_Gen);
79  Ta_Gen.AddChildFrame(Xa_Gen);
80  Tb_Gen.AddChildFrame(Bb_Gen);
81  Tb_Gen.AddChildFrame(Xb_Gen);
82  Xa_Gen.AddChildFrame(La_Gen);
83  Xa_Gen.AddChildFrame(Na_Gen);
84  Xb_Gen.AddChildFrame(Lb_Gen);
85  Xb_Gen.AddChildFrame(Nb_Gen);
86 
87  if(LAB_Gen.InitializeTree())
88  g_Log << LogInfo << "...Successfully initialized generator tree" << LogEnd;
89  else
90  g_Log << LogError << "...Failed initializing generator tree" << LogEnd;
91 
92  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
93 
94  if(NMX < 1) return;
95 
96  // non-resonant stop production
97  CM_Gen.SetVariableMass();
98  // set stop masses
99  Ta_Gen.SetMass(mT); Tb_Gen.SetMass(mT);
100  // set chargino masses
101  Xa_Gen.SetMass(0.5*(mT+mN)); Xb_Gen.SetMass(0.5*(mT+mN));
102  // set B masses
103  Ba_Gen.SetMass(mB); Bb_Gen.SetMass(mB);
104  // set : masses
105  La_Gen.SetMass(mL); Lb_Gen.SetMass(mL);
106  // set sneutrino masses
107  Na_Gen.SetMass(mN); Nb_Gen.SetMass(mN);
108 
109  // set b-jet/lepton pT/eta cuts
110  Ba_Gen.SetPtCut(20.); Bb_Gen.SetPtCut(20.);
111  Ba_Gen.SetEtaCut(2.5); Bb_Gen.SetEtaCut(2.5);
112  La_Gen.SetPtCut(15.); Lb_Gen.SetPtCut(15.);
113  La_Gen.SetEtaCut(2.5); Lb_Gen.SetEtaCut(2.5);
114 
115  if(LAB_Gen.InitializeAnalysis())
116  g_Log << LogInfo << "...Successfully initialized generator analysis" << LogEnd;
117  else
118  g_Log << LogError << "...Failed initializing generator analysis" << LogEnd;
122  g_Log << LogInfo << "Initializing reconstruction frames and tree..." << LogEnd;
124 
125  LabRecoFrame LAB("LAB","LAB");
126  DecayRecoFrame CM("CM","CM");
127  DecayRecoFrame Ta("Ta","#tilde{t}_{a}");
128  DecayRecoFrame Tb("Tb","#tilde{t}_{b}");
129  DecayRecoFrame Xa("Xa","#tilde{#chi}_{a}^{ #pm}");
130  DecayRecoFrame Xb("Xb","#tilde{#chi}_{b}^{ #mp}");
131  VisibleRecoFrame Ba("Ba","b_{a}");
132  VisibleRecoFrame La("La","#it{l}_{a}");
133  InvisibleRecoFrame Na("Na","#tilde{#nu}_{a}");
134  VisibleRecoFrame Bb("Bb","b_{b}");
135  VisibleRecoFrame Lb("Lb","#it{l}_{b}");
136  InvisibleRecoFrame Nb("Nb","#tilde{#nu}_{b}");
137 
138  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
139 
140  LAB.SetChildFrame(CM);
141  CM.AddChildFrame(Ta);
142  CM.AddChildFrame(Tb);
143  Ta.AddChildFrame(Ba);
144  Ta.AddChildFrame(Xa);
145  Tb.AddChildFrame(Bb);
146  Tb.AddChildFrame(Xb);
147  Xa.AddChildFrame(La);
148  Xa.AddChildFrame(Na);
149  Xb.AddChildFrame(Lb);
150  Xb.AddChildFrame(Nb);
151 
152  if(LAB.InitializeTree())
153  g_Log << LogInfo << "...Successfully initialized reconstruction tree" << LogEnd;
154  else
155  g_Log << LogError << "...Failed initializing reconstruction tree" << LogEnd;
156 
157  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
158 
160  std::string group_name;
161 
162  // Invisible Group
163  group_name = "#tilde{#nu} #tilde{#nu} Jigsaws";
164  InvisibleGroup INV("INV", group_name);
165  INV.AddFrame(Na);
166  INV.AddFrame(Nb);
167  CombinatoricGroup B("VIS","b-jet Jigsaws");
168  B.AddFrame(Ba);
169  B.AddFrame(Bb);
170  B.SetNElementsForFrame(Ba, 1);
171  B.SetNElementsForFrame(Bb, 1);
172 
174  std::string jigsaw_name;
175 
176  // Minimize difference Mt jigsaws
177  jigsaw_name = "M_{#tilde{#nu} #tilde{#nu}} ~ m_{#it{l} #it{l}}";
178  SetMassInvJigsaw NuNuM("NuNuM", jigsaw_name);
179  INV.AddJigsaw(NuNuM);
180 
181  jigsaw_name = "#eta_{#tilde{#nu} #tilde{#nu}} = #eta_{b #it{l} b #it{l}}";
182  SetRapidityInvJigsaw NuNuR("NuNuR", jigsaw_name);
183  INV.AddJigsaw(NuNuR);
184  NuNuR.AddVisibleFrames(CM.GetListVisibleFrames());
185 
186  jigsaw_name = "min ( M_{#tilde{t} a}- M_{#tilde{t} b} )^{2}";
187  //MinMassDiffInvJigsaw MinDeltaMt("MinDeltaMt", jigsaw_name, 2);
188  ContraBoostInvJigsaw MinDeltaMt("MinDeltaMt", jigsaw_name);
189  INV.AddJigsaw(MinDeltaMt);
190  MinDeltaMt.AddInvisibleFrame(Na, 0);
191  MinDeltaMt.AddInvisibleFrame(Nb, 1);
192  MinDeltaMt.AddVisibleFrames(La+Ba, 0);
193  MinDeltaMt.AddVisibleFrames(Lb+Bb, 1);
194  MinDeltaMt.AddMassFrame(La, 0);
195  MinDeltaMt.AddMassFrame(Lb, 1);
196 
197  // b-jet combinatoric jigsaws for all trees
198  jigsaw_name = "Minimize M(b #it{l} )_{a} , M(b #it{l} )_{b}";
199 
200  MinMassesCombJigsaw MinBL("MinBL", jigsaw_name);
201  B.AddJigsaw(MinBL);
202  MinBL.AddFrames(La+Ba,0);
203  MinBL.AddFrames(Lb+Bb,1);
204 
205  if(LAB.InitializeAnalysis())
206  g_Log << LogInfo << "...Successfully initialized analysis" << LogEnd;
207  else
208  g_Log << LogError << "...Failed initializing analysis" << LogEnd;
209 
212 
213  TreePlot* treePlot = new TreePlot("TreePlot","TreePlot");
214 
215  treePlot->SetTree(LAB_Gen);
216  treePlot->Draw("GenTree", "Generator Tree", true);
217 
218  treePlot->SetTree(LAB);
219  treePlot->Draw("RecoTree", "Reconstruction Tree");
220 
221  treePlot->SetTree(B);
222  treePlot->Draw("VisTree", "b-jet Jigsaws", true);
223 
224  treePlot->SetTree(INV);
225  treePlot->Draw("InvTree", "Inivisible Jigsaws", true);
226 
227  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
228 
229  std::string plot_title =
230  "#tilde{t} #tilde{t} #rightarrow 2x [ b #tilde{#chi}^{ #pm}(#it{l} #tilde{#nu}) ]";
231  char smasses[200];
232  sprintf(smasses, "m_{#tilde{t}}= %.0f, m_{#tilde{#nu}}= %.0f", mT, mN);
233  HistPlot* histPlot = new HistPlot("HistPlot", plot_title+" , "+std::string(smasses));
234 
235  std::string R_MXN = "#frac{m_{#tilde{#chi}^{ #pm}} - m_{#tilde{#nu}}}";
236  std::string R_MXD = "{m_{#tilde{t}} - m_{#tilde{#nu}}}";
238  for(int m = 0; m < NMX; m++){
239  char smass[200], scat[50];
240  sprintf(scat, "MX%d_%d", m+1, NMX+1);
241  sprintf(smass, " = #frac{%d}{%d}", m+1, NMX+1);
242 
243  cat_list += histPlot->GetNewCategory(scat, "#scale[0.7]{"+R_MXN+R_MXD+std::string(smass)+"}");
244  }
245 
246  //const HistPlotCategory& cat_Reco = histPlot->GetNewCategory("Reco", smasses);
247  const HistPlotCategory& cat_Reco = histPlot->GetNewCategory("Reco",
248  plot_title+" , "+std::string(smasses));
249 
250  const HistPlotVar& MCM = histPlot->GetNewVar("MCM", "M_{CM} / m_{CM}^{ true}", 0., 2.);
251  const HistPlotVar& Eb_Ta = histPlot->GetNewVar("Eb_Ta",
252  "E_{b a}^{#tilde{t} a} / E_{b a}^{#tilde{t} a gen}", 0., 2.);
253  const HistPlotVar& Eb_Tb = histPlot->GetNewVar("Eb_Tb",
254  "E_{b b}^{#tilde{t} b} / E_{b b}^{#tilde{t} b gen}", 0., 2.);
255  const HistPlotVar& El_Xa = histPlot->GetNewVar("El_Xa",
256  "E_{#it{l} a}^{#tilde{#chi}^{ #pm} a} / E_{#it{l} a}^{#tilde{#chi}^{ #pm} a true}", 0., 2.);
257  const HistPlotVar& El_Xb = histPlot->GetNewVar("El_Xb",
258  "E_{#it{l} b}^{#tilde{#chi}^{ #pm} b} / E_{#it{l} b}^{#tilde{#chi}^{ #pm} b true}", 0., 2.);
259  const HistPlotVar& cosCM = histPlot->GetNewVar("cosCM","cos #theta_{CM}", -1., 1.);
260  const HistPlotVar& cosTa = histPlot->GetNewVar("cosTa","cos #theta_{#tilde{t} a}", -1., 1.);
261  const HistPlotVar& cosTb = histPlot->GetNewVar("cosTb","cos #theta_{#tilde{t} b}", -1., 1.);
262  const HistPlotVar& cosXa = histPlot->GetNewVar("cosXa","cos #theta_{#tilde{#chi}^{ #pm} a}", -1., 1.);
263  const HistPlotVar& cosXb = histPlot->GetNewVar("cosXb","cos #theta_{#tilde{#chi}^{ #pm} b}", -1., 1.);
264  const HistPlotVar& DcosCM = histPlot->GetNewVar("DcosCM","#theta_{CM} - #theta_{CM}^{true}",
265  -acos(-1.)/2., acos(-1.)/2.);
266  const HistPlotVar& DcosTa = histPlot->GetNewVar("Dcosta",
267  "#theta_{#tilde{t} a} - #theta_{#tilde{t} a}^{true}", -acos(-1.)/2., acos(-1.)/2.);
268  const HistPlotVar& DcosTb = histPlot->GetNewVar("Dcostb",
269  "#theta_{#tilde{t} b} - #theta_{#tilde{t} b}^{true}", -acos(-1.)/2., acos(-1.)/2.);
270  const HistPlotVar& DcosXa = histPlot->GetNewVar("DcosXa",
271  "#theta_{#tilde{#chi}^{ #pm} a} - #theta_{#tilde{#chi}^{ #pm} a}^{true}", -acos(-1.)/2., acos(-1.)/2.);
272  const HistPlotVar& DcosXb = histPlot->GetNewVar("DcosXb",
273  "#theta_{#tilde{#chi}^{ #pm} b} - #theta_{#tilde{#chi}^{ #pm} b}^{true}", -acos(-1.)/2., acos(-1.)/2.);
274  const HistPlotVar& RMX = histPlot->GetNewVar("RMX", "#scale[0.6]{"+R_MXN+R_MXD+"}", 0., 1.);
275 
276  histPlot->AddPlot(MCM, cat_list);
277  histPlot->AddPlot(Eb_Ta, cat_list);
278  histPlot->AddPlot(El_Xa, cat_list);
279  histPlot->AddPlot(DcosCM, cat_list);
280  histPlot->AddPlot(DcosTa, cat_list);
281  histPlot->AddPlot(DcosXa, cat_list);
282 
283  histPlot->AddPlot(MCM, Eb_Ta, cat_list[NMX/2]);
284  histPlot->AddPlot(MCM, El_Xa, cat_list[NMX/2]);
285  histPlot->AddPlot(Eb_Ta, Eb_Tb, cat_list[NMX/2]);
286  histPlot->AddPlot(El_Xa, El_Xb, cat_list[NMX/2]);
287  histPlot->AddPlot(Eb_Ta, El_Xa, cat_list[NMX/2]);
288  histPlot->AddPlot(DcosCM, MCM, cat_list[NMX/2]);
289  histPlot->AddPlot(DcosTa, Eb_Ta, cat_list[NMX/2]);
290  histPlot->AddPlot(DcosXa, El_Xa, cat_list[NMX/2]);
291  histPlot->AddPlot(DcosCM, DcosTa, cat_list[NMX/2]);
292  histPlot->AddPlot(DcosTa, DcosTb, cat_list[NMX/2]);
293  histPlot->AddPlot(DcosXa, DcosXb, cat_list[NMX/2]);
294 
295  if(NMX > 1){
296  histPlot->AddPlot(MCM, RMX, cat_Reco);
297  histPlot->AddPlot(Eb_Ta, RMX, cat_Reco);
298  histPlot->AddPlot(El_Xa, RMX, cat_Reco);
299  }
300 
303 
304  for(int m = 0; m < NMX; m++){
305  g_Log << LogInfo << "Generating events for ";
306  g_Log << "(mXp - mN) / (mT - mN)";
307  g_Log << " = " << m+1 << " / " << NMX+1 << LogEnd;
308 
309  RMX = double(m+1)/double(NMX+1);
310  Xa_Gen.SetMass(RMX*mT + (1.-RMX)*mN);
311 
312  if(fix_MXb){
313  Xb_Gen.SetMass(mXb);
314  } else {
315  Xb_Gen.SetMass(RMX*mT + (1.-RMX)*mN);
316  }
317 
318  LAB_Gen.InitializeAnalysis();
319 
320  for(int igen = 0; igen < Ngen; igen++){
321  if(igen%((std::max(Ngen,10))/10) == 0)
322  g_Log << LogInfo << "Generating event " << igen << " of " << Ngen << LogEnd;
323 
324  // generate event
325  LAB_Gen.ClearEvent(); // clear the gen tree
326 
327  LAB_Gen.AnalyzeEvent(); // generate a new event
328 
329  TVector3 MET = LAB_Gen.GetInvisibleMomentum(); // Get the MET from gen tree
330  MET.SetZ(0.);
331 
332  // analyze event
333  LAB.ClearEvent(); // clear the event
334 
335  INV.SetLabFrameThreeVector(MET); // Set MET in reco tree
336 
337  La.SetLabFrameFourVector(La_Gen.GetFourVector()); // Set lepton 4-vectors
338  Lb.SetLabFrameFourVector(Lb_Gen.GetFourVector());
339 
340  B.AddLabFrameFourVector(Ba_Gen.GetFourVector()); // Set b-jet 4-vectors
341  B.AddLabFrameFourVector(Bb_Gen.GetFourVector());
342 
343  LAB.AnalyzeEvent(); // analyze the event
344 
346  // Observable Calculations
348 
349  double MCMgen = CM_Gen.GetMass();
350  double Eb_Tagen = Ba_Gen.GetFourVector(Ta_Gen).E();
351  double Eb_Tbgen = Bb_Gen.GetFourVector(Tb_Gen).E();
352  double El_Xagen = La_Gen.GetFourVector(Xa_Gen).E();
353  double El_Xbgen = Lb_Gen.GetFourVector(Xb_Gen).E();
354  double cosCMgen = CM_Gen.GetCosDecayAngle();
355  double cosTagen = Ta_Gen.GetCosDecayAngle();
356  double cosTbgen = Tb_Gen.GetCosDecayAngle();
357  double cosXagen = Xa_Gen.GetCosDecayAngle();
358  double cosXbgen = Xb_Gen.GetCosDecayAngle();
359 
360  MCM = CM.GetMass()/MCMgen;
361  Eb_Ta = Ba.GetFourVector(Ta).E()/Eb_Tagen;
362  Eb_Tb = Bb.GetFourVector(Tb).E()/Eb_Tbgen;
363  El_Xa = La.GetFourVector(Xa).E()/El_Xagen;
364  El_Xb = Lb.GetFourVector(Xb).E()/El_Xbgen;
365  cosCM = CM.GetCosDecayAngle();
366  cosTa = Ta.GetCosDecayAngle();
367  cosTb = Tb.GetCosDecayAngle();
368  cosXa = Xa.GetCosDecayAngle();
369  cosXb = Xb.GetCosDecayAngle();
370  DcosCM = asin(sqrt(1.-cosCM*cosCM)*cosCMgen-sqrt(1.-cosCMgen*cosCMgen)*cosCM);
371  DcosTa = asin(sqrt(1.-cosTa*cosTa)*cosTagen-sqrt(1.-cosTagen*cosTagen)*cosTa);
372  DcosTb = asin(sqrt(1.-cosTb*cosTb)*cosTbgen-sqrt(1.-cosTbgen*cosTbgen)*cosTb);
373  DcosXa = asin(sqrt(1.-cosXa*cosXa)*cosXagen-sqrt(1.-cosXagen*cosXagen)*cosXa);
374  DcosXb = asin(sqrt(1.-cosXb*cosXb)*cosXbgen-sqrt(1.-cosXbgen*cosXbgen)*cosXb);
375 
376  histPlot->Fill(cat_list[m]);
377  histPlot->Fill(cat_Reco);
378  }
379 
380  LAB_Gen.PrintGeneratorEfficiency();
381  }
382 
383  histPlot->Draw();
384 
385  TFile fout(output_name.c_str(),"RECREATE");
386  fout.Close();
387  histPlot->WriteOutput(output_name);
388  histPlot->WriteHist(output_name);
389  treePlot->WriteOutput(output_name);
390 
391  g_Log << LogInfo << "Finished" << LogEnd;
392 }
393 
394 # ifndef __CINT__ // main function for stand-alone compilation
395 int main(){
396  example_DiStop_to_bXp_bXm_to_blNblN();
397  return 0;
398 }
399 #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::HistPlot::WriteHist
void WriteHist(const std::string &filename)
Stores all histograms in root file.
Definition: HistPlot.cc:581
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::RFList
Definition: RFList.hh:43
RestFrames::HistPlot
Definition: HistPlot.hh:44
RestFrames::RFPlot::WriteOutput
void WriteOutput(const std::string &filename)
Stores all plots in root file.
Definition: RFPlot.cc:59
RestFrames.hh
RestFrames::ContraBoostInvJigsaw
Definition: ContraBoostInvJigsaw.hh:37
RestFrames::VisibleRecoFrame
Definition: VisibleRecoFrame.hh:41
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