Geant4  10.02.p03
WLSTrajectory Class Reference

#include <WLSTrajectory.hh>

Inheritance diagram for WLSTrajectory:
Collaboration diagram for WLSTrajectory:

Public Member Functions

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

Private Attributes

WLSTrajectoryPointContainerfpPointsContainer
 
G4int fTrackID
 
G4int fParentID
 
G4double fPDGCharge
 
G4int fPDGEncoding
 
G4String fParticleName
 
G4ThreeVector fInitialMomentum
 
G4ParticleDefinitionfParticleDefinition
 

Detailed Description

Definition at line 53 of file WLSTrajectory.hh.

Constructor & Destructor Documentation

◆ WLSTrajectory() [1/3]

WLSTrajectory::WLSTrajectory ( )

Definition at line 59 of file WLSTrajectory.cc.

63 {
64  fParticleDefinition = NULL;
65 }
CLHEP::Hep3Vector G4ThreeVector
G4String fParticleName
WLSTrajectoryPointContainer * fpPointsContainer
G4ParticleDefinition * fParticleDefinition
G4double fPDGCharge
G4ThreeVector fInitialMomentum

◆ WLSTrajectory() [2/3]

WLSTrajectory::WLSTrajectory ( const G4Track *  aTrack)

Definition at line 69 of file WLSTrajectory.cc.

70 {
71  fParticleDefinition = aTrack->GetDefinition();
75  fTrackID = aTrack->GetTrackID();
76  fParentID = aTrack->GetParentID();
77  fInitialMomentum = aTrack->GetMomentum();
79  // Following is for the first trajectory point
80  fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack));
81 }
G4String fParticleName
const G4String & GetParticleName() const
WLSTrajectoryPointContainer * fpPointsContainer
G4ParticleDefinition * fParticleDefinition
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
G4double fPDGCharge
G4ThreeVector fInitialMomentum
G4double GetPDGCharge() const
Here is the call graph for this function:

◆ WLSTrajectory() [3/3]

WLSTrajectory::WLSTrajectory ( WLSTrajectory right)

Definition at line 85 of file WLSTrajectory.cc.

85  : G4VTrajectory()
86 {
89  fPDGCharge = right.fPDGCharge;
90  fPDGEncoding = right.fPDGEncoding;
91  fTrackID = right.fTrackID;
92  fParentID = right.fParentID;
95 
96  for(size_t i=0;i<right.fpPointsContainer->size();++i) {
97  WLSTrajectoryPoint* rightPoint
98  = (WLSTrajectoryPoint*)((*(right.fpPointsContainer))[i]);
99  fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint));
100  }
101 }
G4String fParticleName
WLSTrajectoryPointContainer * fpPointsContainer
G4ParticleDefinition * fParticleDefinition
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
G4double fPDGCharge
G4ThreeVector fInitialMomentum

◆ ~WLSTrajectory()

WLSTrajectory::~WLSTrajectory ( )
virtual

Definition at line 105 of file WLSTrajectory.cc.

106 {
107  for(size_t i=0;i<fpPointsContainer->size();++i){
108  delete (*fpPointsContainer)[i];
109  }
110  fpPointsContainer->clear();
111 
112  delete fpPointsContainer;
113 }
WLSTrajectoryPointContainer * fpPointsContainer

Member Function Documentation

◆ AppendStep()

void WLSTrajectory::AppendStep ( const G4Step *  aStep)
virtual

Implements G4VTrajectory.

Definition at line 126 of file WLSTrajectory.cc.

127 {
128  fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep));
129 }
WLSTrajectoryPointContainer * fpPointsContainer
Here is the caller graph for this function:

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 198 of file WLSTrajectory.cc.

199 {
200  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
201 
202  values->push_back
204 
205  values->push_back
207 
208  values->push_back(G4AttValue("PN",fParticleName,""));
209 
210  values->push_back
212 
213  values->push_back
215 
216  values->push_back
217  (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
218 
219  values->push_back
220  (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
221 
222  values->push_back
224 
225 #ifdef G4ATTDEBUG
226  G4cout << G4AttCheck(values,GetAttDefs());
227 #endif
228  return values;
229 }
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:371
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4String fParticleName
G4GLOB_DLL std::ostream G4cout
double mag() const
virtual int GetPointEntries() const
G4double fPDGCharge
G4ThreeVector fInitialMomentum
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 156 of file WLSTrajectory.cc.

157 {
158  G4bool isNew;
159  std::map<G4String,G4AttDef>* store
160  = G4AttDefStore::GetInstance("Trajectory",isNew);
161 
162  if (isNew) {
163 
164  G4String ID("ID");
165  (*store)[ID] = G4AttDef(ID,"Track ID","Bookkeeping","","G4int");
166 
167  G4String PID("PID");
168  (*store)[PID] = G4AttDef(PID,"Parent ID","Bookkeeping","","G4int");
169 
170  G4String PN("PN");
171  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
172 
173  G4String Ch("Ch");
174  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
175 
176  G4String PDG("PDG");
177  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
178 
179  G4String IMom("IMom");
180  (*store)[IMom] = G4AttDef(IMom,
181  "Momentum of track at start of trajectory",
182  "Physics","G4BestUnit","G4ThreeVector");
183 
184  G4String IMag("IMag");
185  (*store)[IMag] = G4AttDef(IMag,
186  "Magnitude of momentum of track at start of trajectory",
187  "Physics","G4BestUnit","G4double");
188 
189  G4String NTP("NTP");
190  (*store)[NTP] = G4AttDef(NTP,"No. of points","Bookkeeping","","G4int");
191 
192  }
193  return store;
194 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCharge()

virtual G4double WLSTrajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 79 of file WLSTrajectory.hh.

79 { return fPDGCharge; }
G4double fPDGCharge

◆ GetInitialMomentum()

virtual G4ThreeVector WLSTrajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 81 of file WLSTrajectory.hh.

82  {return fInitialMomentum;}
G4ThreeVector fInitialMomentum
Here is the call graph for this function:

◆ GetParentID()

virtual G4int WLSTrajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 77 of file WLSTrajectory.hh.

77 { return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * WLSTrajectory::GetParticleDefinition ( )

Definition at line 133 of file WLSTrajectory.cc.

134 {
135  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
136 }
G4String fParticleName
static G4ParticleTable * GetParticleTable()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetParticleName()

virtual G4String WLSTrajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 78 of file WLSTrajectory.hh.

78 { return fParticleName; }
G4String fParticleName

◆ GetPDGEncoding()

virtual G4int WLSTrajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 80 of file WLSTrajectory.hh.

80 { return fPDGEncoding; }

◆ GetPoint()

virtual G4VTrajectoryPoint* WLSTrajectory::GetPoint ( G4int  i) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 94 of file WLSTrajectory.hh.

95  { return (*fpPointsContainer)[i]; }
WLSTrajectoryPointContainer * fpPointsContainer
Here is the call graph for this function:

◆ GetPointEntries()

virtual int WLSTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 92 of file WLSTrajectory.hh.

93  { return fpPointsContainer->size(); }
WLSTrajectoryPointContainer * fpPointsContainer
Here is the caller graph for this function:

◆ GetTrackID()

virtual G4int WLSTrajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 76 of file WLSTrajectory.hh.

76 { return fTrackID; }

◆ MergeTrajectory()

void WLSTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 140 of file WLSTrajectory.cc.

141 {
142  if(!secondTrajectory) return;
143 
144  WLSTrajectory* second = (WLSTrajectory*)secondTrajectory;
145  G4int ent = second->GetPointEntries();
146  // initial point of the second trajectory should not be merged
147  for(G4int i=1; i<ent; ++i) {
148  fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
149  }
150  delete (*second->fpPointsContainer)[0];
151  second->fpPointsContainer->clear();
152 }
int G4int
Definition: G4Types.hh:78
static const double second
Definition: G4SIunits.hh:156
virtual int GetPointEntries() const
WLSTrajectoryPointContainer * fpPointsContainer
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator delete()

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

Definition at line 127 of file WLSTrajectory.hh.

127  {
128  WLSTrajectoryAllocator->FreeSingle((WLSTrajectory*)aTrajectory);
129 }
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator

◆ operator new()

void * WLSTrajectory::operator new ( size_t  )
inline

Definition at line 121 of file WLSTrajectory.hh.

121  {
124  return (void*) WLSTrajectoryAllocator->MallocSingle();
125 }
G4ThreadLocal G4Allocator< WLSTrajectory > * WLSTrajectoryAllocator

◆ operator==()

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

Definition at line 71 of file WLSTrajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 117 of file WLSTrajectory.cc.

118 {
119  // Invoke the default implementation in G4VTrajectory...
121  // ... or override with your own code here.
122 }
virtual void ShowTrajectory(std::ostream &os=G4cout) const
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fInitialMomentum

G4ThreeVector WLSTrajectory::fInitialMomentum
private

Definition at line 113 of file WLSTrajectory.hh.

◆ fParentID

G4int WLSTrajectory::fParentID
private

Definition at line 109 of file WLSTrajectory.hh.

◆ fParticleDefinition

G4ParticleDefinition* WLSTrajectory::fParticleDefinition
private

Definition at line 115 of file WLSTrajectory.hh.

◆ fParticleName

G4String WLSTrajectory::fParticleName
private

Definition at line 112 of file WLSTrajectory.hh.

◆ fPDGCharge

G4double WLSTrajectory::fPDGCharge
private

Definition at line 110 of file WLSTrajectory.hh.

◆ fPDGEncoding

G4int WLSTrajectory::fPDGEncoding
private

Definition at line 111 of file WLSTrajectory.hh.

◆ fpPointsContainer

WLSTrajectoryPointContainer* WLSTrajectory::fpPointsContainer
private

Definition at line 106 of file WLSTrajectory.hh.

◆ fTrackID

G4int WLSTrajectory::fTrackID
private

Definition at line 108 of file WLSTrajectory.hh.


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