Geant4  10.02.p03
G4RichTrajectory Class Reference

#include <G4RichTrajectory.hh>

Inheritance diagram for G4RichTrajectory:
Collaboration diagram for G4RichTrajectory:

Public Member Functions

 G4RichTrajectory ()
 
 G4RichTrajectory (const G4Track *aTrack)
 
 G4RichTrajectory (G4RichTrajectory &)
 
virtual ~G4RichTrajectory ()
 
void * operator 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 ()
 
void * operator 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
 

Private Member Functions

G4RichTrajectoryoperator= (const G4RichTrajectory &)
 

Private Attributes

RichTrajectoryPointsContainerfpRichPointsContainer
 
G4TouchableHandle fpInitialVolume
 
G4TouchableHandle fpInitialNextVolume
 
const G4VProcessfpCreatorProcess
 
G4int fCreatorModelID
 
G4TouchableHandle fpFinalVolume
 
G4TouchableHandle fpFinalNextVolume
 
const G4VProcessfpEndingProcess
 
G4double fFinalKineticEnergy
 

Detailed Description

Definition at line 65 of file G4RichTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectory() [1/3]

G4RichTrajectory::G4RichTrajectory ( )

Definition at line 63 of file G4RichTrajectory.cc.

63  :
66  fCreatorModelID(0),
67  fpEndingProcess(0),
69 {
70 }
G4double fFinalKineticEnergy
RichTrajectoryPointsContainer * fpRichPointsContainer
const G4VProcess * fpEndingProcess
const G4VProcess * fpCreatorProcess

◆ G4RichTrajectory() [2/3]

G4RichTrajectory::G4RichTrajectory ( const G4Track *  aTrack)

Definition at line 72 of file G4RichTrajectory.cc.

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

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory right)

Definition at line 97 of file G4RichTrajectory.cc.

97  :
98  G4Trajectory(right)
99 {
109  for(size_t i=0;i<right.fpRichPointsContainer->size();i++)
110  {
111  G4RichTrajectoryPoint* rightPoint =
113  fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
114  }
115 }
G4TouchableHandle fpInitialNextVolume
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
G4double fFinalKineticEnergy
G4TouchableHandle fpInitialVolume
RichTrajectoryPointsContainer * fpRichPointsContainer
G4TouchableHandle fpFinalVolume
G4TouchableHandle fpFinalNextVolume
const G4VProcess * fpEndingProcess
const G4VProcess * fpCreatorProcess

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
virtual

Definition at line 117 of file G4RichTrajectory.cc.

118 {
119  if (fpRichPointsContainer) {
120  // fpRichPointsContainer->clearAndDestroy();
121  size_t i;
122  for(i=0;i<fpRichPointsContainer->size();i++){
123  delete (*fpRichPointsContainer)[i];
124  }
125  fpRichPointsContainer->clear();
126  delete fpRichPointsContainer;
127  }
128 }
RichTrajectoryPointsContainer * fpRichPointsContainer

Member Function Documentation

◆ AppendStep()

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

Reimplemented from G4Trajectory.

Definition at line 130 of file G4RichTrajectory.cc.

131 {
132  fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
133  // Except for first step, which is a sort of virtual step to start
134  // the track, compute the final values...
135  const G4Track* track = aStep->GetTrack();
136  const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
137  if (track->GetCurrentStepNumber() > 0) {
138  fpFinalVolume = track->GetTouchableHandle();
139  fpFinalNextVolume = track->GetNextTouchableHandle();
140  fpEndingProcess = postStepPoint->GetProcessDefinedStep();
142  aStep->GetPreStepPoint()->GetKineticEnergy() -
143  aStep->GetTotalEnergyDeposit();
144  }
145 }
G4double fFinalKineticEnergy
RichTrajectoryPointsContainer * fpRichPointsContainer
G4TouchableHandle fpFinalVolume
G4TouchableHandle fpFinalNextVolume
const G4VProcess * fpEndingProcess
Here is the caller graph for this function:

◆ CreateAttValues()

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

Reimplemented from G4Trajectory.

Definition at line 249 of file G4RichTrajectory.cc.

250 {
251  // Create base class att values...
252  std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
253 
255  values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
256  } else {
257  values->push_back(G4AttValue("IVPath","None",""));
258  }
259 
261  values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
262  } else {
263  values->push_back(G4AttValue("INVPath","None",""));
264  }
265 
266  if (fpCreatorProcess) {
267  values->push_back
270  values->push_back
271  (G4AttValue("CPTN",G4VProcess::GetProcessTypeName(type),""));
272  values->push_back
274  const G4String& creatorModelName =
276  values->push_back(G4AttValue("CMN",creatorModelName,""));
277  } else {
278  values->push_back(G4AttValue("CPN","None",""));
279  values->push_back(G4AttValue("CPTN","None",""));
280  values->push_back(G4AttValue("CMID","None",""));
281  values->push_back(G4AttValue("CMN","None",""));
282  }
283 
285  values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
286  } else {
287  values->push_back(G4AttValue("FVPath","None",""));
288  }
289 
291  values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
292  } else {
293  values->push_back(G4AttValue("FNVPath","None",""));
294  }
295 
296  if (fpEndingProcess) {
297  values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
299  values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
300  } else {
301  values->push_back(G4AttValue("EPN","None",""));
302  values->push_back(G4AttValue("EPTN","None",""));
303  }
304 
305  values->push_back
306  (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
307 
308 #ifdef G4ATTDEBUG
309  G4cout << G4AttCheck(values,GetAttDefs());
310 #endif
311 
312  return values;
313 }
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:141
virtual std::vector< G4AttValue > * CreateAttValues() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:371
static const G4String & GetModelName(G4int)
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
static G4String Path(const G4TouchableHandle &th)
G4TouchableHandle fpInitialNextVolume
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4GLOB_DLL std::ostream G4cout
G4double fFinalKineticEnergy
G4TouchableHandle fpInitialVolume
G4TouchableHandle fpFinalVolume
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
G4TouchableHandle fpFinalNextVolume
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
const G4VProcess * fpEndingProcess
const G4VProcess * fpCreatorProcess
G4ProcessType
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DrawTrajectory()

void G4RichTrajectory::DrawTrajectory ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 169 of file G4RichTrajectory.cc.

170 {
171  // Invoke the default implementation in G4VTrajectory...
173  // ... or override with your own code here.
174 }
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 > * G4RichTrajectory::GetAttDefs ( ) const
virtual

Reimplemented from G4Trajectory.

Definition at line 176 of file G4RichTrajectory.cc.

177 {
178  G4bool isNew;
179  std::map<G4String,G4AttDef>* store
180  = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
181  if (isNew) {
182 
183  // Get att defs from base class...
184  *store = *(G4Trajectory::GetAttDefs());
185 
186  G4String ID;
187 
188  ID = "IVPath";
189  (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
190  "Physics","","G4String");
191 
192  ID = "INVPath";
193  (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
194  "Physics","","G4String");
195 
196  ID = "CPN";
197  (*store)[ID] = G4AttDef(ID,"Creator Process Name",
198  "Physics","","G4String");
199 
200  ID = "CPTN";
201  (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
202  "Physics","","G4String");
203 
204  ID = "CMID";
205  (*store)[ID] = G4AttDef(ID,"Creator Model ID",
206  "Physics","","G4int");
207 
208  ID = "CMN";
209  (*store)[ID] = G4AttDef(ID,"Creator Model Name",
210  "Physics","","G4String");
211 
212  ID = "FVPath";
213  (*store)[ID] = G4AttDef(ID,"Final Volume Path",
214  "Physics","","G4String");
215 
216  ID = "FNVPath";
217  (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
218  "Physics","","G4String");
219 
220  ID = "EPN";
221  (*store)[ID] = G4AttDef(ID,"Ending Process Name",
222  "Physics","","G4String");
223 
224  ID = "EPTN";
225  (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
226  "Physics","","G4String");
227 
228  ID = "FKE";
229  (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
230  "Physics","G4BestUnit","G4double");
231 
232  }
233 
234  return store;
235 }
bool G4bool
Definition: G4Types.hh:79
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
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:

◆ GetPoint()

G4VTrajectoryPoint* G4RichTrajectory::GetPoint ( G4int  i) const
inlinevirtual

Reimplemented from G4Trajectory.

Definition at line 91 of file G4RichTrajectory.hh.

92  { return (*fpRichPointsContainer)[i]; }
RichTrajectoryPointsContainer * fpRichPointsContainer
Here is the call graph for this function:

◆ GetPointEntries()

int G4RichTrajectory::GetPointEntries ( ) const
inlinevirtual

Reimplemented from G4Trajectory.

Definition at line 90 of file G4RichTrajectory.hh.

90 { return fpRichPointsContainer->size(); }
RichTrajectoryPointsContainer * fpRichPointsContainer
Here is the caller graph for this function:

◆ MergeTrajectory()

void G4RichTrajectory::MergeTrajectory ( G4VTrajectory secondTrajectory)
virtual

Reimplemented from G4Trajectory.

Definition at line 147 of file G4RichTrajectory.cc.

148 {
149  if(!secondTrajectory) return;
150 
151  G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
152  G4int ent = seco->GetPointEntries();
153  for(G4int i=1;i<ent;i++) {
154  // initial point of the second trajectory should not be merged
155  fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
156  // fpRichPointsContainer->push_back(seco->fpRichPointsContainer->removeAt(1));
157  }
158  delete (*seco->fpRichPointsContainer)[0];
159  seco->fpRichPointsContainer->clear();
160 }
int G4int
Definition: G4Types.hh:78
int GetPointEntries() const
RichTrajectoryPointsContainer * fpRichPointsContainer
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator delete()

void G4RichTrajectory::operator delete ( void *  aRichTrajectory)
inline

Definition at line 123 of file G4RichTrajectory.hh.

124 {
125  aRichTrajectoryAllocator->FreeSingle((G4RichTrajectory*)aRichTrajectory);
126 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4RichTrajectory > * aRichTrajectoryAllocator

◆ operator new()

void * G4RichTrajectory::operator new ( size_t  )
inline

Definition at line 116 of file G4RichTrajectory.hh.

117 {
120  return (void*)aRichTrajectoryAllocator->MallocSingle();
121 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4RichTrajectory > * aRichTrajectoryAllocator

◆ operator=()

G4RichTrajectory& G4RichTrajectory::operator= ( const G4RichTrajectory )
private

◆ operator==()

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

Definition at line 82 of file G4RichTrajectory.hh.

83  {return (this==&right);}
Here is the call graph for this function:

◆ ShowTrajectory()

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

Reimplemented from G4Trajectory.

Definition at line 162 of file G4RichTrajectory.cc.

163 {
164  // Invoke the default implementation in G4VTrajectory...
166  // ... or override with your own code here.
167 }
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

◆ fCreatorModelID

G4int G4RichTrajectory::fCreatorModelID
private

Definition at line 105 of file G4RichTrajectory.hh.

◆ fFinalKineticEnergy

G4double G4RichTrajectory::fFinalKineticEnergy
private

Definition at line 109 of file G4RichTrajectory.hh.

◆ fpCreatorProcess

const G4VProcess* G4RichTrajectory::fpCreatorProcess
private

Definition at line 104 of file G4RichTrajectory.hh.

◆ fpEndingProcess

const G4VProcess* G4RichTrajectory::fpEndingProcess
private

Definition at line 108 of file G4RichTrajectory.hh.

◆ fpFinalNextVolume

G4TouchableHandle G4RichTrajectory::fpFinalNextVolume
private

Definition at line 107 of file G4RichTrajectory.hh.

◆ fpFinalVolume

G4TouchableHandle G4RichTrajectory::fpFinalVolume
private

Definition at line 106 of file G4RichTrajectory.hh.

◆ fpInitialNextVolume

G4TouchableHandle G4RichTrajectory::fpInitialNextVolume
private

Definition at line 103 of file G4RichTrajectory.hh.

◆ fpInitialVolume

G4TouchableHandle G4RichTrajectory::fpInitialVolume
private

Definition at line 102 of file G4RichTrajectory.hh.

◆ fpRichPointsContainer

RichTrajectoryPointsContainer* G4RichTrajectory::fpRichPointsContainer
private

Definition at line 101 of file G4RichTrajectory.hh.


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