Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
WLSTrajectoryPoint Class Reference

#include <Doxymodules_optical.h>

Inheritance diagram for WLSTrajectoryPoint:
G4TrajectoryPoint G4VTrajectoryPoint

Public Member Functions

 WLSTrajectoryPoint ()=default
 
 WLSTrajectoryPoint (const G4Track *)
 
 WLSTrajectoryPoint (const G4Step *)
 
 WLSTrajectoryPoint (const WLSTrajectoryPoint &right)
 
 ~WLSTrajectoryPoint () override=default
 
void * operator new (size_t)
 
void operator delete (void *aTrajectoryPoint)
 
G4bool operator== (const WLSTrajectoryPoint &right) const
 
G4double GetTime () const
 
const G4ThreeVector GetMomentum () const
 
G4StepStatus GetStepStatus () const
 
G4String GetVolumeName () const
 
const std::map< G4String, G4AttDef > * GetAttDefs () const override
 
std::vector< G4AttValue > * CreateAttValues () const override
 

Private Attributes

G4double fTime = 0.
 
G4ThreeVector fMomentum
 
G4StepStatus fStepStatus = fUndefined
 
G4String fVolumeName
 

Detailed Description

Definition at line 107 of file Doxymodules_optical.h.

Constructor & Destructor Documentation

◆ WLSTrajectoryPoint() [1/4]

WLSTrajectoryPoint::WLSTrajectoryPoint ( )
default

◆ WLSTrajectoryPoint() [2/4]

WLSTrajectoryPoint::WLSTrajectoryPoint ( const G4Track aTrack)

Definition at line 66 of file WLSTrajectoryPoint.cc.

67 : G4TrajectoryPoint(aTrack->GetPosition())
68{
69 fTime = aTrack->GetGlobalTime();
70 fMomentum = aTrack->GetMomentum();
71 fStepStatus = fUndefined;
72 fVolumeName = aTrack->GetVolume()->GetName();
73}

◆ WLSTrajectoryPoint() [3/4]

WLSTrajectoryPoint::WLSTrajectoryPoint ( const G4Step aStep)

Definition at line 47 of file WLSTrajectoryPoint.cc.

48 : G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition())
49{
50 auto postStepPoint = aStep->GetPostStepPoint();
51 fTime = postStepPoint->GetGlobalTime();
52 fMomentum = postStepPoint->GetMomentum();
53 fStepStatus = postStepPoint->GetStepStatus();
54 if(postStepPoint->GetPhysicalVolume())
55 {
56 fVolumeName = postStepPoint->GetPhysicalVolume()->GetName();
57 }
58 else
59 {
60 fVolumeName = G4String(" ");
61 }
62}

◆ WLSTrajectoryPoint() [4/4]

WLSTrajectoryPoint::WLSTrajectoryPoint ( const WLSTrajectoryPoint right)

Definition at line 77 of file WLSTrajectoryPoint.cc.

79{
80 fTime = right.fTime;
81 fMomentum = right.fMomentum;
82 fStepStatus = right.fStepStatus;
83 fVolumeName = right.fVolumeName;
84}

◆ ~WLSTrajectoryPoint()

WLSTrajectoryPoint::~WLSTrajectoryPoint ( )
overridedefault

Member Function Documentation

◆ operator new()

void * WLSTrajectoryPoint::operator new ( size_t  )
inline

Definition at line 80 of file WLSTrajectoryPoint.hh.

81{
84 return (void*) WLSTrajPointAllocator->MallocSingle();
85}
G4ThreadLocal G4Allocator< WLSTrajectoryPoint > * WLSTrajPointAllocator

◆ operator delete()

void WLSTrajectoryPoint::operator delete ( void *  aTrajectoryPoint)
inline

Definition at line 87 of file WLSTrajectoryPoint.hh.

88{
89 WLSTrajPointAllocator->FreeSingle((WLSTrajectoryPoint*) aTrajectoryPoint);
90}

◆ operator==()

G4bool WLSTrajectoryPoint::operator== ( const WLSTrajectoryPoint right) const
inline

Definition at line 58 of file WLSTrajectoryPoint.hh.

59 {
60 return (this == &right);
61 };

◆ GetTime()

G4double WLSTrajectoryPoint::GetTime ( ) const
inline

Definition at line 63 of file WLSTrajectoryPoint.hh.

63{ return fTime; };

◆ GetMomentum()

const G4ThreeVector WLSTrajectoryPoint::GetMomentum ( ) const
inline

Definition at line 64 of file WLSTrajectoryPoint.hh.

64{ return fMomentum; };

◆ GetStepStatus()

G4StepStatus WLSTrajectoryPoint::GetStepStatus ( ) const
inline

Definition at line 65 of file WLSTrajectoryPoint.hh.

65{ return fStepStatus; };

◆ GetVolumeName()

G4String WLSTrajectoryPoint::GetVolumeName ( ) const
inline

Definition at line 66 of file WLSTrajectoryPoint.hh.

66{ return fVolumeName; };

◆ GetAttDefs()

const std::map< G4String, G4AttDef > * WLSTrajectoryPoint::GetAttDefs ( ) const
override

Definition at line 88 of file WLSTrajectoryPoint.cc.

89{
90 G4bool isNew;
91 std::map<G4String, G4AttDef>* store =
92 G4AttDefStore::GetInstance("TrajectoryPoint", isNew);
93 if(isNew)
94 {
95 G4String Pos("Pos");
96 (*store)[Pos] =
97 G4AttDef(Pos, "Position", "Physics", "G4BestUnit", "G4ThreeVector");
98
99 G4String Time("Time");
100 (*store)[Time] =
101 G4AttDef(Time, "Time", "Physics", "G4BestUnit", "G4double");
102
103 G4String Momentum("Momentum");
104 (*store)[Momentum] =
105 G4AttDef(Momentum, "Momentum", "Physics", "G4BestUnit", "G4ThreeVector");
106
107 G4String StepStatus("StepStatus");
108 (*store)[StepStatus] =
109 G4AttDef(StepStatus, "StepStatus", "Physics", "", "G4StepStatus");
110
111 G4String VolumeName("VolumeName");
112 (*store)[VolumeName] =
113 G4AttDef(VolumeName, "VolumeName", "Physics", "", "G4String");
114 }
115 return store;
116}

◆ CreateAttValues()

std::vector< G4AttValue > * WLSTrajectoryPoint::CreateAttValues ( ) const
override

Definition at line 120 of file WLSTrajectoryPoint.cc.

121{
122 auto values = new std::vector<G4AttValue>;
123
124 values->push_back(G4AttValue("Time", G4BestUnit(fTime, "Time"), ""));
125 values->push_back(
126 G4AttValue("Momentum", G4BestUnit(fMomentum, "Momentum"), ""));
127 values->push_back(G4AttValue("StepStatus", G4UIcommand::ConvertToString(fStepStatus), ""));
128 values->push_back(G4AttValue("VolumeName", fVolumeName, ""));
129
130 return values;
131}

Member Data Documentation

◆ fTime

G4double WLSTrajectoryPoint::fTime = 0.
private

Definition at line 72 of file WLSTrajectoryPoint.hh.

◆ fMomentum

G4ThreeVector WLSTrajectoryPoint::fMomentum
private

Definition at line 73 of file WLSTrajectoryPoint.hh.

◆ fStepStatus

G4StepStatus WLSTrajectoryPoint::fStepStatus = fUndefined
private

Definition at line 74 of file WLSTrajectoryPoint.hh.

◆ fVolumeName

G4String WLSTrajectoryPoint::fVolumeName
private

Definition at line 75 of file WLSTrajectoryPoint.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration