RestFrames
v1.0.1
RestFrames HEP Event Analysis Software Library
|
abstract base class for all Frame objects
The RestFrame class is the base class for all frames objects.
Definition at line 45 of file RestFrame.hh.
Public Member Functions | |
RestFrame (const std::string &sname, const std::string &stitle) | |
Standard constructor. More... | |
RestFrame () | |
Empty constructor. | |
virtual void | Clear () |
Clears RestFrame of all connections to other objects. | |
RestFrame type methods | |
FrameType | GetType () const |
Returns RestFrame (FrameType) type. | |
bool | IsVisibleFrame () const |
Is this a VisibleFrame ? (yes/no) | |
bool | IsInvisibleFrame () const |
Is this an InvisibleFrame ? (yes/no) | |
bool | IsDecayFrame () const |
Is this a DecayFrame ? (yes/no) | |
bool | IsLabFrame () const |
Is this a LabFrame ? (yes/no) | |
bool | IsRecoFrame () const |
Is this an ReconstructionFrame ? (yes/no) | |
bool | IsGenFrame () const |
Is this a GeneratorFrame ? (yes/no) | |
virtual std::string | PrintString (LogType type) const |
String of information about RestFrame. | |
RestFrame tree construction methods | |
Member functions for assembling/disassembling trees of connected RestFrame | |
virtual void | AddChildFrame (RestFrame &frame) |
Adds a child RestFrame to this frame. More... | |
void | AddChildFrames (const RestFrameList &frames) |
Adds a list of children to this frame. More... | |
virtual void | SetParentFrame (RestFrame &frame=RestFrame::Empty()) |
Sets the parent frame for this frame. More... | |
void | RemoveChildFrame (RestFrame &frame) |
Removes a child of this frame. More... | |
void | RemoveChildFrames () |
Removes all the children of this frame. More... | |
RestFrame tree structure methods | |
Member functions which can be used to access RestFrames connected to this frame through parent(s) or children. | |
int | GetNChildren () const |
Returns the number of child frames inheriting from this one. | |
int | GetNDescendants () const |
Returns the number of descendents of this frame. | |
virtual RestFrame const & | GetParentFrame () const |
Returns the parent of this frame. More... | |
virtual RestFrame & | GetChildFrame (int i=0) const |
Get the RestFrame of the i th child. More... | |
RestFrameList const & | GetChildFrames () const |
Returns a list of this frame's child RestFrames. More... | |
virtual RestFrame const & | GetLabFrame () const |
Returns the LabFrame that this frame inherits from. More... | |
virtual RestFrame const & | GetProductionFrame () const |
Returns the production frame of this frame. More... | |
virtual RestFrame const & | GetSiblingFrame () const |
Returns the sibling frame of this frame. More... | |
int | GetFrameDepth (const RestFrame &frame) const |
Returns the depth of frame More... | |
virtual RestFrame const & | GetFrameAtDepth (int depth, const RestFrame &frame) const |
Returns the frame at depth. More... | |
virtual RestFrameList | GetListFrames (FrameType type=kLabFrame) const |
Returns a list of frames inheriting from this one. More... | |
virtual RestFrameList | GetListVisibleFrames () const |
Returns a list of VisibleFrame s inheriting from this. | |
virtual RestFrameList | GetListInvisibleFrames () const |
Returns a list of InvisibleFrame s inheriting from this. | |
RestFrame event analysis functions | |
Member functions which can be used to analyze an event. Each these functions generally requires that the method "AnalyzeEvent()" be successfully called from a LabFrame which is connected to this frame. Otherwise, a trivial value will be returned and a warning message printed. | |
RestFrameList | operator+ (RestFrame &frame) |
Combines RestFrames into RestFrameList. More... | |
RestFrameList | operator+ (const RestFrameList &frames) |
Combines RestFrames into RestFrameList. More... | |
virtual RFCharge | GetCharge () const |
Returns the charge of this frame. More... | |
virtual double | GetMass () const |
Returns the mass of this frame. More... | |
TLorentzVector | GetFourVector (const RestFrame &frame=RestFrame::Empty()) const |
Returns this frame's four-vector in a specified frame. More... | |
TLorentzVector | GetTransverseFourVector (const RestFrame &frame=RestFrame::Empty(), const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &axis_frame=RestFrame::Empty()) const |
Returns this frame's four-vector in a specified frame. More... | |
TLorentzVector | GetVisibleFourVector (const RestFrame &frame=RestFrame::Empty()) const |
Returns four-vector of visible descendants. More... | |
TLorentzVector | GetInvisibleFourVector (const RestFrame &frame=RestFrame::Empty()) const |
Returns four-vector of invisible descendants. More... | |
double | GetEnergy (const RestFrame &frame) const |
Returns energy of this frame in specified reference frame. More... | |
double | GetMomentum (const RestFrame &frame) const |
Returns magnitude of momentum. More... | |
TLorentzVector | GetFourVector (const TLorentzVector &P, const RestFrame &def_frame=RestFrame::Empty()) const |
Returns four vector boosted to different frame. More... | |
TLorentzVector | GetTransverseFourVector (const TLorentzVector &P, const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &axis_frame=RestFrame::Empty()) const |
Returns transverse four vector in this frame. | |
double | GetTransverseMomentum (const RestFrame &frame=RestFrame::Empty(), const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &axis_frame=RestFrame::Empty()) const |
Returns magnitude of transverse momentum. More... | |
double | GetTransverseMomentum (const TLorentzVector &P, const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &axis_frame=RestFrame::Empty()) const |
Returns magnitude of transverse momentum of this frame. More... | |
double | GetVisibleShape () const |
Returns visible shape of frame. More... | |
double | GetSumVisibleMomentum () const |
Returns scalar sum of visible child momenta. More... | |
double | GetSumInvisibleMomentum () const |
Returns scalar sum of invisible child momenta. More... | |
TVector3 | GetBoostInParentFrame () const |
Returns the boost of this frame in it's parent's frame. More... | |
double | GetGammaInParentFrame () const |
Returns the gamma of this frame in its parent's frame. More... | |
TVector3 | GetDecayPlaneNormalVector (const RestFrame &frame=RestFrame::Empty()) const |
Returns the vector normal to the decay plane of this frame. More... | |
double | GetDeltaPhiDecayPlanes (const RestFrame &frame) const |
Returns difference of azimuthal angle between decay planes. More... | |
double | GetCosDecayAngle (const RestFrame &frame=RestFrame::Empty()) const |
Returns the cosine of this frame's decay angle. More... | |
double | GetDeltaPhiDecayAngle (const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &frame=RestFrame::Empty()) const |
Returns difference of azimuthal decay angles in an axis. More... | |
double | GetDeltaPhiBoostVisible (const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &frame=RestFrame::Empty()) const |
Returns difference of azimuthal angles between the decay angles of the visible particles decaying from this frame. More... | |
double | GetDeltaPhiDecayVisible (const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &frame=RestFrame::Empty()) const |
Returns difference of azimuthal angles between the decay angles of the visible particles decaying from this frame. More... | |
double | GetDeltaPhiVisible (const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &frame=RestFrame::Empty()) const |
Returns difference of azimuthal angles between the angles of its visible particles kids. More... | |
Public Member Functions inherited from RestFrames::RFBase | |
RFBase (const std::string &sname, const std::string &stitle, int key) | |
Standard constructor. More... | |
RFBase () | |
Empty constructor. | |
bool | IsEmpty () const |
Checks whether this is default (empty) instance of class. | |
bool | operator! () const |
Tests whether key is the same as this. | |
void | Print (LogType type) const |
Print information associated with object. | |
RFKey | GetKey () const |
gets object identification key | |
std::string | GetName () const |
Returns object name. | |
std::string | GetTitle () const |
Returns object title. | |
bool | IsSame (const RFKey &key) const |
Tests whether key is the same as this. | |
bool | IsSame (const RFBase &obj) const |
Tests whether obj is the same as this. | |
bool | operator== (const RFKey &key) const |
Tests whether key is the same as this. | |
bool | operator== (const RFBase &obj) const |
Tests whether obj is the same as this. | |
bool | operator!= (const RFKey &key) const |
Tests whether key is the same as this. | |
bool | operator!= (const RFBase &obj) const |
Tests whether obj is the same as this. | |
Static Public Member Functions | |
static void | SetAxis (const TVector3 &axis) |
Set axis perpendicular to transverse plane. More... | |
static TVector3 const & | GetAxis () |
Retrieve axis which defines transverse plane. More... | |
static RestFrame & | Empty () |
Returns empty RestFrame. More... | |
static ConstRestFrameList const & | EmptyList () |
Returns empty RestFrameList. More... | |
Static Public Member Functions inherited from RestFrames::RFBase | |
static RFBase & | Empty () |
Returns empty RFBase. More... | |
Protected Member Functions | |
virtual bool | IsSoundBody () const |
TVector3 const & | GetChildBoostVector (RestFrame &frame) const |
TVector3 const & | GetParentBoostVector () const |
virtual bool | InitializeTreeRecursive () |
Recursively initialize this frame's tree. | |
virtual bool | InitializeAnalysisRecursive ()=0 |
Recursively initialize this frame and its children for analysis. | |
virtual bool | AnalyzeEventRecursive ()=0 |
Recursively analyze event in this frame and its children. | |
virtual bool | ClearEventRecursive ()=0 |
Recursively clear event information from this frame and its children. | |
bool | IsCircularTree (std::vector< RFKey > &keys) const |
Check this RestFrame 's tree for circular connections. | |
Protected Member Functions inherited from RestFrames::RFBase | |
bool | SetBody (bool body) const |
bool | SetMind (bool mind) const |
bool | SetSpirit (bool spirit) const |
virtual bool | IsSoundMind () const |
virtual bool | IsSoundSpirit () const |
void | UnSoundBody (const std::string &function) const |
void | UnSoundMind (const std::string &function) const |
void | UnSoundSpirit (const std::string &function) const |
void | AddDependent (RFBase *dep) |
pointer to RFBase object owned by this one | |
Protected Attributes | |
FrameType | m_Type |
AnaType | m_Ana |
Protected Attributes inherited from RestFrames::RFBase | |
RFLog | m_Log |
RFBase * | m_This |
Friends | |
class | ReconstructionFrame |
class | GeneratorFrame |
Additional Inherited Members | |
Static Protected Attributes inherited from RestFrames::RFBase | |
static const TVector3 | m_Empty3Vector |
static const TLorentzVector | m_Empty4Vector |
RestFrames::RestFrame::RestFrame | ( | const std::string & | sname, |
const std::string & | stitle | ||
) |
Standard constructor.
sname | Class instance name used for log statements |
stitle | Class instance title used in figures |
Definition at line 48 of file RestFrame.cc.
|
virtual |
Adds a child RestFrame to this frame.
frame | RestFrame to be added as child |
Method for adding a RestFrame as a child of this frame.
Reimplemented in RestFrames::GeneratorFrame, and RestFrames::ReconstructionFrame.
Definition at line 198 of file RestFrame.cc.
void RestFrames::RestFrame::AddChildFrames | ( | const RestFrameList & | frames | ) |
Adds a list of children to this frame.
frames | RestFrames to be added as children |
Method for adding a list of RestFrames as childs of this frame. Calls AddChildFrame() for each RestFrame on the list
Definition at line 225 of file RestFrame.cc.
|
static |
|
static |
|
static |
Retrieve axis which defines transverse plane.
Definition at line 76 of file RestFrame.cc.
TVector3 RestFrames::RestFrame::GetBoostInParentFrame | ( | ) | const |
Returns the boost of this frame in it's parent's frame.
Definition at line 889 of file RestFrame.cc.
|
virtual |
Returns the charge of this frame.
Definition at line 427 of file RestFrame.cc.
|
virtual |
Get the RestFrame of the i th child.
Reimplemented in RestFrames::ReconstructionFrame, and RestFrames::GeneratorFrame.
Definition at line 235 of file RestFrame.cc.
RestFrameList const & RestFrames::RestFrame::GetChildFrames | ( | ) | const |
Returns a list of this frame's child RestFrames.
Definition at line 252 of file RestFrame.cc.
double RestFrames::RestFrame::GetCosDecayAngle | ( | const RestFrame & | frame = RestFrame::Empty() | ) | const |
Returns the cosine of this frame's decay angle.
frame | (optional) Frame defining child axis |
Definition at line 912 of file RestFrame.cc.
TVector3 RestFrames::RestFrame::GetDecayPlaneNormalVector | ( | const RestFrame & | frame = RestFrame::Empty() | ) | const |
Returns the vector normal to the decay plane of this frame.
frame | (optional) Frame defining child axis |
\[ \hat{n}_{\perp}\ = \frac{ \vec{p}_{C}^{~P} \times \vec{p}_{F}^{~P} } {\left|\vec{p}_{C}^{~P} \times \vec{p}_{F}^{~P}\right|}~, \]
where \(\vec{p}_{F}^{~P}\) is the momentum of this frame evaluated in it's parent's rest frame and \(\vec{p}_{C}^{~P}\) is the momentum of it's first child frame evaluated in the same frame. If this frame is of type LabFrame then \(\hat{n}_{\perp}\) is defined alternatively as:\[ \hat{n}_{\perp}\ = \frac{ \vec{p}_{C}^{~F} \times \vec{n}_{\parallel} } {\left|\vec{p}_{C}^{~F} \times \vec{n}_{\parallel}\right|}~, \]
where \(\vec{p}_{C}^{~F}\) is the momentum of the child frame evaluated in this frame and \(\vec{n}_{\parallel}\) is the vector defining the transverse plane as returned by RestFrame::GetAxis().Definition at line 933 of file RestFrame.cc.
double RestFrames::RestFrame::GetDeltaPhiBoostVisible | ( | const TVector3 & | axis = RestFrame::GetAxis() , |
const RestFrame & | frame = RestFrame::Empty() |
||
) | const |
Returns difference of azimuthal angles between the decay angles of the visible particles decaying from this frame.
axis | Axis to define coordinates (default is z-axis) |
frame | (optional) frame defining child axis |
Method boosts to visible particles frame and calculates difference of azimuthal angles between the decay angles of the visible particles decaying from this frame
Definition at line 1026 of file RestFrame.cc.
double RestFrames::RestFrame::GetDeltaPhiDecayAngle | ( | const TVector3 & | axis = RestFrame::GetAxis() , |
const RestFrame & | frame = RestFrame::Empty() |
||
) | const |
Returns difference of azimuthal decay angles in an axis.
axis | Axis to define coordinates (default is z-axis) |
frame | (optional) Frame defining child axis |
Definition at line 997 of file RestFrame.cc.
double RestFrames::RestFrame::GetDeltaPhiDecayPlanes | ( | const RestFrame & | frame | ) | const |
Returns difference of azimuthal angle between decay planes.
frame | Frame corresponding to other decay plane |
Definition at line 977 of file RestFrame.cc.
double RestFrames::RestFrame::GetDeltaPhiDecayVisible | ( | const TVector3 & | axis = RestFrame::GetAxis() , |
const RestFrame & | frame = RestFrame::Empty() |
||
) | const |
Returns difference of azimuthal angles between the decay angles of the visible particles decaying from this frame.
axis | Axis to define coordinates (default is z-axis) |
frame | (optional) frame defining child axis |
Definition at line 1052 of file RestFrame.cc.
double RestFrames::RestFrame::GetDeltaPhiVisible | ( | const TVector3 & | axis = RestFrame::GetAxis() , |
const RestFrame & | frame = RestFrame::Empty() |
||
) | const |
Returns difference of azimuthal angles between the angles of its visible particles kids.
axis | Axis to define coordinates (default is z-axis) |
frame | (optional) frame defining child axis |
Definition at line 1084 of file RestFrame.cc.
double RestFrames::RestFrame::GetEnergy | ( | const RestFrame & | frame | ) | const |
Returns energy of this frame in specified reference frame.
frame | Rest frame in which to evaluate energy |
Definition at line 544 of file RestFrame.cc.
TLorentzVector RestFrames::RestFrame::GetFourVector | ( | const RestFrame & | frame = RestFrame::Empty() | ) | const |
Returns this frame's four-vector in a specified frame.
frame | Rest frame in which to evaluate four-vector |
Definition at line 451 of file RestFrame.cc.
TLorentzVector RestFrames::RestFrame::GetFourVector | ( | const TLorentzVector & | P, |
const RestFrame & | def_frame = RestFrame::Empty() |
||
) | const |
Returns four vector boosted to different frame.
P | Four-vector to boost |
<br> |
Definition at line 595 of file RestFrame.cc.
|
virtual |
Returns the frame at depth.
depth | Depth of desired frame |
frame | Descendent frame that defines path |
Definition at line 806 of file RestFrame.cc.
int RestFrames::RestFrame::GetFrameDepth | ( | const RestFrame & | frame | ) | const |
Returns the depth of frame
frame | Frame whose depth is returned |
Definition at line 790 of file RestFrame.cc.
double RestFrames::RestFrame::GetGammaInParentFrame | ( | ) | const |
Returns the gamma of this frame in its parent's frame.
Definition at line 901 of file RestFrame.cc.
TLorentzVector RestFrames::RestFrame::GetInvisibleFourVector | ( | const RestFrame & | frame = RestFrame::Empty() | ) | const |
Returns four-vector of invisible descendants.
frame | Rest frame in which to evaluate four-vector |
Definition at line 517 of file RestFrame.cc.
|
virtual |
Returns the LabFrame that this frame inherits from.
Definition at line 256 of file RestFrame.cc.
|
virtual |
Returns a list of frames inheriting from this one.
Definition at line 338 of file RestFrame.cc.
|
virtual |
Returns the mass of this frame.
Reimplemented in RestFrames::GeneratorFrame.
Definition at line 443 of file RestFrame.cc.
double RestFrames::RestFrame::GetMomentum | ( | const RestFrame & | frame | ) | const |
Returns magnitude of momentum.
frame | Rest frame in which to evaluate momentum |
Definition at line 552 of file RestFrame.cc.
|
virtual |
Returns the parent of this frame.
Reimplemented in RestFrames::ReconstructionFrame, and RestFrames::GeneratorFrame.
Definition at line 245 of file RestFrame.cc.
|
virtual |
Returns the production frame of this frame.
Definition at line 877 of file RestFrame.cc.
|
virtual |
Returns the sibling frame of this frame.
Definition at line 269 of file RestFrame.cc.
double RestFrames::RestFrame::GetSumInvisibleMomentum | ( | ) | const |
Returns scalar sum of invisible child momenta.
Definition at line 863 of file RestFrame.cc.
double RestFrames::RestFrame::GetSumVisibleMomentum | ( | ) | const |
Returns scalar sum of visible child momenta.
Definition at line 849 of file RestFrame.cc.
TLorentzVector RestFrames::RestFrame::GetTransverseFourVector | ( | const RestFrame & | frame = RestFrame::Empty() , |
const TVector3 & | axis = RestFrame::GetAxis() , |
||
const RestFrame & | axis_frame = RestFrame::Empty() |
||
) | const |
Returns this frame's four-vector in a specified frame.
frame | Rest frame in which to evaluate four-vector |
axis | Longitudinal axis |
axis_frame | Rest frame in which axis is defined |
Definition at line 633 of file RestFrame.cc.
double RestFrames::RestFrame::GetTransverseMomentum | ( | const RestFrame & | frame = RestFrame::Empty() , |
const TVector3 & | axis = RestFrame::GetAxis() , |
||
const RestFrame & | axis_frame = RestFrame::Empty() |
||
) | const |
Returns magnitude of transverse momentum.
frame | Rest frame in which to evaluate momentum |
axis | Longitudinal axis |
axis_frame | Rest frame in which axis is defined |
Definition at line 560 of file RestFrame.cc.
double RestFrames::RestFrame::GetTransverseMomentum | ( | const TLorentzVector & | P, |
const TVector3 & | axis = RestFrame::GetAxis() , |
||
const RestFrame & | axis_frame = RestFrame::Empty() |
||
) | const |
Returns magnitude of transverse momentum of this frame.
P | Four-vector of frame |
axis | Longitudinal axis |
axis_frame | Rest frame in which axis is defined |
Definition at line 733 of file RestFrame.cc.
TLorentzVector RestFrames::RestFrame::GetVisibleFourVector | ( | const RestFrame & | frame = RestFrame::Empty() | ) | const |
Returns four-vector of visible descendants.
frame | Rest frame in which to evaluate four-vector |
Definition at line 490 of file RestFrame.cc.
double RestFrames::RestFrame::GetVisibleShape | ( | ) | const |
Returns visible shape of frame.
\[ \mathrm{visible~shape} \equiv \frac{ \sqrt{ \sum_{i}^{N}\sum_{j < i}^{N} 2\left(\left|\vec{p}_{i}\right|\left|\vec{p}_{j}\right|+ \vec{p}_{i}\cdot\vec{p}_{j}\right) } } { \sum_{i}^{N}\left|\vec{p}_{i}\right| }~. \]
Definition at line 826 of file RestFrame.cc.
RestFrameList RestFrames::RestFrame::operator+ | ( | const RestFrameList & | frames | ) |
Combines RestFrames into RestFrameList.
frames | List of additional RestFrames to add in list |
Definition at line 87 of file RestFrame.cc.
RestFrameList RestFrames::RestFrame::operator+ | ( | RestFrame & | frame | ) |
Combines RestFrames into RestFrameList.
frame | Additional RestFrame to add in list |
Definition at line 80 of file RestFrame.cc.
void RestFrames::RestFrame::RemoveChildFrame | ( | RestFrame & | frame | ) |
Removes a child of this frame.
frame | Child frame to be removed |
Method for removing a child RestFrame from the list of children of this frame
Definition at line 164 of file RestFrame.cc.
void RestFrames::RestFrame::RemoveChildFrames | ( | ) |
Removes all the children of this frame.
Method for removing all the children of this frame. No child left behind.
Definition at line 173 of file RestFrame.cc.
|
static |
Set axis perpendicular to transverse plane.
axis | Input axis |
Method sets axis to define the "transverse plane", perpendicular to the axis.
Definition at line 72 of file RestFrame.cc.
|
virtual |
Sets the parent frame for this frame.
frame | Frame to be set as parent |
Method for connecting a child frame to its parent frame Empty default sets parent frame to none
Reimplemented in RestFrames::ReconstructionFrame, and RestFrames::GeneratorFrame.
Definition at line 181 of file RestFrame.cc.
#include "RestFrames/RestFrame.hh"