LOGO

RestFrames  v1.0.1
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 
49  GeneratorFrame(const std::string& sname, const std::string& stitle);
50 
53  virtual ~GeneratorFrame();
54 
56  virtual void Clear();
57 
68  virtual void AddChildFrame(RestFrame& frame);
69 
76  virtual void SetParentFrame(RestFrame& frame =
78 
85  virtual GeneratorFrame const& GetParentFrame() const;
86 
92  virtual GeneratorFrame& GetChildFrame(int i = 0) const;
93 
102  void SetPCut(double cut);
103 
112  void SetPtCut(double cut);
113 
120  void SetEtaCut(double cut);
121 
133  void SetMassWindowCut(double min, double max);
134 
136  void RemovePCut();
137 
139  void RemovePtCut();
140 
142  void RemoveEtaCut();
143 
145  void RemoveMassWindowCut();
146 
148  void PrintGeneratorEfficiency() const;
149 
151  virtual double GetMass() const;
152 
161 
163  bool IsVariableMassMCMC() const;
164 
166  virtual double GetMinimumMassMCMC() const;
167 
169  virtual void GenerateMassMCMC(double& mass, double& prob,
170  double max = -1.) const;
171 
173  virtual double GetProbMCMC(double mass = -1.) const;
174 
178  double GetRandom() const;
179 
184  double GetGaus(double mu, double sig) const;
185 
187 
191  static GeneratorFrame& Empty();
192 
193  protected:
194  double m_Mass;
195 
197  bool AnalyzeEventRecursive();
198  bool ClearEventRecursive();
199 
200  virtual void ResetGenFrame() = 0;
201  virtual bool GenerateFrame() = 0;
202 
203  void SetChildren(const std::vector<TLorentzVector>& P_children);
204  virtual bool InitializeGenAnalysis();
205 
206  virtual bool IterateMCMC();
207  bool IterateRecursiveMCMC();
208 
209  void SetVariableMassMCMC(bool var = true);
210  virtual void SetMassMCMC(double mass);
211  void SetMassMCMC(double mass, GeneratorFrame& frame) const;
212 
213  bool EventInAcceptance() const;
214 
215  private:
216  TRandom *m_Random;
217  bool m_VarMassMCMC;
218 
219  mutable long m_Ngen;
220  mutable long m_Npass;
221 
222  double m_PCut;
223  double m_PtCut;
224  double m_EtaCut;
225  double m_minMassCut;
226  double m_maxMassCut;
227 
228  bool m_doCuts;
229  bool m_doPCut;
230  bool m_doPtCut;
231  bool m_doEtaCut;
232  bool m_dominMassCut;
233  bool m_domaxMassCut;
234 
235  };
236 
237 }
238 
239 #endif
RestFrames::GeneratorFrame::AddChildFrame
virtual void AddChildFrame(RestFrame &frame)
Adds a child RestFrame to this frame.
Definition: GeneratorFrame.cc:92
RestFrames::GeneratorFrame::GetParentFrame
virtual GeneratorFrame const & GetParentFrame() const
Returns the parent of this frame.
Definition: GeneratorFrame.cc:104
RestFrames::GeneratorFrame::PrintGeneratorEfficiency
void PrintGeneratorEfficiency() const
Print generator efficiency information.
Definition: GeneratorFrame.cc:273
RestFrames::RestFrame
abstract base class for all Frame objects
Definition: RestFrame.hh:45
RestFrames::GeneratorFrame::GenerateMassMCMC
virtual void GenerateMassMCMC(double &mass, double &prob, double max=-1.) const
Generates mass for Markov Chain MonteCarlo event generation.
Definition: GeneratorFrame.cc:247
RestFrames::GeneratorFrame::GetGaus
double GetGaus(double mu, double sig) const
Returns gaussian random number.
Definition: GeneratorFrame.cc:197
RestFrames::GeneratorFrame
Definition: GeneratorFrame.hh:42
RestFrames::GeneratorFrame::SetEtaCut
void SetEtaCut(double cut)
Sets pseudorapidity cut.
Definition: GeneratorFrame.cc:372
RestFrames::GeneratorFrame::GetRandom
double GetRandom() const
Returns random value.
Definition: GeneratorFrame.cc:193
RestFrames::GeneratorFrame::RemovePtCut
void RemovePtCut()
Removes transverse momentum cut.
Definition: GeneratorFrame.cc:401
RestFrames::GeneratorFrame::Empty
static GeneratorFrame & Empty()
Returns empty GeneratorFrame.
Definition: GeneratorFrame.cc:84
RestFrames::GeneratorFrame::GetMinimumMassMCMC
virtual double GetMinimumMassMCMC() const
Returns minimun mass of Markov Chain MonteCarlo event generation.
Definition: GeneratorFrame.cc:230
RestFrames::GeneratorFrame::RemoveEtaCut
void RemoveEtaCut()
Removes pseudorapidity cut.
Definition: GeneratorFrame.cc:407
RestFrames::GeneratorFrame::AnalyzeEventRecursive
bool AnalyzeEventRecursive()
Recursively analyze event in this frame and its children.
Definition: GeneratorFrame.cc:133
RestFrames::GeneratorFrame::IsVariableMassMCMC
bool IsVariableMassMCMC() const
Is this frame capable having a variable mass? (true/false)
Definition: GeneratorFrame.cc:226
RestFrames::GeneratorFrame::RemovePCut
void RemovePCut()
Removes momentum cut.
Definition: GeneratorFrame.cc:395
RestFrames::GeneratorFrame::GetChildFrame
virtual GeneratorFrame & GetChildFrame(int i=0) const
Returns the frame of the i th child.
Definition: GeneratorFrame.cc:112
RestFrames::GeneratorFrame::SetMassWindowCut
void SetMassWindowCut(double min, double max)
Sets mass frame cut.
Definition: GeneratorFrame.cc:378
RestFrames::GeneratorFrame::GetProbMCMC
virtual double GetProbMCMC(double mass=-1.) const
Evaluates probability of the state is in right now or the probablility of a state it could go in.
Definition: GeneratorFrame.cc:269
RestFrames::GeneratorFrame::GeneratorFrame
GeneratorFrame()
Empty constructor.
Definition: GeneratorFrame.cc:72
RestFrames::GeneratorFrame::SetPCut
void SetPCut(double cut)
Sets momentum cut.
Definition: GeneratorFrame.cc:358
RestFrames::RestFrame::Empty
static RestFrame & Empty()
Returns empty RestFrame.
Definition: RestFrame.cc:62
RestFrames::GeneratorFrame::SetPtCut
void SetPtCut(double cut)
Sets transverse momentum cut.
Definition: GeneratorFrame.cc:365
RestFrames::GeneratorFrame::InitializeAnalysisRecursive
bool InitializeAnalysisRecursive()
Recursively initialize this frame and its children for analysis.
Definition: GeneratorFrame.cc:169
RestFrames::GeneratorFrame::ClearEventRecursive
bool ClearEventRecursive()
Recursively clear event information from this frame and its children.
Definition: GeneratorFrame.cc:120
RestFrames::GeneratorFrame::GetMass
virtual double GetMass() const
Get the mass of this frame.
Definition: GeneratorFrame.cc:88
RestFrames::GeneratorFrame::SetParentFrame
virtual void SetParentFrame(RestFrame &frame=RestFrame::Empty())
Set the parent frame for this frame.
Definition: GeneratorFrame.cc:98
RestFrames::GeneratorFrame::Clear
virtual void Clear()
Clears GeneratorFrame of all connections to other objects.
Definition: GeneratorFrame.cc:79
RestFrames::GeneratorFrame::RemoveMassWindowCut
void RemoveMassWindowCut()
Removes mass frame cuts.
Definition: GeneratorFrame.cc:413
RestFrame.hh