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

#include <Doxymodules_field.h>

Inheritance diagram for F04Trajectory:
G4VTrajectory

Public Member Functions

 F04Trajectory ()
 
 F04Trajectory (const G4Track *aTrack)
 
 F04Trajectory (F04Trajectory &)
 
 ~F04Trajectory () override
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const F04Trajectory &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

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

Detailed Description

Definition at line 123 of file Doxymodules_field.h.

Constructor & Destructor Documentation

◆ F04Trajectory() [1/3]

F04Trajectory::F04Trajectory ( )

Definition at line 51 of file F04Trajectory.cc.

52{}

◆ F04Trajectory() [2/3]

F04Trajectory::F04Trajectory ( const G4Track aTrack)

Definition at line 56 of file F04Trajectory.cc.

57{
58 G4ParticleDefinition * particleDefinition = aTrack->GetDefinition();
59 fParticleName = particleDefinition->GetParticleName();
60 fPDGCharge = particleDefinition->GetPDGCharge();
61 fPDGEncoding = particleDefinition->GetPDGEncoding();
62 fTrackID = aTrack->GetTrackID();
63 fParentID = aTrack->GetParentID();
64 fInitialMomentum = aTrack->GetMomentum();
66 // Following is for the first trajectory point
67 fpPointsContainer->push_back(new F04TrajectoryPoint(aTrack));
68}
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
G4String fParticleName
G4ThreeVector fInitialMomentum
G4double fPDGCharge
TrajectoryPointContainer * fpPointsContainer

◆ F04Trajectory() [3/3]

F04Trajectory::F04Trajectory ( F04Trajectory right)

Definition at line 72 of file F04Trajectory.cc.

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

◆ ~F04Trajectory()

F04Trajectory::~F04Trajectory ( )
override

Definition at line 91 of file F04Trajectory.cc.

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

Member Function Documentation

◆ operator new()

void * F04Trajectory::operator new ( size_t  )
inline

Definition at line 114 of file F04Trajectory.hh.

114 {
117 return (void*) F04TrajectoryAllocator->MallocSingle();
118}
G4ThreadLocal G4Allocator< F04Trajectory > * F04TrajectoryAllocator

◆ operator delete()

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

Definition at line 120 of file F04Trajectory.hh.

120 {
121 F04TrajectoryAllocator->FreeSingle((F04Trajectory*)aTrajectory);
122}

◆ operator==()

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

Definition at line 67 of file F04Trajectory.hh.

68 { return (this==&right); }

◆ GetTrackID()

G4int F04Trajectory::GetTrackID ( ) const
inlineoverride

Definition at line 72 of file F04Trajectory.hh.

72{ return fTrackID; }

◆ GetParentID()

G4int F04Trajectory::GetParentID ( ) const
inlineoverride

Definition at line 73 of file F04Trajectory.hh.

73{ return fParentID; }

◆ GetParticleName()

G4String F04Trajectory::GetParticleName ( ) const
inlineoverride

Definition at line 74 of file F04Trajectory.hh.

74{ return fParticleName; }

◆ GetCharge()

G4double F04Trajectory::GetCharge ( ) const
inlineoverride

Definition at line 75 of file F04Trajectory.hh.

75{ return fPDGCharge; }

◆ GetPDGEncoding()

G4int F04Trajectory::GetPDGEncoding ( ) const
inlineoverride

Definition at line 76 of file F04Trajectory.hh.

76{ return fPDGEncoding; }

◆ GetInitialMomentum()

G4ThreeVector F04Trajectory::GetInitialMomentum ( ) const
inlineoverride

Definition at line 77 of file F04Trajectory.hh.

77{return fInitialMomentum;}

◆ ShowTrajectory()

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

Definition at line 103 of file F04Trajectory.cc.

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

◆ AppendStep()

void F04Trajectory::AppendStep ( const G4Step aStep)
override

Definition at line 112 of file F04Trajectory.cc.

113{
114 fpPointsContainer->push_back(new F04TrajectoryPoint(aStep));
115}

◆ MergeTrajectory()

void F04Trajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
override

Definition at line 126 of file F04Trajectory.cc.

127{
128 if(!secondTrajectory) return;
129
130 auto second = (F04Trajectory*)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 fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
135 }
136 delete (*second->fpPointsContainer)[0];
137 second->fpPointsContainer->clear();
138}

◆ GetParticleDefinition()

G4ParticleDefinition * F04Trajectory::GetParticleDefinition ( )

Definition at line 119 of file F04Trajectory.cc.

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

◆ GetPointEntries()

int F04Trajectory::GetPointEntries ( ) const
inlineoverride

Definition at line 87 of file F04Trajectory.hh.

88 { return fpPointsContainer->size(); }

◆ GetPoint()

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

Definition at line 89 of file F04Trajectory.hh.

90 { return (*fpPointsContainer)[i]; }

◆ GetAttDefs()

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

Definition at line 142 of file F04Trajectory.cc.

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

◆ CreateAttValues()

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

Definition at line 184 of file F04Trajectory.cc.

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

Member Data Documentation

◆ fpPointsContainer

TrajectoryPointContainer* F04Trajectory::fpPointsContainer = nullptr
private

Definition at line 101 of file F04Trajectory.hh.

◆ fTrackID

G4int F04Trajectory::fTrackID = 0
private

Definition at line 103 of file F04Trajectory.hh.

◆ fParentID

G4int F04Trajectory::fParentID = 0
private

Definition at line 104 of file F04Trajectory.hh.

◆ fPDGCharge

G4double F04Trajectory::fPDGCharge = 0.
private

Definition at line 105 of file F04Trajectory.hh.

◆ fPDGEncoding

G4int F04Trajectory::fPDGEncoding = 0
private

Definition at line 106 of file F04Trajectory.hh.

◆ fParticleName

G4String F04Trajectory::fParticleName
private

Definition at line 107 of file F04Trajectory.hh.

◆ fInitialMomentum

G4ThreeVector F04Trajectory::fInitialMomentum
private

Definition at line 108 of file F04Trajectory.hh.


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

Applications | User Support | Publications | Collaboration