Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4RichTrajectory Class Reference

#include <G4RichTrajectory.hh>

Inheritance diagram for G4RichTrajectory:
G4Trajectory G4VTrajectory

Public Member Functions

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

Detailed Description

Definition at line 65 of file G4RichTrajectory.hh.

Constructor & Destructor Documentation

G4RichTrajectory::G4RichTrajectory ( )

Definition at line 62 of file G4RichTrajectory.cc.

62  :
63  fpRichPointsContainer(0),
64  fpCreatorProcess(0),
65  fpEndingProcess(0),
66  fFinalKineticEnergy(0.)
67 {
68 }
G4RichTrajectory::G4RichTrajectory ( const G4Track aTrack)

Definition at line 70 of file G4RichTrajectory.cc.

References G4Track::GetCreatorProcess(), G4Track::GetKineticEnergy(), G4Track::GetNextTouchableHandle(), and G4Track::GetTouchableHandle().

70  :
71  G4Trajectory(aTrack) // Note: this initialises the base class data
72  // members and, unfortunately but never mind,
73  // creates a G4TrajectoryPoint in
74  // TrajectoryPointContainer that we cannot
75  // access because it's private. We store the
76  // same information (plus more) in a
77  // G4RichTrajectoryPoint in the
78  // RichTrajectoryPointsContainer
79 {
80  fpInitialVolume = aTrack->GetTouchableHandle();
81  fpInitialNextVolume = aTrack->GetNextTouchableHandle();
82  fpCreatorProcess = aTrack->GetCreatorProcess();
83  // On construction, set final values to initial values.
84  // Final values are updated at the addition of every step - see AppendStep.
85  fpFinalVolume = aTrack->GetTouchableHandle();
86  fpFinalNextVolume = aTrack->GetNextTouchableHandle();
87  fpEndingProcess = aTrack->GetCreatorProcess();
88  fFinalKineticEnergy = aTrack->GetKineticEnergy();
89  // Insert the first rich trajectory point (see note above)...
90  fpRichPointsContainer = new RichTrajectoryPointsContainer;
91  fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
92 }
const G4VProcess * GetCreatorProcess() const
G4double GetKineticEnergy() const
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
const G4TouchableHandle & GetNextTouchableHandle() const
const G4TouchableHandle & GetTouchableHandle() const
G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory right)

Definition at line 94 of file G4RichTrajectory.cc.

94  :
95  G4Trajectory(right)
96 {
97  fpInitialVolume = right.fpInitialVolume;
98  fpInitialNextVolume = right.fpInitialNextVolume;
99  fpCreatorProcess = right.fpCreatorProcess;
100  fpFinalVolume = right.fpFinalVolume;
101  fpFinalNextVolume = right.fpFinalNextVolume;
102  fpEndingProcess = right.fpEndingProcess;
103  fFinalKineticEnergy = right.fFinalKineticEnergy;
104  fpRichPointsContainer = new RichTrajectoryPointsContainer;
105  for(size_t i=0;i<right.fpRichPointsContainer->size();i++)
106  {
107  G4RichTrajectoryPoint* rightPoint =
108  (G4RichTrajectoryPoint*)((*(right.fpRichPointsContainer))[i]);
109  fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
110  }
111 }
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
G4RichTrajectory::~G4RichTrajectory ( )
virtual

Definition at line 113 of file G4RichTrajectory.cc.

114 {
115  if (fpRichPointsContainer) {
116  // fpRichPointsContainer->clearAndDestroy();
117  size_t i;
118  for(i=0;i<fpRichPointsContainer->size();i++){
119  delete (*fpRichPointsContainer)[i];
120  }
121  fpRichPointsContainer->clear();
122  delete fpRichPointsContainer;
123  }
124 }

Member Function Documentation

void G4RichTrajectory::AppendStep ( const G4Step aStep)
virtual

Reimplemented from G4Trajectory.

Definition at line 126 of file G4RichTrajectory.cc.

References G4Track::GetCurrentStepNumber(), G4StepPoint::GetKineticEnergy(), G4Track::GetNextTouchableHandle(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetProcessDefinedStep(), G4Step::GetTotalEnergyDeposit(), G4Track::GetTouchableHandle(), and G4Step::GetTrack().

127 {
128  fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
129  // Except for first step, which is a sort of virtual step to start
130  // the track, compute the final values...
131  const G4Track* track = aStep->GetTrack();
132  const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
133  if (track->GetCurrentStepNumber() > 0) {
134  fpFinalVolume = track->GetTouchableHandle();
135  fpFinalNextVolume = track->GetNextTouchableHandle();
136  fpEndingProcess = postStepPoint->GetProcessDefinedStep();
137  fFinalKineticEnergy =
138  aStep->GetPreStepPoint()->GetKineticEnergy() -
139  aStep->GetTotalEnergyDeposit();
140  }
141 }
G4StepPoint * GetPreStepPoint() const
G4int GetCurrentStepNumber() const
const G4TouchableHandle & GetNextTouchableHandle() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetTotalEnergyDeposit() const
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
std::vector< G4AttValue > * G4RichTrajectory::CreateAttValues ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 237 of file G4RichTrajectory.cc.

References G4Trajectory::CreateAttValues(), G4BestUnit, G4cout, GetAttDefs(), G4VProcess::GetProcessName(), G4VProcess::GetProcessType(), G4VProcess::GetProcessTypeName(), and G4VTouchable::GetVolume().

238 {
239  // Create base class att values...
240  std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
241 
242  if (fpInitialVolume && fpInitialVolume->GetVolume()) {
243  values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
244  } else {
245  values->push_back(G4AttValue("IVPath","None",""));
246  }
247 
248  if (fpInitialNextVolume && fpInitialNextVolume->GetVolume()) {
249  values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
250  } else {
251  values->push_back(G4AttValue("INVPath","None",""));
252  }
253 
254  if (fpCreatorProcess) {
255  values->push_back(G4AttValue("CPN",fpCreatorProcess->GetProcessName(),""));
256  G4ProcessType type = fpCreatorProcess->GetProcessType();
257  values->push_back(G4AttValue("CPTN",G4VProcess::GetProcessTypeName(type),""));
258  } else {
259  values->push_back(G4AttValue("CPN","None",""));
260  values->push_back(G4AttValue("CPTN","None",""));
261  }
262 
263  if (fpFinalVolume && fpFinalVolume->GetVolume()) {
264  values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
265  } else {
266  values->push_back(G4AttValue("FVPath","None",""));
267  }
268 
269  if (fpFinalNextVolume && fpFinalNextVolume->GetVolume()) {
270  values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
271  } else {
272  values->push_back(G4AttValue("FNVPath","None",""));
273  }
274 
275  if (fpEndingProcess) {
276  values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
277  G4ProcessType type = fpEndingProcess->GetProcessType();
278  values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
279  } else {
280  values->push_back(G4AttValue("EPN","None",""));
281  values->push_back(G4AttValue("EPTN","None",""));
282  }
283 
284  values->push_back
285  (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
286 
287 #ifdef G4ATTDEBUG
288  G4cout << G4AttCheck(values,GetAttDefs());
289 #endif
290 
291  return values;
292 }
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:141
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
G4GLOB_DLL std::ostream G4cout
virtual std::vector< G4AttValue > * CreateAttValues() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4ProcessType
void G4RichTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 165 of file G4RichTrajectory.cc.

References G4VTrajectory::DrawTrajectory().

166 {
167  // Invoke the default implementation in G4VTrajectory...
169  // ... or override with your own code here.
170 }
virtual void DrawTrajectory() const
const std::map< G4String, G4AttDef > * G4RichTrajectory::GetAttDefs ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 172 of file G4RichTrajectory.cc.

References G4Trajectory::GetAttDefs(), and G4AttDefStore::GetInstance().

Referenced by CreateAttValues(), G4VisCommandList::SetNewValue(), and G4VisCommandSceneAddTrajectories::SetNewValue().

173 {
174  G4bool isNew;
175  std::map<G4String,G4AttDef>* store
176  = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
177  if (isNew) {
178 
179  // Get att defs from base class...
180  *store = *(G4Trajectory::GetAttDefs());
181 
182  G4String ID;
183 
184  ID = "IVPath";
185  (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
186  "Physics","","G4String");
187 
188  ID = "INVPath";
189  (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
190  "Physics","","G4String");
191 
192  ID = "CPN";
193  (*store)[ID] = G4AttDef(ID,"Creator Process Name",
194  "Physics","","G4String");
195 
196  ID = "CPTN";
197  (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
198  "Physics","","G4String");
199 
200  ID = "FVPath";
201  (*store)[ID] = G4AttDef(ID,"Final Volume Path",
202  "Physics","","G4String");
203 
204  ID = "FNVPath";
205  (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
206  "Physics","","G4String");
207 
208  ID = "EPN";
209  (*store)[ID] = G4AttDef(ID,"Ending Process Name",
210  "Physics","","G4String");
211 
212  ID = "EPTN";
213  (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
214  "Physics","","G4String");
215 
216  ID = "FKE";
217  (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
218  "Physics","G4BestUnit","G4double");
219 
220  }
221 
222  return store;
223 }
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
G4VTrajectoryPoint* G4RichTrajectory::GetPoint ( G4int  i) const
inlinevirtual

Reimplemented from G4Trajectory.

Definition at line 91 of file G4RichTrajectory.hh.

92  { return (*fpRichPointsContainer)[i]; }
int G4RichTrajectory::GetPointEntries ( ) const
inlinevirtual

Reimplemented from G4Trajectory.

Definition at line 90 of file G4RichTrajectory.hh.

Referenced by MergeTrajectory().

90 { return fpRichPointsContainer->size(); }
void G4RichTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Reimplemented from G4Trajectory.

Definition at line 143 of file G4RichTrajectory.cc.

References GetPointEntries().

144 {
145  if(!secondTrajectory) return;
146 
147  G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
148  G4int ent = seco->GetPointEntries();
149  for(G4int i=1;i<ent;i++) {
150  // initial point of the second trajectory should not be merged
151  fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
152  // fpRichPointsContainer->push_back(seco->fpRichPointsContainer->removeAt(1));
153  }
154  delete (*seco->fpRichPointsContainer)[0];
155  seco->fpRichPointsContainer->clear();
156 }
int GetPointEntries() const
int G4int
Definition: G4Types.hh:78
void G4RichTrajectory::operator delete ( void aRichTrajectory)
inline

Definition at line 122 of file G4RichTrajectory.hh.

References aRichTrajectoryAllocator.

123 {
124  aRichTrajectoryAllocator->FreeSingle((G4RichTrajectory*)aRichTrajectory);
125 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4RichTrajectory > * aRichTrajectoryAllocator
void * G4RichTrajectory::operator new ( size_t  )
inline

Definition at line 115 of file G4RichTrajectory.hh.

References aRichTrajectoryAllocator.

116 {
119  return (void*)aRichTrajectoryAllocator->MallocSingle();
120 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4RichTrajectory > * aRichTrajectoryAllocator
int G4RichTrajectory::operator== ( const G4RichTrajectory right) const
inline

Definition at line 82 of file G4RichTrajectory.hh.

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

Reimplemented from G4Trajectory.

Definition at line 158 of file G4RichTrajectory.cc.

References G4VTrajectory::ShowTrajectory().

159 {
160  // Invoke the default implementation in G4VTrajectory...
162  // ... or override with your own code here.
163 }
virtual void ShowTrajectory(std::ostream &os=G4cout) const

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