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

#include <Doxymodules_optical.h>

Inheritance diagram for WLSTrajectory:
G4VTrajectory

Public Member Functions

 WLSTrajectory ()=default
 
 WLSTrajectory (const G4Track *)
 
 WLSTrajectory (WLSTrajectory &)
 
 ~WLSTrajectory () override
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const WLSTrajectory &right) const
 
G4int GetTrackID () const override
 
G4int GetParentID () const override
 
G4String GetParticleName () const override
 
G4double GetCharge () const override
 
G4int GetPDGEncoding () const override
 
G4ThreeVector GetInitialMomentum () const override
 
void ShowTrajectory (std::ostream &os=G4cout) const override
 
void AppendStep (const G4Step *aStep) override
 
void MergeTrajectory (G4VTrajectory *secondTrajectory) override
 
G4ParticleDefinitionGetParticleDefinition ()
 
int GetPointEntries () const override
 
G4VTrajectoryPointGetPoint (G4int i) const override
 
const std::map< G4String, G4AttDef > * GetAttDefs () const override
 
std::vector< G4AttValue > * CreateAttValues () const override
 

Private Attributes

WLSTrajectoryPointContainerfpPointsContainer = nullptr
 
G4int fTrackID = 0
 
G4int fParentID = 0
 
G4double fPDGCharge = 0.
 
G4int fPDGEncoding = 0
 
G4String fParticleName
 
G4ThreeVector fInitialMomentum
 
G4ParticleDefinitionfParticleDefinition = nullptr
 

Detailed Description

Definition at line 106 of file Doxymodules_optical.h.

Constructor & Destructor Documentation

◆ WLSTrajectory() [1/3]

WLSTrajectory::WLSTrajectory ( )
default

◆ WLSTrajectory() [2/3]

WLSTrajectory::WLSTrajectory ( const G4Track aTrack)

Definition at line 53 of file WLSTrajectory.cc.

54{
55 fParticleDefinition = aTrack->GetDefinition();
56 fParticleName = fParticleDefinition->GetParticleName();
57 fPDGCharge = fParticleDefinition->GetPDGCharge();
58 fPDGEncoding = fParticleDefinition->GetPDGEncoding();
59 fTrackID = aTrack->GetTrackID();
60 fParentID = aTrack->GetParentID();
61 fInitialMomentum = aTrack->GetMomentum();
63 // Following is for the first trajectory point
64 fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack));
65}
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
G4ThreeVector fInitialMomentum
WLSTrajectoryPointContainer * fpPointsContainer
G4ParticleDefinition * fParticleDefinition
G4String fParticleName
G4double fPDGCharge

◆ WLSTrajectory() [3/3]

WLSTrajectory::WLSTrajectory ( WLSTrajectory right)

Definition at line 69 of file WLSTrajectory.cc.

71{
72 fParticleDefinition = right.fParticleDefinition;
73 fParticleName = right.fParticleName;
74 fPDGCharge = right.fPDGCharge;
75 fPDGEncoding = right.fPDGEncoding;
76 fTrackID = right.fTrackID;
77 fParentID = right.fParentID;
78 fInitialMomentum = right.fInitialMomentum;
80
81 for(size_t i = 0; i < right.fpPointsContainer->size(); ++i)
82 {
83 auto rightPoint = (WLSTrajectoryPoint*) ((*(right.fpPointsContainer))[i]);
84 fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint));
85 }
86}

◆ ~WLSTrajectory()

WLSTrajectory::~WLSTrajectory ( )
override

Definition at line 90 of file WLSTrajectory.cc.

91{
92 for(size_t i = 0; i < fpPointsContainer->size(); ++i)
93 {
94 delete(*fpPointsContainer)[i];
95 }
96 fpPointsContainer->clear();
97 delete fpPointsContainer;
98}

Member Function Documentation

◆ operator new()

void * WLSTrajectory::operator new ( size_t  )
inline

Definition at line 104 of file WLSTrajectory.hh.

105{
108 return (void*) WLSTrajectoryAllocator->MallocSingle();
109}
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator

◆ operator delete()

void WLSTrajectory::operator delete ( void *  aTrajectory)
inline

Definition at line 111 of file WLSTrajectory.hh.

112{
113 WLSTrajectoryAllocator->FreeSingle((WLSTrajectory*) aTrajectory);
114}

◆ operator==()

int WLSTrajectory::operator== ( const WLSTrajectory right) const
inline

Definition at line 59 of file WLSTrajectory.hh.

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

◆ GetTrackID()

G4int WLSTrajectory::GetTrackID ( ) const
inlineoverride

Definition at line 64 of file WLSTrajectory.hh.

64{ return fTrackID; }

◆ GetParentID()

G4int WLSTrajectory::GetParentID ( ) const
inlineoverride

Definition at line 65 of file WLSTrajectory.hh.

65{ return fParentID; }

◆ GetParticleName()

G4String WLSTrajectory::GetParticleName ( ) const
inlineoverride

Definition at line 66 of file WLSTrajectory.hh.

66{ return fParticleName; }

◆ GetCharge()

G4double WLSTrajectory::GetCharge ( ) const
inlineoverride

Definition at line 67 of file WLSTrajectory.hh.

67{ return fPDGCharge; }

◆ GetPDGEncoding()

G4int WLSTrajectory::GetPDGEncoding ( ) const
inlineoverride

Definition at line 68 of file WLSTrajectory.hh.

68{ return fPDGEncoding; }

◆ GetInitialMomentum()

G4ThreeVector WLSTrajectory::GetInitialMomentum ( ) const
inlineoverride

Definition at line 69 of file WLSTrajectory.hh.

70 {
71 return fInitialMomentum;
72 }

◆ ShowTrajectory()

void WLSTrajectory::ShowTrajectory ( std::ostream &  os = G4cout) const
override

Definition at line 102 of file WLSTrajectory.cc.

103{
104 // Invoke the default implementation in G4VTrajectory...
105 G4VTrajectory::ShowTrajectory(os);
106 // ... or override with your own code here.
107}

◆ AppendStep()

void WLSTrajectory::AppendStep ( const G4Step aStep)
override

Definition at line 111 of file WLSTrajectory.cc.

112{
113 fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep));
114}

◆ MergeTrajectory()

void WLSTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
override

Definition at line 125 of file WLSTrajectory.cc.

126{
127 if(!secondTrajectory)
128 return;
129
130 auto second = (WLSTrajectory*) secondTrajectory;
131 G4int ent = second->GetPointEntries();
132 // initial point of the second trajectory should not be merged
133 for(G4int i = 1; i < ent; ++i)
134 {
135 fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
136 }
137 delete(*second->fpPointsContainer)[0];
138 second->fpPointsContainer->clear();
139}

◆ GetParticleDefinition()

G4ParticleDefinition * WLSTrajectory::GetParticleDefinition ( )

Definition at line 118 of file WLSTrajectory.cc.

119{
120 return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
121}

◆ GetPointEntries()

int WLSTrajectory::GetPointEntries ( ) const
inlineoverride

Definition at line 80 of file WLSTrajectory.hh.

80{ return fpPointsContainer->size(); }

◆ GetPoint()

G4VTrajectoryPoint * WLSTrajectory::GetPoint ( G4int  i) const
inlineoverride

Definition at line 81 of file WLSTrajectory.hh.

82 {
83 return (*fpPointsContainer)[i];
84 }

◆ GetAttDefs()

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

Definition at line 143 of file WLSTrajectory.cc.

144{
145 G4bool isNew;
146 std::map<G4String, G4AttDef>* store =
147 G4AttDefStore::GetInstance("Trajectory", isNew);
148
149 if(isNew)
150 {
151 G4String ID("ID");
152 (*store)[ID] = G4AttDef(ID, "Track ID", "Bookkeeping", "", "G4int");
153
154 G4String PID("PID");
155 (*store)[PID] = G4AttDef(PID, "Parent ID", "Bookkeeping", "", "G4int");
156
157 G4String PN("PN");
158 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
159
160 G4String Ch("Ch");
161 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
162
163 G4String PDG("PDG");
164 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
165
166 G4String IMom("IMom");
167 (*store)[IMom] = G4AttDef(IMom, "Momentum of track at start of trajectory",
168 "Physics", "G4BestUnit", "G4ThreeVector");
169
170 G4String IMag("IMag");
171 (*store)[IMag] =
172 G4AttDef(IMag, "Magnitude of momentum of track at start of trajectory",
173 "Physics", "G4BestUnit", "G4double");
174
175 G4String NTP("NTP");
176 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Bookkeeping", "", "G4int");
177 }
178 return store;
179}

◆ CreateAttValues()

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

Definition at line 183 of file WLSTrajectory.cc.

184{
185 auto values = new std::vector<G4AttValue>;
186
187 values->push_back(
188 G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
189
190 values->push_back(
191 G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
192
193 values->push_back(G4AttValue("PN", fParticleName, ""));
194
195 values->push_back(
196 G4AttValue("Ch", G4UIcommand::ConvertToString(fPDGCharge), ""));
197
198 values->push_back(
199 G4AttValue("PDG", G4UIcommand::ConvertToString(fPDGEncoding), ""));
200
201 values->push_back(
202 G4AttValue("IMom", G4BestUnit(fInitialMomentum, "Energy"), ""));
203
204 values->push_back(
205 G4AttValue("IMag", G4BestUnit(fInitialMomentum.mag(), "Energy"), ""));
206
207 values->push_back(
208 G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
209
210 return values;
211}
int GetPointEntries() const override

Member Data Documentation

◆ fpPointsContainer

WLSTrajectoryPointContainer* WLSTrajectory::fpPointsContainer = nullptr
private

Definition at line 90 of file WLSTrajectory.hh.

◆ fTrackID

G4int WLSTrajectory::fTrackID = 0
private

Definition at line 92 of file WLSTrajectory.hh.

◆ fParentID

G4int WLSTrajectory::fParentID = 0
private

Definition at line 93 of file WLSTrajectory.hh.

◆ fPDGCharge

G4double WLSTrajectory::fPDGCharge = 0.
private

Definition at line 94 of file WLSTrajectory.hh.

◆ fPDGEncoding

G4int WLSTrajectory::fPDGEncoding = 0
private

Definition at line 95 of file WLSTrajectory.hh.

◆ fParticleName

G4String WLSTrajectory::fParticleName
private

Definition at line 96 of file WLSTrajectory.hh.

◆ fInitialMomentum

G4ThreeVector WLSTrajectory::fInitialMomentum
private

Definition at line 97 of file WLSTrajectory.hh.

◆ fParticleDefinition

G4ParticleDefinition* WLSTrajectory::fParticleDefinition = nullptr
private

Definition at line 99 of file WLSTrajectory.hh.


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

Applications | User Support | Publications | Collaboration