LOGO

RestFrames  v1.0.1
RestFrames HEP Event Analysis Software Library
RestFrame.hh
Go to the documentation of this file.
1 // RestFrames: particle physics event analysis library
3 // --------------------------------------------------------------------
4 // Copyright (c) 2014-2018, 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 RestFrame_HH
31 #define RestFrame_HH
32 
33 #include "RestFrames/RFBase.hh"
34 #include "RestFrames/RFCharge.hh"
35 
36 namespace RestFrames {
37 
39  enum FrameType { kVanillaFrame, kVisibleFrame, kInvisibleFrame,
40  kDecayFrame, kLabFrame};
41 
43  enum AnaType { kRecoFrame, kGenFrame };
44 
45  class RestFrame : public RFBase {
46  public:
47 
52  RestFrame(const std::string& sname, const std::string& stitle);
53 
55  RestFrame();
56 
57  virtual ~RestFrame();
58 
65  static void SetAxis(const TVector3& axis);
66 
70  static TVector3 const& GetAxis();
71 
73  virtual void Clear();
74 
83 
85  FrameType GetType() const;
86 
88  bool IsVisibleFrame() const;
89 
91  bool IsInvisibleFrame() const;
92 
94  bool IsDecayFrame() const;
95 
97  bool IsLabFrame() const;
98 
100  bool IsRecoFrame() const;
101 
103  bool IsGenFrame() const;
104 
106  virtual std::string PrintString(LogType type) const;
107 
109 
117 
127  virtual void AddChildFrame(RestFrame& frame);
128 
136  void AddChildFrames(const RestFrameList& frames);
137 
144  virtual void SetParentFrame(RestFrame& frame =
145  RestFrame::Empty());
146 
155  void RemoveChildFrame(RestFrame& frame);
156 
161  void RemoveChildFrames();
162 
170 
172  int GetNChildren() const;
173 
175  int GetNDescendants() const;
176 
183  virtual RestFrame const& GetParentFrame() const;
184 
188  virtual RestFrame& GetChildFrame(int i = 0) const;
189 
193  RestFrameList const& GetChildFrames() const;
194 
201  virtual RestFrame const& GetLabFrame() const;
202 
206  virtual RestFrame const& GetProductionFrame() const;
207 
214  virtual RestFrame const& GetSiblingFrame() const;
215 
224  int GetFrameDepth(const RestFrame& frame) const;
225 
237  virtual RestFrame const& GetFrameAtDepth(int depth, const RestFrame& frame) const;
238 
245  virtual RestFrameList GetListFrames(FrameType type = kLabFrame) const;
246 
248  virtual RestFrameList GetListVisibleFrames() const;
249 
251  virtual RestFrameList GetListInvisibleFrames() const;
252 
254 
265 
272 
278  RestFrameList operator + (const RestFrameList& frames);
279 
283  virtual RFCharge GetCharge() const;
284 
288  virtual double GetMass() const;
289 
300  TLorentzVector GetFourVector(const RestFrame& frame =
301  RestFrame::Empty()) const;
302 
314  TLorentzVector GetTransverseFourVector(const RestFrame& frame =
316  const TVector3& axis =
318  const RestFrame& axis_frame =
319  RestFrame::Empty()) const;
320 
330  TLorentzVector GetVisibleFourVector(const RestFrame& frame =
331  RestFrame::Empty()) const;
332 
342  TLorentzVector GetInvisibleFourVector(const RestFrame& frame =
343  RestFrame::Empty()) const;
344 
354  double GetEnergy(const RestFrame& frame) const;
355 
365  double GetMomentum(const RestFrame& frame) const;
366 
371  TLorentzVector GetFourVector(const TLorentzVector& P,
372  const RestFrame& def_frame =
373  RestFrame::Empty()) const;
374 
377  TLorentzVector GetTransverseFourVector(const TLorentzVector& P,
378  const TVector3& axis =
380  const RestFrame& axis_frame =
381  RestFrame::Empty()) const;
382 
394  double GetTransverseMomentum(const RestFrame& frame =
396  const TVector3& axis =
398  const RestFrame& axis_frame =
399  RestFrame::Empty()) const;
400 
412  double GetTransverseMomentum(const TLorentzVector& P,
413  const TVector3& axis =
415  const RestFrame& axis_frame =
416  RestFrame::Empty()) const;
417 
432  double GetVisibleShape() const;
433 
440  double GetSumVisibleMomentum() const;
441 
448  double GetSumInvisibleMomentum() const;
449 
456  TVector3 GetBoostInParentFrame() const;
457 
461  double GetGammaInParentFrame() const;
462 
485  TVector3 GetDecayPlaneNormalVector(const RestFrame& frame =
486  RestFrame::Empty()) const;
487 
495  double GetDeltaPhiDecayPlanes(const RestFrame& frame) const;
496 
502  double GetCosDecayAngle(const RestFrame& frame =
503  RestFrame::Empty()) const;
504 
511  double GetDeltaPhiDecayAngle(const TVector3& axis =
513  const RestFrame& frame =
514  RestFrame::Empty()) const;
515 
525  double GetDeltaPhiBoostVisible(const TVector3& axis =
527  const RestFrame& frame =
528  RestFrame::Empty()) const;
529 
536  double GetDeltaPhiDecayVisible(const TVector3& axis =
538  const RestFrame& frame =
539  RestFrame::Empty()) const;
540 
547  double GetDeltaPhiVisible(const TVector3& axis =
549  const RestFrame& frame =
550  RestFrame::Empty()) const;
552 
556  static RestFrame& Empty();
557 
561  static ConstRestFrameList const& EmptyList();
562 
563  protected:
564  FrameType m_Type;
565  AnaType m_Ana;
566 
567  virtual bool IsSoundBody() const;
568 
569  TVector3 const& GetChildBoostVector(RestFrame& frame) const;
570  TVector3 const& GetParentBoostVector() const;
571 
573  virtual bool InitializeTreeRecursive();
574 
576  virtual bool InitializeAnalysisRecursive() = 0;
577 
579  virtual bool AnalyzeEventRecursive() = 0;
580 
582  virtual bool ClearEventRecursive() = 0;
583 
585  bool IsCircularTree(std::vector<RFKey>& keys) const;
586 
587  private:
589  static int m_class_key;
590 
592  static TVector3 m_Axis;
593 
594  // 4-vector of this state in the frame it's initialized
595  TLorentzVector m_P;
596 
597  // the reference frame where this four-vector is initialized
598  const RestFrame* m_ProdFramePtr;
599 
600  // list of child frames and boosts
601  RestFrameList m_ChildFrames;
602  mutable std::map<const RestFrame*, TVector3> m_ChildBoosts;
603 
604  // parent frame and boost
605  RestFrame* m_ParentFramePtr;
606  TVector3 m_ParentBoost;
607 
608  void SetFourVector(const TLorentzVector& V, const RestFrame& frame);
609  void SetChildBoostVector(RestFrame& frame, const TVector3& boost);
610  void SetParentBoostVector(const TVector3& boost);
611 
612  // Recursively get lists of frames
613  void FillListFramesRecursive(RestFrameList& frames,
614  FrameType type = kLabFrame) const;
615 
616  bool FindPathToFrame(const RestFrame& dest_frame,
617  const RestFrame& prev_frame,
618  std::vector<const TVector3*>& boosts) const;
619 
620  static const ConstRestFrameList m_EmptyList;
621 
622  friend class ReconstructionFrame;
623  friend class GeneratorFrame;
624 
625  };
626 
627 }
628 
629 #endif
RestFrames::RestFrame::GetMomentum
double GetMomentum(const RestFrame &frame) const
Returns magnitude of momentum.
Definition: RestFrame.cc:552
RestFrames::RestFrame::GetLabFrame
virtual RestFrame const & GetLabFrame() const
Returns the LabFrame that this frame inherits from.
Definition: RestFrame.cc:256
RestFrames::RestFrame::GetVisibleShape
double GetVisibleShape() const
Returns visible shape of frame.
Definition: RestFrame.cc:826
RestFrames::RestFrame::SetParentFrame
virtual void SetParentFrame(RestFrame &frame=RestFrame::Empty())
Sets the parent frame for this frame.
Definition: RestFrame.cc:181
RestFrames::RestFrame::GetVisibleFourVector
TLorentzVector GetVisibleFourVector(const RestFrame &frame=RestFrame::Empty()) const
Returns four-vector of visible descendants.
Definition: RestFrame.cc:490
RestFrames::RestFrame::GetChildFrame
virtual RestFrame & GetChildFrame(int i=0) const
Get the RestFrame of the i th child.
Definition: RestFrame.cc:235
RestFrames::RestFrame::RemoveChildFrames
void RemoveChildFrames()
Removes all the children of this frame.
Definition: RestFrame.cc:173
RestFrames::RestFrame::GetListInvisibleFrames
virtual RestFrameList GetListInvisibleFrames() const
Returns a list of InvisibleFrame s inheriting from this.
Definition: RestFrame.cc:348
RestFrames::RestFrame::GetNDescendants
int GetNDescendants() const
Returns the number of descendents of this frame.
Definition: RestFrame.cc:288
RestFrames::RestFrame::InitializeTreeRecursive
virtual bool InitializeTreeRecursive()
Recursively initialize this frame's tree.
Definition: RestFrame.cc:150
RestFrames::RestFrame::GetType
FrameType GetType() const
Returns RestFrame (FrameType) type.
Definition: RestFrame.cc:93
RestFrames::AnaType
AnaType
Type of RestFrame, with respect to its analysis capabilities.
Definition: RestFrame.hh:43
RFCharge.hh
RestFrames::RFCharge
Definition: RFCharge.hh:40
RestFrames::RestFrame::GetDecayPlaneNormalVector
TVector3 GetDecayPlaneNormalVector(const RestFrame &frame=RestFrame::Empty()) const
Returns the vector normal to the decay plane of this frame.
Definition: RestFrame.cc:933
RestFrames::RestFrame::GetTransverseMomentum
double GetTransverseMomentum(const RestFrame &frame=RestFrame::Empty(), const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &axis_frame=RestFrame::Empty()) const
Returns magnitude of transverse momentum.
Definition: RestFrame.cc:560
RestFrames::RestFrame::GetDeltaPhiVisible
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.
Definition: RestFrame.cc:1084
RestFrames::RestFrame::GetDeltaPhiDecayVisible
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 fro...
Definition: RestFrame.cc:1052
RestFrames::RestFrame::IsRecoFrame
bool IsRecoFrame() const
Is this an ReconstructionFrame ? (yes/no)
Definition: RestFrame.cc:113
RestFrames::RestFrame::IsDecayFrame
bool IsDecayFrame() const
Is this a DecayFrame ? (yes/no)
Definition: RestFrame.cc:105
RestFrames::RestFrame
abstract base class for all Frame objects
Definition: RestFrame.hh:45
RestFrames::RestFrame::AddChildFrame
virtual void AddChildFrame(RestFrame &frame)
Adds a child RestFrame to this frame.
Definition: RestFrame.cc:198
RestFrames::RestFrame::GetFourVector
TLorentzVector GetFourVector(const RestFrame &frame=RestFrame::Empty()) const
Returns this frame's four-vector in a specified frame.
Definition: RestFrame.cc:451
RestFrames::RestFrame::IsGenFrame
bool IsGenFrame() const
Is this a GeneratorFrame ? (yes/no)
Definition: RestFrame.cc:117
RestFrames::RestFrame::GetCosDecayAngle
double GetCosDecayAngle(const RestFrame &frame=RestFrame::Empty()) const
Returns the cosine of this frame's decay angle.
Definition: RestFrame.cc:912
RestFrames::RFList< RestFrame >
RestFrames::RestFrame::GetSumVisibleMomentum
double GetSumVisibleMomentum() const
Returns scalar sum of visible child momenta.
Definition: RestFrame.cc:849
RestFrames::RestFrame::SetAxis
static void SetAxis(const TVector3 &axis)
Set axis perpendicular to transverse plane.
Definition: RestFrame.cc:72
RestFrames::RFBase
Base class for all RestFrame package objects.
Definition: RFBase.hh:53
RestFrames::RestFrame::GetFrameAtDepth
virtual RestFrame const & GetFrameAtDepth(int depth, const RestFrame &frame) const
Returns the frame at depth.
Definition: RestFrame.cc:806
RestFrames::RestFrame::InitializeAnalysisRecursive
virtual bool InitializeAnalysisRecursive()=0
Recursively initialize this frame and its children for analysis.
RestFrames::RestFrame::GetSumInvisibleMomentum
double GetSumInvisibleMomentum() const
Returns scalar sum of invisible child momenta.
Definition: RestFrame.cc:863
RestFrames::RestFrame::GetProductionFrame
virtual RestFrame const & GetProductionFrame() const
Returns the production frame of this frame.
Definition: RestFrame.cc:877
RestFrames::GeneratorFrame
Definition: GeneratorFrame.hh:42
RestFrames::RestFrame::GetNChildren
int GetNChildren() const
Returns the number of child frames inheriting from this one.
Definition: RestFrame.cc:231
RestFrames::RestFrame::PrintString
virtual std::string PrintString(LogType type) const
String of information about RestFrame.
Definition: RestFrame.cc:121
RestFrames::RestFrame::AddChildFrames
void AddChildFrames(const RestFrameList &frames)
Adds a list of children to this frame.
Definition: RestFrame.cc:225
RestFrames::RestFrame::GetChildFrames
RestFrameList const & GetChildFrames() const
Returns a list of this frame's child RestFrames.
Definition: RestFrame.cc:252
RestFrames::RestFrame::IsCircularTree
bool IsCircularTree(std::vector< RFKey > &keys) const
Check this RestFrame 's tree for circular connections.
Definition: RestFrame.cc:360
RestFrames::LogType
LogType
Type of Log Message.
Definition: RFLog.hh:45
RFBase.hh
RestFrames::RestFrame::GetParentFrame
virtual RestFrame const & GetParentFrame() const
Returns the parent of this frame.
Definition: RestFrame.cc:245
RestFrames::RestFrame::IsInvisibleFrame
bool IsInvisibleFrame() const
Is this an InvisibleFrame ? (yes/no)
Definition: RestFrame.cc:101
RestFrames::RestFrame::GetTransverseFourVector
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.
Definition: RestFrame.cc:633
RestFrames::RestFrame::GetCharge
virtual RFCharge GetCharge() const
Returns the charge of this frame.
Definition: RestFrame.cc:427
RestFrames::RestFrame::GetSiblingFrame
virtual RestFrame const & GetSiblingFrame() const
Returns the sibling frame of this frame.
Definition: RestFrame.cc:269
RestFrames::RestFrame::GetInvisibleFourVector
TLorentzVector GetInvisibleFourVector(const RestFrame &frame=RestFrame::Empty()) const
Returns four-vector of invisible descendants.
Definition: RestFrame.cc:517
RestFrames::RestFrame::Clear
virtual void Clear()
Clears RestFrame of all connections to other objects.
Definition: RestFrame.cc:57
RestFrames::RestFrame::GetListVisibleFrames
virtual RestFrameList GetListVisibleFrames() const
Returns a list of VisibleFrame s inheriting from this.
Definition: RestFrame.cc:344
RestFrames::RestFrame::GetMass
virtual double GetMass() const
Returns the mass of this frame.
Definition: RestFrame.cc:443
RestFrames::RestFrame::GetListFrames
virtual RestFrameList GetListFrames(FrameType type=kLabFrame) const
Returns a list of frames inheriting from this one.
Definition: RestFrame.cc:338
RestFrames::RestFrame::GetGammaInParentFrame
double GetGammaInParentFrame() const
Returns the gamma of this frame in its parent's frame.
Definition: RestFrame.cc:901
RestFrames::FrameType
FrameType
Type of RestFrame, with respect to its decays.
Definition: RestFrame.hh:39
RestFrames::RestFrame::IsVisibleFrame
bool IsVisibleFrame() const
Is this a VisibleFrame ? (yes/no)
Definition: RestFrame.cc:97
RestFrames::RestFrame::GetAxis
static TVector3 const & GetAxis()
Retrieve axis which defines transverse plane.
Definition: RestFrame.cc:76
RestFrames::RestFrame::RestFrame
RestFrame()
Empty constructor.
Definition: RestFrame.cc:39
RestFrames::RestFrame::EmptyList
static ConstRestFrameList const & EmptyList()
Returns empty RestFrameList.
Definition: RestFrame.cc:66
RestFrames::RestFrame::ClearEventRecursive
virtual bool ClearEventRecursive()=0
Recursively clear event information from this frame and its children.
RestFrames::RestFrame::GetEnergy
double GetEnergy(const RestFrame &frame) const
Returns energy of this frame in specified reference frame.
Definition: RestFrame.cc:544
RestFrames::RestFrame::GetFrameDepth
int GetFrameDepth(const RestFrame &frame) const
Returns the depth of frame
Definition: RestFrame.cc:790
RestFrames::RestFrame::IsLabFrame
bool IsLabFrame() const
Is this a LabFrame ? (yes/no)
Definition: RestFrame.cc:109
RestFrames::RestFrame::Empty
static RestFrame & Empty()
Returns empty RestFrame.
Definition: RestFrame.cc:62
RestFrames::RestFrame::RemoveChildFrame
void RemoveChildFrame(RestFrame &frame)
Removes a child of this frame.
Definition: RestFrame.cc:164
RestFrames::RestFrame::GetBoostInParentFrame
TVector3 GetBoostInParentFrame() const
Returns the boost of this frame in it's parent's frame.
Definition: RestFrame.cc:889
RestFrames::RestFrame::AnalyzeEventRecursive
virtual bool AnalyzeEventRecursive()=0
Recursively analyze event in this frame and its children.
RestFrames::RestFrame::GetDeltaPhiDecayPlanes
double GetDeltaPhiDecayPlanes(const RestFrame &frame) const
Returns difference of azimuthal angle between decay planes.
Definition: RestFrame.cc:977
RestFrames::ReconstructionFrame
Definition: ReconstructionFrame.hh:41
RestFrames::RestFrame::GetDeltaPhiDecayAngle
double GetDeltaPhiDecayAngle(const TVector3 &axis=RestFrame::GetAxis(), const RestFrame &frame=RestFrame::Empty()) const
Returns difference of azimuthal decay angles in an axis.
Definition: RestFrame.cc:997
RestFrames::RestFrame::GetDeltaPhiBoostVisible
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 fro...
Definition: RestFrame.cc:1026
RestFrames::RestFrame::operator+
RestFrameList operator+(RestFrame &frame)
Combines RestFrames into RestFrameList.
Definition: RestFrame.cc:80