LOGO

RestFrames  v1.0.1
RestFrames HEP Event Analysis Software Library
example_X2X2_to_ZllXHggX.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 using std::cout;
41 using std::endl;
42 
43 void example_X2X2_to_ZllXHggX(std::string output_name =
44  "output_X2X2_to_ZllXHggX.root"){
45 
46  // set particle masses and widths [GeV]
47  double mX2 = 300.;
48  double mZ = 91.19; // PDG 2016
49  double wZ = 2.50;
50  double mH = 125.;
51  double wH = 0.04;
52 
53  // number of different neutralino masses to evaluate
54  int NmX1 = 1;
55  std::vector<double> mX1;
56  mX1.push_back(250.); // lightest X1 mass to evaluate
57 
58  // Number of events to generate
59  int Ngen = 100000;
60 
62  g_Log << LogInfo << "Initializing generator frames and tree..." << LogEnd;
64  ppLabGenFrame LAB_Gen("LAB_Gen","LAB");
65  DecayGenFrame X2X2_Gen("X2X2_Gen","#tilde{#chi}^{ 0}_{2} #tilde{#chi}^{ 0}_{2}");
66  DecayGenFrame X2a_Gen("X2a_Gen","#tilde{#chi}^{ 0}_{2 a}");
67  DecayGenFrame X2b_Gen("X2b_Gen","#tilde{#chi}^{ 0}_{2 b}");
68  ResonanceGenFrame Za_Gen("Za_Gen","Z_{a}");
69  ResonanceGenFrame Hb_Gen("Hb_Gen","H_{b}");
70  VisibleGenFrame L1_Gen("L1_Gen","#it{l}_{1}");
71  VisibleGenFrame L2_Gen("L2_Gen","#it{l}_{2}");
72  VisibleGenFrame G1_Gen("G1_Gen","#gamma_{1}");
73  VisibleGenFrame G2_Gen("G2_Gen","#gamma_{2}");
74  InvisibleGenFrame X1a_Gen("X1a_Gen","#tilde{#chi}^{ 0}_{1 a}");
75  InvisibleGenFrame X1b_Gen("X1b_Gen","#tilde{#chi}^{ 0}_{1 b}");
76 
77  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
78 
79  LAB_Gen.SetChildFrame(X2X2_Gen);
80  X2X2_Gen.AddChildFrame(X2a_Gen);
81  X2X2_Gen.AddChildFrame(X2b_Gen);
82  X2a_Gen.AddChildFrame(Za_Gen);
83  X2a_Gen.AddChildFrame(X1a_Gen);
84  X2b_Gen.AddChildFrame(Hb_Gen);
85  X2b_Gen.AddChildFrame(X1b_Gen);
86  Za_Gen.AddChildFrame(L1_Gen);
87  Za_Gen.AddChildFrame(L2_Gen);
88  Hb_Gen.AddChildFrame(G1_Gen);
89  Hb_Gen.AddChildFrame(G2_Gen);
90 
91  if(LAB_Gen.InitializeTree())
92  g_Log << LogInfo << "...Successfully initialized generator tree" << LogEnd;
93  else
94  g_Log << LogError << "...Failed initializing generator tree" << LogEnd;
95 
96  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
97 
98  if(NmX1 < 1 || mX1.size() < 1)
99  return;
100  if(mX1[0] >= mX2)
101  return;
102 
103  X2X2_Gen.SetVariableMass();
104  X2a_Gen.SetMass(mX2); X2b_Gen.SetMass(mX2);
105  X1a_Gen.SetMass(mX1[0]); X1b_Gen.SetMass(mX1[0]);
106  Za_Gen.SetMass(mZ);
107  Za_Gen.SetWidth(wZ);
108  Hb_Gen.SetMass(mH);
109  Hb_Gen.SetWidth(wH);
110 
111  L1_Gen.SetPtCut(8.); L1_Gen.SetEtaCut(2.5);
112  L2_Gen.SetPtCut(8.); L2_Gen.SetEtaCut(2.5);
113  G1_Gen.SetPtCut(20.); G1_Gen.SetEtaCut(3.0);
114  G2_Gen.SetPtCut(20.); G2_Gen.SetEtaCut(3.0);
115 
116  if(LAB_Gen.InitializeAnalysis())
117  g_Log << LogInfo << "...Successfully initialized generator analysis" << LogEnd;
118  else
119  g_Log << LogError << "...Failed initializing generator analysis" << LogEnd;
122 
124  g_Log << LogInfo << "Initializing reconstruction frames and tree..." << LogEnd;
126  LabRecoFrame LAB("LAB","LAB");
127  DecayRecoFrame X2X2("X2X2","#tilde{#chi}^{ 0}_{2} #tilde{#chi}^{ 0}_{2}");
128  DecayRecoFrame X2a("X2a","#tilde{#chi}^{ 0}_{2 a}");
129  DecayRecoFrame X2b("X2b","#tilde{#chi}^{ 0}_{2 b}");
130  DecayRecoFrame Za("Za","Z_{a}");
131  DecayRecoFrame Hb("Hb","H_{b}");
132  VisibleRecoFrame L1("L1","#it{l}_{1}");
133  VisibleRecoFrame L2("L2","#it{l}_{2}");
134  VisibleRecoFrame G1("G1","#gamma_{1}");
135  VisibleRecoFrame G2("G2","#gamma_{2}");
136  InvisibleRecoFrame X1a("X1a","#tilde{#chi}^{ 0}_{1 a}");
137  InvisibleRecoFrame X1b("X1b","#tilde{#chi}^{ 0}_{1 b}");
138 
139  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
140 
141  LAB.SetChildFrame(X2X2);
142  X2X2.AddChildFrame(X2a);
143  X2X2.AddChildFrame(X2b);
144  X2a.AddChildFrame(Za);
145  X2a.AddChildFrame(X1a);
146  X2b.AddChildFrame(Hb);
147  X2b.AddChildFrame(X1b);
148  Za.AddChildFrame(L1);
149  Za.AddChildFrame(L2);
150  Hb.AddChildFrame(G1);
151  Hb.AddChildFrame(G2);
152 
153  if(LAB.InitializeTree())
154  g_Log << LogInfo << "...Successfully initialized reconstruction trees" << LogEnd;
155  else
156  g_Log << LogError << "...Failed initializing reconstruction trees" << LogEnd;
157 
158  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
159 
160  // Invisible Groups
161  InvisibleGroup INV("INV","#tilde{#chi}_{1}^{ 0} Jigsaws");
162  INV.AddFrames(X1a+X1b);
163 
164  // Set di-LSP mass to minimum Lorentz-invariant expression
165  // SetMassInvJigsaw X1_mass("X1_mass","Set M_{#tilde{#chi}_{1}^{ 0} #tilde{#chi}_{1}^{ 0}} to minimum");
166  // INV.AddJigsaw(X1_mass);
167 
168 
169  // Set di-LSP rapidity to that of visible particles
170  SetRapidityInvJigsaw X1_eta("X1_eta","#eta_{#tilde{#chi}_{1}^{ 0} #tilde{#chi}_{1}^{ 0}} = #eta_{2#gamma+2#it{l}}");
171  INV.AddJigsaw(X1_eta);
172  X1_eta.AddVisibleFrames(X2X2.GetListVisibleFrames());
173 
174  //ContraBoostInvJigsaw X1X1_minM2("X1X1_contra","Contraboost invariant Jigsaw");
175  //MinMassDiffInvJigsaw X1X1_contra("MinMh_R","min M_{h}, M_{h}^{ a}= M_{h}^{ b}",2);
176  MinMassesSqInvJigsaw X1X1_contra("MinMWa_R","min M_{W}, M_{W}^{a,a}= M_{W}^{a,b}", 2);
177  INV.AddJigsaw(X1X1_contra);
178  X1X1_contra.AddVisibleFrames(X2a.GetListVisibleFrames(), 0);
179  X1X1_contra.AddVisibleFrames(X2b.GetListVisibleFrames(), 1);
180  X1X1_contra.AddInvisibleFrame(X1a, 0);
181  X1X1_contra.AddInvisibleFrame(X1b, 1);
182 
183  if(LAB.InitializeAnalysis())
184  g_Log << LogInfo << "...Successfully initialized analysis" << std::endl << LogEnd;
185  else
186  g_Log << LogError << "...Failed initializing analysis" << LogEnd;
189 
190  TreePlot* treePlot = new TreePlot("TreePlot","TreePlot");
191 
192  treePlot->SetTree(LAB_Gen);
193  treePlot->Draw("GenTree", "Generator Tree", true);
194 
195  treePlot->SetTree(LAB);
196  treePlot->Draw("RecoTree", "Reconstruction Tree");
197 
198  treePlot->SetTree(INV);
199  treePlot->Draw("InvTree", "Inivisible Jigsaws", true);
200 
201  //-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//-//
202 
203  // Declare observables for histogram booking
204  HistPlot* histPlot = new HistPlot("Plots",
205  std::string("#tilde{#chi}_{2}^{ 0} #tilde{#chi}_{2}^{ 0}") +
206  "#rightarrow Z(#it{l}#it{l}) #tilde{#chi}_{1}^{ 0}"+
207  "H(#gamma #gamma) #tilde{#chi}_{1}^{ 0}");
208 
209  histPlot->SetRebin(1);
210 
212  char smassX2[200];
213  sprintf(smassX2, "m_{#tilde{#chi}^{0}_{2}}= %.0f", mX2);
214  for(int m = 0; m < NmX1; m++){
215  if(m > 0)
216  mX1.push_back(mX1[0] + m*(mX2-mX1[0])/double(NmX1));
217 
218  char smassX1[200], scat[50];
219  sprintf(scat, "MX1_%.0f", mX1[m]);
220  sprintf(smassX1, "m_{#tilde{#chi}^{0}_{1}}= %.0f", mX1[m]);
221 
222  cat_list += histPlot->GetNewCategory(scat, std::string(smassX2)+
223  " , "+std::string(smassX1));
224  }
225 
226  const HistPlotVar& MCM = histPlot->GetNewVar("MCM",
227  "M_{#tilde{#chi}_{2}^{ 0} #tilde{#chi}_{2}^{ 0}}",
228  0., 2.);
229  const HistPlotVar& EZX2a = histPlot->GetNewVar("EZX2a", "E_{Z}^{ #tilde{#chi}_{2 a}^{ 0}}", 0., 2.);
230  const HistPlotVar& EHX2b = histPlot->GetNewVar("EHX2b", "E_{H}^{ #tilde{#chi}_{2 b}^{ 0}}", 0., 2.);
231  const HistPlotVar& cosX2a = histPlot->GetNewVar("cosX2a","cos #theta_{#tilde{#chi}_{2 a}^{ 0}}", -1., 1.);
232  const HistPlotVar& cosX2b = histPlot->GetNewVar("cosX2b","cos #theta_{#tilde{#chi}_{2 b}^{ 0}}", -1., 1.);
233  const HistPlotVar& cosZ = histPlot->GetNewVar("cosZ","cos #theta_{Z}", -1., 1.);
234  const HistPlotVar& cosH = histPlot->GetNewVar("cosH","cos #theta_{H}", -1., 1.);
235  const HistPlotVar& DcosZ = histPlot->GetNewVar("DcosZ","#theta_{Z} - #theta_{Z}^{gen}", -1., 1.);
236  const HistPlotVar& DcosH = histPlot->GetNewVar("DcosH","#theta_{H} - #theta_{H}^{gen}", -1., 1.);
237  const HistPlotVar& DcosX2a = histPlot->GetNewVar("DcosX2a","#theta_{X2a} - #theta_{X2a}^{gen}", -1., 1.);
238  const HistPlotVar& DcosX2b = histPlot->GetNewVar("DcosX2b","#theta_{X2b} - #theta_{X2b}^{gen}", -1., 1.);
239  const HistPlotVar& RISR = histPlot->GetNewVar("RISR","R_{ISR}", 0.3, 1.);
240  const HistPlotVar& PTISR = histPlot->GetNewVar("PTISR","p_{T}^{ ISR}", 0., 1000., "[GeV]");
241 
242  const HistPlotVar& MI2 = histPlot->GetNewVar("MI2","m_{I}^{2}", -2000., 2000., "[GeV^{2}]");
243  const HistPlotVar& MI = histPlot->GetNewVar("MI","m_{I}", 0., 800., "[GeV]");
244  const HistPlotVar& MP = histPlot->GetNewVar("MP","m_{P}", 0., 800., "[GeV]");
245  const HistPlotVar& RISR2 = histPlot->GetNewVar("RISR2","R_{ISR2}", 0.3, 1.);
246 
247  /*
248  histPlot->AddPlot(MCM, cat_list);
249  histPlot->AddPlot(EZX2a, cat_list);
250  histPlot->AddPlot(EHX2b, cat_list);
251  histPlot->AddPlot(cosZ, cat_list);
252  histPlot->AddPlot(cosH, cat_list);
253  histPlot->AddPlot(DcosZ, cat_list);
254  histPlot->AddPlot(DcosH, cat_list);
255  histPlot->AddPlot(DcosX2a, cat_list);
256  histPlot->AddPlot(DcosX2b, cat_list);
257  histPlot->AddPlot(MCM, EZX2a, cat_list[NmX1/2]);
258  histPlot->AddPlot(MCM, EHX2b, cat_list[NmX1/2]);
259  histPlot->AddPlot(EZX2a, EHX2b, cat_list[NmX1/2]);
260  histPlot->AddPlot(EZX2a, DcosX2a, cat_list[NmX1/2]);
261  histPlot->AddPlot(EHX2b, DcosX2b, cat_list[NmX1/2]);
262  histPlot->AddPlot(cosZ, DcosZ, cat_list[NmX1/2]);
263  histPlot->AddPlot(cosH, DcosH, cat_list[NmX1/2]);
264  */
265 
266  histPlot->AddPlot(RISR, cat_list);
267  histPlot->AddPlot(RISR, PTISR, cat_list[NmX1/2]);
268  //histPlot->AddPlot(RISR, EZX2a, cat_list[NmX1/2]);
269 
270  histPlot->AddPlot(MI, cat_list);
271  histPlot->AddPlot(MP, cat_list);
272  histPlot->AddPlot(RISR2, cat_list);
273  histPlot->AddPlot(MI, PTISR, cat_list);
274  histPlot->AddPlot(MP, PTISR, cat_list);
275  histPlot->AddPlot(MI, RISR, cat_list);
276  histPlot->AddPlot(MP, RISR, cat_list);
277  histPlot->AddPlot(MI, RISR2, cat_list);
278  histPlot->AddPlot(MP, RISR2, cat_list);
279  histPlot->AddPlot(MP, MI, cat_list);
280  histPlot->AddPlot(RISR2, PTISR, cat_list[NmX1/2]);
281  histPlot->AddPlot(RISR2, RISR, cat_list[NmX1/2]);
282 
283  for(int m = 0; m < NmX1; m++){
284  g_Log << LogInfo << "Generating events for ";
285  g_Log << "mX2 = " << mX2 << " , ";
286  g_Log << "mX1 = " << mX1[m] << LogEnd;
287 
288  X1a_Gen.SetMass(mX1[m]);
289  X1b_Gen.SetMass(mX1[m]);
290 
291  LAB_Gen.InitializeAnalysis();
292 
293  for(int igen = 0; igen < Ngen; igen++){
294  if(igen%((std::max(Ngen,10))/10) == 0)
295  g_Log << LogInfo << "Generating event " << igen << " of " << Ngen << LogEnd;
296 
297  // generate event
298  LAB_Gen.ClearEvent(); // clear the gen tree
299 
300  //LAB_Gen.SetTransverseMomentum(600.); // give X2X2 some Pt
301  LAB_Gen.SetTransverseMomentum(800.*gRandom->Rndm()); // give X2X2 some Pt
302 
303  LAB_Gen.AnalyzeEvent(); // generate a new event
304 
305  // testing
306  TLorentzVector P_V = X2X2_Gen.GetListVisibleFrames().GetFourVector();
307  TLorentzVector P_I = X2X2_Gen.GetListInvisibleFrames().GetFourVector();
308 
309  // TVector3 BetaZ = (P_V+P_I).BoostVector();
310  // BetaZ.SetX(0.);
311  // BetaZ.SetY(0.);
312 
313  // P_V.Boost(-BetaZ);
314  // P_V.Boost(-BetaZ);
315 
316  P_V.SetPxPyPzE(P_V.Px(),P_V.Py(),0.,sqrt(P_V.Pt()*P_V.Pt() + P_V.M2()));
317  P_I.SetPxPyPzE(P_I.Px(),P_I.Py(),0.,sqrt(P_I.Pt()*P_I.Pt() + P_I.M2()));
318 
319  // double Ei = P_V.E()*(P_I.P()*P_I.P() + P_V.Vect().Dot(P_I.Vect()))/(P_V.P()*P_V.P() + P_V.Vect().Dot(P_I.Vect()));
320 
321  // if(Ei < 0){
322  // MI2 = 0.;
323  // MI = 0;
324  // MP = 0;
325  // continue;
326  // }
327 
328  // MI2 = (pow(P_V.E()*(P_I.P()*P_I.P() + P_V.Vect().Dot(P_I.Vect()))/(P_V.P()*P_V.P() + P_V.Vect().Dot(P_I.Vect())),2.) -
329  // P_I.P()*P_I.P() - P_V.M2() +
330  // pow(X2a_Gen.GetListVisibleFrames().GetFourVector().M()+X2b_Gen.GetListVisibleFrames().GetFourVector().M(),2.)) / 4.;
331 
332  // MI = sqrt( MI2 );
333 
334  double a = ((P_V+P_I).P()*(P_V+P_I).P()*(P_V+P_I).P()*(P_V+P_I).P()-P_I.Vect().Dot((P_V+P_I).Vect())*P_I.Vect().Dot((P_V+P_I).Vect()));
335  double b = -P_V.E()*P_I.Vect().Dot((P_V+P_I).Vect())*P_I.Vect().Dot((P_V+P_I).Vect());
336  double c1 = - P_V.E()*P_V.E()*P_I.Vect().Dot((P_V+P_I).Vect())*P_I.Vect().Dot((P_V+P_I).Vect());
337  double c2 = (4.-1.)*(P_V.Vect().Cross(P_I.Vect()).Mag2())*(P_V+P_I).P()*(P_V+P_I).P();
338  double Ei = (-b + sqrt(b*b - a*(c1+c2))) / a;
339  if(b*b - a*(c1+c2) < 0.)
340  cout << " wtf " << endl;
341  if(Ei < 0){
342  cout << "HERE" << endl;
343  MI2 = 0.;
344  MI = 0;
345  MP = 0;
346  continue;
347  }
348 
349  MI2 = (Ei*Ei - P_I.P()*P_I.P() - P_V.M2() +
350  pow(X2a_Gen.GetListVisibleFrames().GetFourVector().M()+X2b_Gen.GetListVisibleFrames().GetFourVector().M(),2.)) / 4.;
351 
352  MI = sqrt( MI2 );
353 
354  //MP = MI / RISR;
355 
356  double Minv = sqrt(Ei*Ei - P_I.P()*P_I.P());
357 
358  // analyze event
359  LAB.ClearEvent(); // clear the reco tree
360  L1.SetLabFrameFourVector(L1_Gen.GetFourVector()); // Set lepton 4-vectors
361  L2.SetLabFrameFourVector(L2_Gen.GetFourVector());
362  G1.SetLabFrameFourVector(G1_Gen.GetFourVector()); // Set photon 4-vectors
363  G2.SetLabFrameFourVector(G2_Gen.GetFourVector());
364 
365  TVector3 vMET = LAB_Gen.GetInvisibleMomentum();
366  vMET.SetZ(0.);
367  TLorentzVector MET;
368  //MET.SetVectM(vMET, Minv);
369  MET.SetVectM(vMET, Minv); // Get the MET from gen tree
370  INV.SetLabFrameFourVector(MET); // Set the MET in reco tree
371 
372  X1a.SetMinimumMass(MI);
373  X1b.SetMinimumMass(MI);
374 
375  LAB.AnalyzeEvent(); //analyze the event
376 
377  // calculate observables
378  MCM = X2X2.GetMass() / X2X2_Gen.GetMass();
379 
380  EZX2a = Za.GetEnergy(X2a) / Za_Gen.GetEnergy(X2a_Gen);
381  EHX2b = Hb.GetEnergy(X2b) / Hb_Gen.GetEnergy(X2b_Gen);
382 
383  cosX2a = X2a.GetCosDecayAngle();
384  double cosX2agen = X2a_Gen.GetCosDecayAngle();
385  cosX2b = X2b.GetCosDecayAngle();
386  double cosX2bgen = X2b_Gen.GetCosDecayAngle();
387  cosZ = Za.GetCosDecayAngle();
388  double cosZgen = Za_Gen.GetCosDecayAngle();
389  cosH = Hb.GetCosDecayAngle();
390  double cosHgen = Hb_Gen.GetCosDecayAngle();
391 
392  DcosX2a = asin(sqrt(1.-cosX2a*cosX2a)*cosX2agen-sqrt(1.-cosX2agen*cosX2agen)*cosX2a);
393  DcosX2b = asin(sqrt(1.-cosX2b*cosX2b)*cosX2bgen-sqrt(1.-cosX2bgen*cosX2bgen)*cosX2b);
394  DcosZ = asin(sqrt(1.-cosZ*cosZ)*cosZgen-sqrt(1.-cosZgen*cosZgen)*cosZ);
395  DcosH = asin(sqrt(1.-cosH*cosH)*cosHgen-sqrt(1.-cosHgen*cosHgen)*cosH);
396 
397  TVector3 vP_ISR = X2X2.GetFourVector(LAB).Vect();
398  TVector3 vP_I = X2X2.GetListInvisibleFrames().GetFourVector(LAB).Vect();
399  vP_ISR.SetZ(0.);
400  vP_I.SetZ(0.);
401 
402  RISR = fabs(vP_I.Dot(vP_ISR.Unit())) / vP_ISR.Mag();
403  PTISR = vP_ISR.Mag();
404 
405 
406  MP = X2a.GetMass();
407 
408  RISR2 = MI/MP;
409 
410  histPlot->Fill(cat_list[m]);
411  }
412 
413  LAB_Gen.PrintGeneratorEfficiency();
414  }
415 
416  histPlot->Draw();
417 
418  TFile fout(output_name.c_str(),"RECREATE");
419  fout.Close();
420  histPlot->WriteOutput(output_name);
421  histPlot->WriteHist(output_name);
422  treePlot->WriteOutput(output_name);
423 
424  g_Log << LogInfo << "Finished" << LogEnd;
425 }
426 
427 # ifndef __CINT__ // main function for stand-alone compilation
428  int main(){
429  example_X2X2_to_ZllXHggX();
430  return 0;
431  }
432 #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::HistPlot::SetRebin
void SetRebin(int rebin=4)
Reduces the number of bins by a factor.
Definition: HistPlot.cc:576
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::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::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::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