Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 ()
 
voidoperator 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
 

Detailed Description

Definition at line 53 of file WLSTrajectory.hh.

Constructor & Destructor Documentation

WLSTrajectory::WLSTrajectory ( )

Definition at line 59 of file WLSTrajectory.cc.

60  : fpPointsContainer(0), fTrackID(0), fParentID(0),
61  fPDGCharge(0.0), fPDGEncoding(0), fParticleName(""),
62  fInitialMomentum(G4ThreeVector())
63 {
64  fParticleDefinition = NULL;
65 }
CLHEP::Hep3Vector G4ThreeVector
WLSTrajectory::WLSTrajectory ( const G4Track aTrack)

Definition at line 69 of file WLSTrajectory.cc.

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

Here is the call graph for this function:

WLSTrajectory::WLSTrajectory ( WLSTrajectory right)

Definition at line 85 of file WLSTrajectory.cc.

85  : G4VTrajectory()
86 {
87  fParticleDefinition=right.fParticleDefinition;
88  fParticleName = right.fParticleName;
89  fPDGCharge = right.fPDGCharge;
90  fPDGEncoding = right.fPDGEncoding;
91  fTrackID = right.fTrackID;
92  fParentID = right.fParentID;
93  fInitialMomentum = right.fInitialMomentum;
94  fpPointsContainer = new WLSTrajectoryPointContainer();
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 }
std::vector< G4VTrajectoryPoint * > WLSTrajectoryPointContainer
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 }

Member Function Documentation

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 }
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
203  (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
204 
205  values->push_back
206  (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
207 
208  values->push_back(G4AttValue("PN",fParticleName,""));
209 
210  values->push_back
211  (G4AttValue("Ch",G4UIcommand::ConvertToString(fPDGCharge),""));
212 
213  values->push_back
214  (G4AttValue("PDG",G4UIcommand::ConvertToString(fPDGEncoding),""));
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:372
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual int GetPointEntries() const
G4GLOB_DLL std::ostream G4cout
double mag() const

Here is the call graph for this function:

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:

virtual G4double WLSTrajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 79 of file WLSTrajectory.hh.

79 { return fPDGCharge; }
virtual G4ThreeVector WLSTrajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 81 of file WLSTrajectory.hh.

82  {return fInitialMomentum;}
virtual G4int WLSTrajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 77 of file WLSTrajectory.hh.

77 { return fParentID; }
G4ParticleDefinition * WLSTrajectory::GetParticleDefinition ( )

Definition at line 133 of file WLSTrajectory.cc.

134 {
135  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
136 }
static G4ParticleTable * GetParticleTable()

Here is the call graph for this function:

virtual G4String WLSTrajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 78 of file WLSTrajectory.hh.

78 { return fParticleName; }
virtual G4int WLSTrajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 80 of file WLSTrajectory.hh.

80 { return fPDGEncoding; }
virtual G4VTrajectoryPoint* WLSTrajectory::GetPoint ( G4int  i) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 94 of file WLSTrajectory.hh.

95  { return (*fpPointsContainer)[i]; }
virtual int WLSTrajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 92 of file WLSTrajectory.hh.

93  { return fpPointsContainer->size(); }

Here is the caller graph for this function:

virtual G4int WLSTrajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 76 of file WLSTrajectory.hh.

76 { return fTrackID; }
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 }
static constexpr double second
Definition: G4SIunits.hh:157
int G4int
Definition: G4Types.hh:78
virtual int GetPointEntries() const

Here is the call graph for this function:

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
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
int WLSTrajectory::operator== ( const WLSTrajectory right) const
inline

Definition at line 71 of file WLSTrajectory.hh.

72  { return (this==&right); }
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:


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