Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RE01Trajectory Class Reference

#include <RE01Trajectory.hh>

Inheritance diagram for RE01Trajectory:
Collaboration diagram for RE01Trajectory:

Public Member Functions

 RE01Trajectory (const G4Track *aTrack)
 
virtual ~RE01Trajectory ()
 
virtual void ShowTrajectory (std::ostream &os=G4cout) const
 
virtual void DrawTrajectory () const
 
virtual const std::map
< G4String, G4AttDef > * 
GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
virtual void AppendStep (const G4Step *aStep)
 
virtual void MergeTrajectory (G4VTrajectory *secondTrajectory)
 
voidoperator new (size_t)
 
void operator delete (void *)
 
int operator== (const RE01Trajectory &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 int GetPointEntries () const
 
virtual G4VTrajectoryPointGetPoint (G4int i) const
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()
 
virtual ~G4VTrajectory ()
 
G4bool operator== (const G4VTrajectory &right) const
 

Detailed Description

Definition at line 54 of file RE01Trajectory.hh.

Constructor & Destructor Documentation

RE01Trajectory::RE01Trajectory ( const G4Track aTrack)

Definition at line 53 of file RE01Trajectory.cc.

54 :G4VTrajectory(),
55  fPositionRecord(0),fParticleDefinition(0)
56 {
57  fParticleDefinition = aTrack->GetDefinition();
58  fParticleName = fParticleDefinition->GetParticleName();
59  fPDGCharge = fParticleDefinition->GetPDGCharge();
60  fPDGEncoding = fParticleDefinition->GetPDGEncoding();
61  if(fParticleName=="unknown")
62  {
64  if(pp)
65  {
66  if(pp->GetCharge()<DBL_MAX) fPDGCharge = pp->GetCharge();
67  fPDGEncoding = pp->GetPDGcode();
68  if(pp->GetG4code()!=0)
69  {
70  fParticleName += " : ";
71  fParticleName += pp->GetG4code()->GetParticleName();
72  }
73  }
74  }
75  fTrackID = aTrack->GetTrackID();
76  RE01TrackInformation* trackInfo
78  fTrackStatus = trackInfo->GetTrackingStatus();
79  if(fTrackStatus == 1)
80  { fParentID = aTrack->GetParentID(); }
81  else if(fTrackStatus == 2)
82  { fParentID = trackInfo->GetSourceTrackID(); }
83  else
84  { fParentID = -1; }
85  fPositionRecord = new RE01TrajectoryPointContainer();
86  fPositionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
87  fMomentum = aTrack->GetMomentum();
88  fVertexPosition = aTrack->GetPosition();
89  fGlobalTime = aTrack->GetGlobalTime();
90 }
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetPosition() const
const G4String & GetParticleName() const
G4ParticleDefinition * GetG4code() const
G4VUserTrackInformation * GetUserInformation() const
G4int GetTrackingStatus() const
G4int GetTrackID() const
G4double GetGlobalTime() const
G4PrimaryParticle * GetPrimaryParticle() const
G4int GetPDGcode() const
G4ThreeVector GetMomentum() const
std::vector< G4VTrajectoryPoint * > RE01TrajectoryPointContainer
G4double GetCharge() const
G4double GetPDGCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4int GetSourceTrackID() const

Here is the call graph for this function:

RE01Trajectory::~RE01Trajectory ( )
virtual

Definition at line 93 of file RE01Trajectory.cc.

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

Member Function Documentation

void RE01Trajectory::AppendStep ( const G4Step aStep)
virtual

Implements G4VTrajectory.

Definition at line 254 of file RE01Trajectory.cc.

255 {
256  fPositionRecord->push_back(
258 }
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const

Here is the call graph for this function:

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

Reimplemented from G4VTrajectory.

Definition at line 217 of file RE01Trajectory.cc.

218 {
219  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
220 
221  values->push_back
222  (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
223 
224  values->push_back
225  (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
226 
227  values->push_back
228  (G4AttValue("Status",G4UIcommand::ConvertToString(fTrackStatus),""));
229 
230  values->push_back(G4AttValue("PN",fParticleName,""));
231 
232  values->push_back
233  (G4AttValue("Ch",G4UIcommand::ConvertToString(fPDGCharge),""));
234 
235  values->push_back
236  (G4AttValue("PDG",G4UIcommand::ConvertToString(fPDGEncoding),""));
237 
238  values->push_back
239  (G4AttValue("IMom",G4BestUnit(fMomentum,"Energy"),""));
240 
241  values->push_back
242  (G4AttValue("IMag",G4BestUnit(fMomentum.mag(),"Energy"),""));
243 
244  values->push_back
245  (G4AttValue("VtxPos",G4BestUnit(fVertexPosition,"Length"),""));
246 
247  values->push_back
249 
250  return values;
251 }
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:372
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual int GetPointEntries() const
double mag() const

Here is the call graph for this function:

void RE01Trajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 127 of file RE01Trajectory.cc.

128 {
129 
132 
133  G4Polyline pPolyline;
134  for (size_t i = 0; i < fPositionRecord->size() ; i++) {
135  G4TrajectoryPoint* aTrajectoryPoint =
136  (G4TrajectoryPoint*)((*fPositionRecord)[i]);
137  pos = aTrajectoryPoint->GetPosition();
138  pPolyline.push_back( pos );
139  }
140 
141  G4Colour colour(0.2,0.2,0.2);
142  if(fParticleDefinition==G4Gamma::GammaDefinition())
143  colour = G4Colour(0.,0.,1.);
144  else if(fParticleDefinition==G4Electron::ElectronDefinition()
145  ||fParticleDefinition==G4Positron::PositronDefinition())
146  colour = G4Colour(1.,1.,0.);
147  else if(fParticleDefinition==G4MuonMinus::MuonMinusDefinition()
148  ||fParticleDefinition==G4MuonPlus::MuonPlusDefinition())
149  colour = G4Colour(0.,1.,0.);
150  else if(fParticleDefinition->GetParticleType()=="meson")
151  {
152  if(fPDGCharge!=0.)
153  colour = G4Colour(1.,0.,0.);
154  else
155  colour = G4Colour(0.5,0.,0.);
156  }
157  else if(fParticleDefinition->GetParticleType()=="baryon")
158  {
159  if(fPDGCharge!=0.)
160  colour = G4Colour(0.,1.,1.);
161  else
162  colour = G4Colour(0.,0.5,0.5);
163  }
164 
165  G4VisAttributes attribs(colour);
166  pPolyline.SetVisAttributes(attribs);
167  if(pVVisManager) pVVisManager->Draw(pPolyline);
168 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:89
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
static G4VVisManager * GetConcreteInstance()
const G4ThreeVector GetPosition() const
const G4String & GetParticleType() const
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
static G4Positron * PositronDefinition()
Definition: G4Positron.cc:89
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
static const G4double pos
static G4Gamma * GammaDefinition()
Definition: G4Gamma.cc:81

Here is the call graph for this function:

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

Reimplemented from G4VTrajectory.

Definition at line 171 of file RE01Trajectory.cc.

172 {
173  G4bool isNew;
174  std::map<G4String,G4AttDef>* store
175  = G4AttDefStore::GetInstance("RE01Trajectory",isNew);
176  if (isNew) {
177 
178  G4String id("ID");
179  (*store)[id] = G4AttDef(id,"Track ID","Bookkeeping","","G4int");
180 
181  G4String pid("PID");
182  (*store)[pid] = G4AttDef(pid,"Parent ID","Bookkeeping","","G4int");
183 
184  G4String status("Status");
185  (*store)[status] = G4AttDef(status,"Track Status","Bookkeeping","","G4int");
186 
187  G4String pn("PN");
188  (*store)[pn] = G4AttDef(pn,"Particle Name","Bookkeeping","","G4String");
189 
190  G4String ch("Ch");
191  (*store)[ch] = G4AttDef(ch,"Charge","Physics","e+","G4double");
192 
193  G4String pdg("PDG");
194  (*store)[pdg] = G4AttDef(pdg,"PDG Encoding","Bookkeeping","","G4int");
195 
196  G4String imom("IMom");
197  (*store)[imom] = G4AttDef(imom, "Momentum of track at start of trajectory",
198  "Physics","G4BestUnit","G4ThreeVector");
199 
200  G4String imag("IMag");
201  (*store)[imag] =
202  G4AttDef(imag, "Magnitude of momentum of track at start of trajectory",
203  "Physics","G4BestUnit","G4double");
204 
205  G4String vtxPos("VtxPos");
206  (*store)[vtxPos] = G4AttDef(vtxPos, "Vertex position",
207  "Physics","G4BestUnit","G4ThreeVector");
208 
209  G4String ntp("NTP");
210  (*store)[ntp] = G4AttDef(ntp,"No. of points","Bookkeeping","","G4int");
211 
212  }
213  return store;
214 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Here is the call graph for this function:

virtual G4double RE01Trajectory::GetCharge ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 75 of file RE01Trajectory.hh.

75 { return fPDGCharge; }
virtual G4ThreeVector RE01Trajectory::GetInitialMomentum ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 77 of file RE01Trajectory.hh.

77 { return fMomentum; }
virtual G4int RE01Trajectory::GetParentID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 73 of file RE01Trajectory.hh.

73 { return fParentID; }
virtual G4String RE01Trajectory::GetParticleName ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 74 of file RE01Trajectory.hh.

74 { return fParticleName; }
virtual G4int RE01Trajectory::GetPDGEncoding ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 76 of file RE01Trajectory.hh.

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

Implements G4VTrajectory.

Definition at line 79 of file RE01Trajectory.hh.

80  { return (*fPositionRecord)[i]; }
virtual int RE01Trajectory::GetPointEntries ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 78 of file RE01Trajectory.hh.

78 { return fPositionRecord->size(); }

Here is the caller graph for this function:

virtual G4int RE01Trajectory::GetTrackID ( ) const
inlinevirtual

Implements G4VTrajectory.

Definition at line 72 of file RE01Trajectory.hh.

72 { return fTrackID; }
void RE01Trajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Implements G4VTrajectory.

Definition at line 261 of file RE01Trajectory.cc.

262 {
263  if(!secondTrajectory) return;
264 
265  RE01Trajectory* seco = (RE01Trajectory*)secondTrajectory;
266  G4int ent = seco->GetPointEntries();
267  //
268  // initial point of the second trajectory should not be merged
269  for(int i=1;i<ent;i++)
270  {
271  fPositionRecord->push_back((*(seco->fPositionRecord))[i]);
272  }
273  delete (*seco->fPositionRecord)[0];
274  seco->fPositionRecord->clear();
275 
276 }
int G4int
Definition: G4Types.hh:78
virtual int GetPointEntries() const

Here is the call graph for this function:

void RE01Trajectory::operator delete ( void aTrajectory)
inline

Definition at line 106 of file RE01Trajectory.hh.

107 {
108  myTrajectoryAllocator->FreeSingle((RE01Trajectory*)aTrajectory);
109 }
G4ThreadLocal G4Allocator< RE01Trajectory > * myTrajectoryAllocator
void * RE01Trajectory::operator new ( size_t  )
inline

Definition at line 99 of file RE01Trajectory.hh.

100 {
103  return (void*)myTrajectoryAllocator->MallocSingle();
104 }
G4ThreadLocal G4Allocator< RE01Trajectory > * myTrajectoryAllocator
int RE01Trajectory::operator== ( const RE01Trajectory right) const
inline

Definition at line 69 of file RE01Trajectory.hh.

70  {return (this==&right);}
void RE01Trajectory::ShowTrajectory ( std::ostream &  os = G4cout) const
virtual

Reimplemented from G4VTrajectory.

Definition at line 105 of file RE01Trajectory.cc.

106 {
107  os << G4endl << "TrackID =" << fTrackID
108  << " : ParentID=" << fParentID << " : TrackStatus=" << fTrackStatus << G4endl;
109  os << "Particle name : " << fParticleName << " PDG code : " << fPDGEncoding
110  << " Charge : " << fPDGCharge << G4endl;
111  os << "Original momentum : " <<
112  G4BestUnit(fMomentum,"Energy") << G4endl;
113  os << "Vertex : " << G4BestUnit(fVertexPosition,"Length")
114  << " Global time : " << G4BestUnit(fGlobalTime,"Time") << G4endl;
115  os << " Current trajectory has " << fPositionRecord->size()
116  << " points." << G4endl;
117 
118  for( size_t i=0 ; i < fPositionRecord->size() ; i++){
119  G4TrajectoryPoint* aTrajectoryPoint =
120  (G4TrajectoryPoint*)((*fPositionRecord)[i]);
121  os << "Point[" << i << "]"
122  << " Position= " << aTrajectoryPoint->GetPosition() << G4endl;
123  }
124 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
const G4ThreeVector GetPosition() const
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

Here is the caller graph for this function:


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