LOGO

RestFrames  v1.0.0
RestFrames HEP Event Analysis Software Library
GeneratorFrame.hh
Go to the documentation of this file.
1 // RestFrames: particle physics event analysis library
3 // --------------------------------------------------------------------
4 // Copyright (c) 2014-2016, Christopher Rogan
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 #ifndef GeneratorFrame_HH
31 #define GeneratorFrame_HH
32 
33 #include <TRandom.h>
34 
35 #include "RestFrames/RestFrame.hh"
36 
37 namespace RestFrames {
38 
40  // GeneratorFrame class
42  class GeneratorFrame : public RestFrame {
43  public:
44  GeneratorFrame(const std::string& sname, const std::string& stitle);
46  virtual ~GeneratorFrame();
47 
49  virtual void Clear();
50 
58  virtual void AddChildFrame(RestFrame& frame);
59 
66  virtual void SetParentFrame(RestFrame& frame =
67  RestFrame::Empty());
68 
74  virtual GeneratorFrame const& GetParentFrame() const;
75 
77  virtual GeneratorFrame& GetChildFrame(int i = 0) const;
78 
79  void SetPCut(double cut);
80  void SetPtCut(double cut);
81  void SetEtaCut(double cut);
82  void SetMassWindowCut(double min, double max);
83 
84  void RemovePCut();
85  void RemovePtCut();
86  void RemoveEtaCut();
87  void RemoveMassWindowCut();
88 
90  void PrintGeneratorEfficiency() const;
91 
93  virtual double GetMass() const;
94 
96  bool IsVariableMassMCMC() const;
97 
98  virtual double GetMinimumMassMCMC() const;
99  virtual void GenerateMassMCMC(double& mass, double& prob,
100  double max = -1.) const;
101  virtual double GetProbMCMC(double mass = -1.) const;
102 
103  double GetRandom() const;
104  double GetGaus(double mu, double sig) const;
105 
106  static GeneratorFrame& Empty();
107 
108  protected:
109  double m_Mass;
110 
111  bool InitializeAnalysisRecursive();
112  bool AnalyzeEventRecursive();
113  bool ClearEventRecursive();
114 
115  virtual void ResetGenFrame() = 0;
116  virtual bool GenerateFrame() = 0;
117 
118  void SetChildren(const std::vector<TLorentzVector>& P_children);
119  virtual bool InitializeGenAnalysis();
120 
121  virtual bool IterateMCMC();
122  bool IterateRecursiveMCMC();
123 
124  void SetVariableMassMCMC(bool var = true);
125  virtual void SetMassMCMC(double mass);
126  void SetMassMCMC(double mass, GeneratorFrame& frame) const;
127 
128  bool EventInAcceptance() const;
129 
130  private:
131  TRandom *m_Random;
132  bool m_VarMassMCMC;
133 
134  mutable long m_Ngen;
135  mutable long m_Npass;
136 
137  double m_PCut;
138  double m_PtCut;
139  double m_EtaCut;
140  double m_minMassCut;
141  double m_maxMassCut;
142 
143  bool m_doCuts;
144  bool m_doPCut;
145  bool m_doPtCut;
146  bool m_doEtaCut;
147  bool m_dominMassCut;
148  bool m_domaxMassCut;
149 
150  };
151 
152 }
153 
154 #endif
virtual double GetMass() const
Get the mass of this frame.
virtual GeneratorFrame const & GetParentFrame() const
Returns the parent of this frame.
bool IsVariableMassMCMC() const
Frame is capable having a variable mass? (true/false)
virtual void AddChildFrame(RestFrame &frame)
Add a child RestFrame to this frame.
virtual void SetParentFrame(RestFrame &frame=RestFrame::Empty())
Set the parent frame for this frame.
virtual void Clear()
Clears GeneratorFrame of all connections to other objects.
abstract base class for all Frame objects
static GeneratorFrame & Empty()
Returns empty instance of class.
virtual GeneratorFrame & GetChildFrame(int i=0) const
Get the frame of the i th child.
void PrintGeneratorEfficiency() const
Print generator efficiency information.