Geant4  10.02.p03
RE04Trajectory Class Reference

#include <RE04Trajectory.hh>

Inheritance diagram for RE04Trajectory:
Collaboration diagram for RE04Trajectory:

Public Member Functions

 RE04Trajectory ()
 
 RE04Trajectory (const G4Track *aTrack)
 
 RE04Trajectory (RE04Trajectory &)
 
virtual ~RE04Trajectory ()
 
void * operator new (size_t)
 
void operator delete (void *)
 
int operator== (const RE04Trajectory &right) const
 
virtual G4int GetTrackID () const
 
virtual G4int GetParentID () const
 
virtual G4String GetParticleName () const
 
virtual G4double GetCharge () const
 
virtual G4int GetPDGEncoding () const
 
virtual G4double GetInitialKineticEnergy () const
 
virtual G4ThreeVector GetInitialMomentum () const
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory () const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
G4ParticleDefinitionGetParticleDefinition ()
 
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
 

Private Attributes

TrajectoryPointContainerfPositionRecord
 
G4int fTrackID
 
G4int fParentID
 
G4int fPDGEncoding
 
G4double fPDGCharge
 
G4String fParticleName
 
G4double fInitialKineticEnergy
 
G4ThreeVector fInitialMomentum
 

Detailed Description

User trajectory class

Definition at line 90 of file RE04Trajectory.hh.

Constructor & Destructor Documentation

◆ RE04Trajectory() [1/3]

RE04Trajectory::RE04Trajectory ( )

Definition at line 48 of file RE04Trajectory.cc.

49 : G4VTrajectory(),
51  fPDGEncoding( 0 ), fPDGCharge(0.0), fParticleName(""),
53 {;}
TrajectoryPointContainer * fPositionRecord
CLHEP::Hep3Vector G4ThreeVector
G4double fInitialKineticEnergy
G4String fParticleName
G4ThreeVector fInitialMomentum

◆ RE04Trajectory() [2/3]

RE04Trajectory::RE04Trajectory ( const G4Track *  aTrack)

Definition at line 56 of file RE04Trajectory.cc.

57 {
58  G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
59  fParticleName = fpParticleDefinition->GetParticleName();
60  fPDGCharge = fpParticleDefinition->GetPDGCharge();
61  fPDGEncoding = fpParticleDefinition->GetPDGEncoding();
62  fTrackID = aTrack->GetTrackID();
63  fParentID = aTrack->GetParentID();
64  fInitialKineticEnergy = aTrack->GetKineticEnergy();
65  fInitialMomentum = aTrack->GetMomentum();
67  // Following is for the first trajectory point
68  fPositionRecord->push_back(new RE04TrajectoryPoint(
69  aTrack->GetPosition(),aTrack->GetMaterial()));
70 }
TrajectoryPointContainer * fPositionRecord
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
G4double fInitialKineticEnergy
const G4String & GetParticleName() const
G4String fParticleName
G4double GetPDGCharge() const
G4ThreeVector fInitialMomentum
Here is the call graph for this function:

◆ RE04Trajectory() [3/3]

RE04Trajectory::RE04Trajectory ( RE04Trajectory right)

Definition at line 73 of file RE04Trajectory.cc.

73  :G4VTrajectory()
74 {
76  fPDGCharge = right.fPDGCharge;
77  fPDGEncoding = right.fPDGEncoding;
78  fTrackID = right.fTrackID;
79  fParentID = right.fParentID;
83 
84  for(size_t i=0;i<right.fPositionRecord->size();i++)
85  {
86  RE04TrajectoryPoint* rightPoint
87  = (RE04TrajectoryPoint*)((*(right.fPositionRecord))[i]);
88  fPositionRecord->push_back(new RE04TrajectoryPoint(*rightPoint));
89  }
90 }
TrajectoryPointContainer * fPositionRecord
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
G4double fInitialKineticEnergy
G4String fParticleName
G4ThreeVector fInitialMomentum

◆ ~RE04Trajectory()

RE04Trajectory::~RE04Trajectory ( )
virtual

Definition at line 93 of file RE04Trajectory.cc.

94 {
95  if (fPositionRecord) {
96  // fPositionRecord->clearAndDestroy();
97  size_t i;
98  for(i=0;i<fPositionRecord->size();i++){
99  delete (*fPositionRecord)[i];
100  }
101  fPositionRecord->clear();
102  delete fPositionRecord;
103  }
104 }
TrajectoryPointContainer * fPositionRecord

Member Function Documentation

◆ AppendStep()

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

Implements G4VTrajectory.

Definition at line 205 of file RE04Trajectory.cc.

206 {
207  fPositionRecord->push_back( new RE04TrajectoryPoint(
208  aStep->GetPostStepPoint()->GetPosition(),
209  aStep->GetPreStepPoint()->GetMaterial() ));
210 }
TrajectoryPointContainer * fPositionRecord
Here is the caller graph for this function:

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 167 of file RE04Trajectory.cc.

168 {
169  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
170 
171  values->push_back
173 
174  values->push_back
176 
177  values->push_back(G4AttValue("PN",fParticleName,""));
178 
179  values->push_back
181 
182  values->push_back
184 
185  values->push_back
186  (G4AttValue("IKE",G4BestUnit(fInitialKineticEnergy,"Energy"),""));
187 
188  values->push_back
189  (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
190 
191  values->push_back
192  (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
193 
194  values->push_back
196 
197 #ifdef G4ATTDEBUG
198  G4cout << G4AttCheck(values,GetAttDefs());
199 #endif
200 
201  return values;
202 }
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:371
G4double fInitialKineticEnergy
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual int GetPointEntries() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4GLOB_DLL std::ostream G4cout
double mag() const
G4String fParticleName
G4ThreeVector fInitialMomentum
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DrawTrajectory()

void RE04Trajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 115 of file RE04Trajectory.cc.

116 {
117  // Invoke the default implementation in G4VTrajectory...
119  // ... or override with your own code here.
120 }
virtual void DrawTrajectory() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 123 of file RE04Trajectory.cc.

124 {
125  G4bool isNew;
126  std::map<G4String,G4AttDef>* store
127  = G4AttDefStore::GetInstance("RE04Trajectory",isNew);
128  if (isNew) {
129 
130  G4String id("ID");
131  (*store)[id] = G4AttDef(id,"Track ID","Physics","","G4int");
132 
133  G4String pid("PID");
134  (*store)[pid] = G4AttDef(pid,"Parent ID","Physics","","G4int");
135 
136  G4String pn("PN");
137  (*store)[pn] = G4AttDef(pn,"Particle Name","Physics","","G4String");
138 
139  G4String ch("Ch");
140  (*store)[ch] = G4AttDef(ch,"Charge","Physics","e+","G4double");
141 
142  G4String pdg("PDG");
143  (*store)[pdg] = G4AttDef(pdg,"PDG Encoding","Physics","","G4int");
144 
145  G4String ike("IKE");
146  (*store)[ike] =
147  G4AttDef(ike, "Initial kinetic energy",
148  "Physics","G4BestUnit","G4double");
149 
150  G4String iMom("IMom");
151  (*store)[iMom] = G4AttDef(iMom, "Initial momentum",
152  "Physics","G4BestUnit","G4ThreeVector");
153 
154  G4String iMag("IMag");
155  (*store)[iMag] =
156  G4AttDef(iMag, "Initial momentum magnitude",
157  "Physics","G4BestUnit","G4double");
158 
159  G4String ntp("NTP");
160  (*store)[ntp] = G4AttDef(ntp,"No. of points","Physics","","G4int");
161 
162  }
163  return store;
164 }
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 RE04Trajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 119 of file RE04Trajectory.hh.

120  { return fPDGCharge; }

◆ GetInitialKineticEnergy()

virtual G4double RE04Trajectory::GetInitialKineticEnergy ( ) const
inlinevirtual

Definition at line 123 of file RE04Trajectory.hh.

124  { return fInitialKineticEnergy; }
G4double fInitialKineticEnergy

◆ GetInitialMomentum()

virtual G4ThreeVector RE04Trajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 125 of file RE04Trajectory.hh.

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

◆ GetParentID()

virtual G4int RE04Trajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 115 of file RE04Trajectory.hh.

116  { return fParentID; }

◆ GetParticleDefinition()

G4ParticleDefinition * RE04Trajectory::GetParticleDefinition ( )

Definition at line 213 of file RE04Trajectory.cc.

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

◆ GetParticleName()

virtual G4String RE04Trajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 117 of file RE04Trajectory.hh.

118  { return fParticleName; }
G4String fParticleName

◆ GetPDGEncoding()

virtual G4int RE04Trajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 121 of file RE04Trajectory.hh.

122  { return fPDGEncoding; }

◆ GetPoint()

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

Implements G4VTrajectory.

Definition at line 134 of file RE04Trajectory.hh.

135  { return (*fPositionRecord)[i]; }
TrajectoryPointContainer * fPositionRecord
Here is the call graph for this function:

◆ GetPointEntries()

virtual int RE04Trajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 133 of file RE04Trajectory.hh.

133 { return fPositionRecord->size(); }
TrajectoryPointContainer * fPositionRecord
Here is the caller graph for this function:

◆ GetTrackID()

virtual G4int RE04Trajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 113 of file RE04Trajectory.hh.

114  { return fTrackID; }

◆ MergeTrajectory()

void RE04Trajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 219 of file RE04Trajectory.cc.

220 {
221  if(!secondTrajectory) return;
222 
223  RE04Trajectory* seco = (RE04Trajectory*)secondTrajectory;
224  G4int ent = seco->GetPointEntries();
225  for(G4int i=1;i<ent;i++) // initial point of the second trajectory
226  // should not be merged
227  {
228  fPositionRecord->push_back((*(seco->fPositionRecord))[i]);
229  // fPositionRecord->push_back(seco->fPositionRecord->removeAt(1));
230  }
231  delete (*seco->fPositionRecord)[0];
232  seco->fPositionRecord->clear();
233 }
TrajectoryPointContainer * fPositionRecord
virtual int GetPointEntries() const
int G4int
Definition: G4Types.hh:78
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator delete()

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

Definition at line 166 of file RE04Trajectory.hh.

167 {
168  faTrajAllocator->FreeSingle((RE04Trajectory*)aTrajectory);
169 }
G4ThreadLocal G4Allocator< RE04Trajectory > * faTrajAllocator

◆ operator new()

void * RE04Trajectory::operator new ( size_t  )
inline

Definition at line 160 of file RE04Trajectory.hh.

161 {
163  return (void*)faTrajAllocator->MallocSingle();
164 }
G4ThreadLocal G4Allocator< RE04Trajectory > * faTrajAllocator

◆ operator==()

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

Definition at line 109 of file RE04Trajectory.hh.

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

◆ ShowTrajectory()

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

Reimplemented from G4VTrajectory.

Definition at line 107 of file RE04Trajectory.cc.

108 {
109  // Invoke the default implementation in G4VTrajectory...
111  // ... or override with your own code here.
112 }
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

◆ fInitialKineticEnergy

G4double RE04Trajectory::fInitialKineticEnergy
private

Definition at line 153 of file RE04Trajectory.hh.

◆ fInitialMomentum

G4ThreeVector RE04Trajectory::fInitialMomentum
private

Definition at line 154 of file RE04Trajectory.hh.

◆ fParentID

G4int RE04Trajectory::fParentID
private

Definition at line 149 of file RE04Trajectory.hh.

◆ fParticleName

G4String RE04Trajectory::fParticleName
private

Definition at line 152 of file RE04Trajectory.hh.

◆ fPDGCharge

G4double RE04Trajectory::fPDGCharge
private

Definition at line 151 of file RE04Trajectory.hh.

◆ fPDGEncoding

G4int RE04Trajectory::fPDGEncoding
private

Definition at line 150 of file RE04Trajectory.hh.

◆ fPositionRecord

TrajectoryPointContainer* RE04Trajectory::fPositionRecord
private

Definition at line 147 of file RE04Trajectory.hh.

◆ fTrackID

G4int RE04Trajectory::fTrackID
private

Definition at line 148 of file RE04Trajectory.hh.


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