RestFrames
v1.0.1
RestFrames HEP Event Analysis Software Library
|
Go to the documentation of this file.
30 #ifndef GeneratorFrame_HH
31 #define GeneratorFrame_HH
37 namespace RestFrames {
49 GeneratorFrame(
const std::string& sname,
const std::string& stitle);
151 virtual double GetMass()
const;
170 double max = -1.)
const;
173 virtual double GetProbMCMC(
double mass = -1.)
const;
184 double GetGaus(
double mu,
double sig)
const;
200 virtual void ResetGenFrame() = 0;
201 virtual bool GenerateFrame() = 0;
203 void SetChildren(
const std::vector<TLorentzVector>& P_children);
204 virtual bool InitializeGenAnalysis();
206 virtual bool IterateMCMC();
207 bool IterateRecursiveMCMC();
209 void SetVariableMassMCMC(
bool var =
true);
210 virtual void SetMassMCMC(
double mass);
213 bool EventInAcceptance()
const;
220 mutable long m_Npass;
virtual void AddChildFrame(RestFrame &frame)
Adds a child RestFrame to this frame.
virtual GeneratorFrame const & GetParentFrame() const
Returns the parent of this frame.
void PrintGeneratorEfficiency() const
Print generator efficiency information.
abstract base class for all Frame objects
virtual void GenerateMassMCMC(double &mass, double &prob, double max=-1.) const
Generates mass for Markov Chain MonteCarlo event generation.
double GetGaus(double mu, double sig) const
Returns gaussian random number.
void SetEtaCut(double cut)
Sets pseudorapidity cut.
double GetRandom() const
Returns random value.
void RemovePtCut()
Removes transverse momentum cut.
static GeneratorFrame & Empty()
Returns empty GeneratorFrame.
virtual double GetMinimumMassMCMC() const
Returns minimun mass of Markov Chain MonteCarlo event generation.
void RemoveEtaCut()
Removes pseudorapidity cut.
bool AnalyzeEventRecursive()
Recursively analyze event in this frame and its children.
bool IsVariableMassMCMC() const
Is this frame capable having a variable mass? (true/false)
void RemovePCut()
Removes momentum cut.
virtual GeneratorFrame & GetChildFrame(int i=0) const
Returns the frame of the i th child.
void SetMassWindowCut(double min, double max)
Sets mass frame cut.
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.
GeneratorFrame()
Empty constructor.
void SetPCut(double cut)
Sets momentum cut.
static RestFrame & Empty()
Returns empty RestFrame.
void SetPtCut(double cut)
Sets transverse momentum cut.
bool InitializeAnalysisRecursive()
Recursively initialize this frame and its children for analysis.
bool ClearEventRecursive()
Recursively clear event information from this frame and its children.
virtual double GetMass() const
Get the mass of 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.
void RemoveMassWindowCut()
Removes mass frame cuts.