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

#include <G4RichTrajectoryPoint.hh>

Inheritance diagram for G4RichTrajectoryPoint:
Collaboration diagram for G4RichTrajectoryPoint:

Public Member Functions

 G4RichTrajectoryPoint ()
 
 G4RichTrajectoryPoint (const G4Track *)
 
 G4RichTrajectoryPoint (const G4Step *)
 
 G4RichTrajectoryPoint (const G4RichTrajectoryPoint &right)
 
virtual ~G4RichTrajectoryPoint ()
 
const std::vector
< G4ThreeVector > * 
GetAuxiliaryPoints () const
 
voidoperator new (size_t)
 
void operator delete (void *aRichTrajectoryPoint)
 
int operator== (const G4RichTrajectoryPoint &right) const
 
virtual const std::map
< G4String, G4AttDef > * 
GetAttDefs () const
 
virtual std::vector< G4AttValue > * CreateAttValues () const
 
- Public Member Functions inherited from G4TrajectoryPoint
 G4TrajectoryPoint ()
 
 G4TrajectoryPoint (G4ThreeVector pos)
 
 G4TrajectoryPoint (const G4TrajectoryPoint &right)
 
virtual ~G4TrajectoryPoint ()
 
voidoperator new (size_t)
 
void operator delete (void *aTrajectoryPoint)
 
int operator== (const G4TrajectoryPoint &right) const
 
const G4ThreeVector GetPosition () const
 
- Public Member Functions inherited from G4VTrajectoryPoint
 G4VTrajectoryPoint ()
 
virtual ~G4VTrajectoryPoint ()
 
G4bool operator== (const G4VTrajectoryPoint &right) const
 

Detailed Description

Definition at line 72 of file G4RichTrajectoryPoint.hh.

Constructor & Destructor Documentation

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( )

Definition at line 66 of file G4RichTrajectoryPoint.cc.

66  :
67  fpAuxiliaryPointVector(0),
68  fTotEDep(0.),
69  fRemainingEnergy(0.),
70  fpProcess(0),
71  fPreStepPointStatus(fUndefined),
72  fPostStepPointStatus(fUndefined),
73  fPreStepPointGlobalTime(0),
74  fPostStepPointGlobalTime(0),
75  fPreStepPointWeight(1.),
76  fPostStepPointWeight(1.)
77 {}
G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Track aTrack)

Definition at line 79 of file G4RichTrajectoryPoint.cc.

79  :
80  G4TrajectoryPoint(aTrack->GetPosition()),
81  fpAuxiliaryPointVector(0),
82  fTotEDep(0.),
83  fRemainingEnergy(aTrack->GetKineticEnergy()),
84  fpProcess(0),
85  fPreStepPointStatus(fUndefined),
86  fPostStepPointStatus(fUndefined),
87  fPreStepPointGlobalTime(aTrack->GetGlobalTime()),
88  fPostStepPointGlobalTime(aTrack->GetGlobalTime()),
89  fpPreStepPointVolume(aTrack->GetTouchableHandle()),
90  fpPostStepPointVolume(aTrack->GetNextTouchableHandle()),
91  fPreStepPointWeight(aTrack->GetWeight()),
92  fPostStepPointWeight(aTrack->GetWeight())
93 {}
const G4ThreeVector & GetPosition() const
G4double GetKineticEnergy() const
const G4TouchableHandle & GetNextTouchableHandle() const
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetWeight() const
G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Step aStep)

Definition at line 95 of file G4RichTrajectoryPoint.cc.

95  :
97  fpAuxiliaryPointVector(aStep->GetPointerToVectorOfAuxiliaryPoints()),
98  fTotEDep(aStep->GetTotalEnergyDeposit())
99 {
100  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
101  G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
102  if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) { // First step
103  fRemainingEnergy = aStep->GetTrack()->GetKineticEnergy();
104  } else {
105  fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep;
106  }
107  fpProcess = postStepPoint->GetProcessDefinedStep();
108  fPreStepPointStatus = preStepPoint->GetStepStatus();
109  fPostStepPointStatus = postStepPoint->GetStepStatus();
110  fPreStepPointGlobalTime = preStepPoint->GetGlobalTime();
111  fPostStepPointGlobalTime = postStepPoint->GetGlobalTime();
112  fpPreStepPointVolume = preStepPoint->GetTouchableHandle();
113  fpPostStepPointVolume = postStepPoint->GetTouchableHandle();
114  fPreStepPointWeight = preStepPoint->GetWeight();
115  fPostStepPointWeight = postStepPoint->GetWeight();
116 }
G4double GetWeight() const
G4StepStatus GetStepStatus() const
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
G4int GetCurrentStepNumber() const
const G4ThreeVector & GetPosition() const
G4double GetTotalEnergyDeposit() const
const G4VProcess * GetProcessDefinedStep() const
G4StepPoint * GetPostStepPoint() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
Definition: G4Step.hh:242
G4double GetGlobalTime() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
const G4TouchableHandle & GetTouchableHandle() const

Here is the call graph for this function:

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4RichTrajectoryPoint right)

Definition at line 119 of file G4RichTrajectoryPoint.cc.

119  :
120  G4TrajectoryPoint(right),
121  fpAuxiliaryPointVector(right.fpAuxiliaryPointVector),
122  fTotEDep(right.fTotEDep),
123  fRemainingEnergy(right.fRemainingEnergy),
124  fpProcess(right.fpProcess),
125  fPreStepPointStatus(right.fPreStepPointStatus),
126  fPostStepPointStatus(right.fPostStepPointStatus),
127  fPreStepPointGlobalTime(right.fPreStepPointGlobalTime),
128  fPostStepPointGlobalTime(right.fPostStepPointGlobalTime),
129  fpPreStepPointVolume(right.fpPreStepPointVolume),
130  fpPostStepPointVolume(right.fpPostStepPointVolume),
131  fPreStepPointWeight(right.fPreStepPointWeight),
132  fPostStepPointWeight(right.fPostStepPointWeight)
133 {}
G4RichTrajectoryPoint::~G4RichTrajectoryPoint ( )
virtual

Definition at line 135 of file G4RichTrajectoryPoint.cc.

136 {
137  if(fpAuxiliaryPointVector) {
138  /*
139  G4cout << "Deleting fpAuxiliaryPointVector at "
140  << (void*) fpAuxiliaryPointVector
141  << G4endl;
142  */
143  delete fpAuxiliaryPointVector;
144  }
145 }

Member Function Documentation

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

Reimplemented from G4TrajectoryPoint.

Definition at line 232 of file G4RichTrajectoryPoint.cc.

233 {
234  // Create base class att values...
235  std::vector<G4AttValue>* values = G4TrajectoryPoint::CreateAttValues();
236 
237  if (fpAuxiliaryPointVector) {
238  std::vector<G4ThreeVector>::iterator iAux;
239  for (iAux = fpAuxiliaryPointVector->begin();
240  iAux != fpAuxiliaryPointVector->end(); ++iAux) {
241  values->push_back(G4AttValue("Aux",G4BestUnit(*iAux,"Length"),""));
242  }
243  }
244 
245  values->push_back(G4AttValue("TED",G4BestUnit(fTotEDep,"Energy"),""));
246 
247  values->push_back(G4AttValue("RE",G4BestUnit(fRemainingEnergy,"Energy"),""));
248 
249  if (fpProcess) {
250  values->push_back
251  (G4AttValue("PDS",fpProcess->GetProcessName(),""));
252  values->push_back
253  (G4AttValue
254  ("PTDS",G4VProcess::GetProcessTypeName(fpProcess->GetProcessType()),
255  ""));
256  } else {
257  values->push_back(G4AttValue("PDS","None",""));
258  values->push_back(G4AttValue("PTDS","None",""));
259  }
260 
261  values->push_back
262  (G4AttValue("PreStatus",Status(fPreStepPointStatus),""));
263 
264  values->push_back
265  (G4AttValue("PostStatus",Status(fPostStepPointStatus),""));
266 
267  values->push_back
268  (G4AttValue("PreT",G4BestUnit(fPreStepPointGlobalTime,"Time"),""));
269 
270  values->push_back
271  (G4AttValue("PostT",G4BestUnit(fPostStepPointGlobalTime,"Time"),""));
272 
273  if (fpPreStepPointVolume && fpPreStepPointVolume->GetVolume()) {
274  values->push_back(G4AttValue("PreVPath",Path(fpPreStepPointVolume),""));
275  } else {
276  values->push_back(G4AttValue("PreVPath","None",""));
277  }
278 
279  if (fpPostStepPointVolume && fpPostStepPointVolume->GetVolume()) {
280  values->push_back(G4AttValue("PostVPath",Path(fpPostStepPointVolume),""));
281  } else {
282  values->push_back(G4AttValue("PostVPath","None",""));
283  }
284 
285  {
286  std::ostringstream oss;
287  oss << fPreStepPointWeight;
288  values->push_back
289  (G4AttValue("PreW",oss.str(),""));
290  }
291 
292  {
293  std::ostringstream oss;
294  oss << fPostStepPointWeight;
295  values->push_back
296  (G4AttValue("PostW",oss.str(),""));
297  }
298 
299 #ifdef G4ATTDEBUG
300  G4cout << G4AttCheck(values,GetAttDefs());
301 #endif
302 
303  return values;
304 }
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:141
virtual std::vector< G4AttValue > * CreateAttValues() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:414
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
static G4String Path(const G4TouchableHandle &th)
static G4String Status(G4StepStatus stps)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const

Here is the call graph for this function:

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

Reimplemented from G4TrajectoryPoint.

Definition at line 148 of file G4RichTrajectoryPoint.cc.

149 {
150  G4bool isNew;
151  std::map<G4String,G4AttDef>* store
152  = G4AttDefStore::GetInstance("G4RichTrajectoryPoint",isNew);
153  if (isNew) {
154 
155  // Copy base class att defs...
156  *store = *(G4TrajectoryPoint::GetAttDefs());
157 
158  G4String ID;
159 
160  ID = "Aux";
161  (*store)[ID] = G4AttDef(ID, "Auxiliary Point Position",
162  "Physics","G4BestUnit","G4ThreeVector");
163  ID = "TED";
164  (*store)[ID] = G4AttDef(ID,"Total Energy Deposit",
165  "Physics","G4BestUnit","G4double");
166  ID = "RE";
167  (*store)[ID] = G4AttDef(ID,"Remaining Energy",
168  "Physics","G4BestUnit","G4double");
169  ID = "PDS";
170  (*store)[ID] = G4AttDef(ID,"Process Defined Step",
171  "Physics","","G4String");
172  ID = "PTDS";
173  (*store)[ID] = G4AttDef(ID,"Process Type Defined Step",
174  "Physics","","G4String");
175  ID = "PreStatus";
176  (*store)[ID] = G4AttDef(ID,"Pre-step-point status",
177  "Physics","","G4String");
178  ID = "PostStatus";
179  (*store)[ID] = G4AttDef(ID,"Post-step-point status",
180  "Physics","","G4String");
181  ID = "PreT";
182  (*store)[ID] = G4AttDef(ID,"Pre-step-point global time",
183  "Physics","G4BestUnit","G4double");
184  ID = "PostT";
185  (*store)[ID] = G4AttDef(ID,"Post-step-point global time",
186  "Physics","G4BestUnit","G4double");
187  ID = "PreVPath";
188  (*store)[ID] = G4AttDef(ID,"Pre-step Volume Path",
189  "Physics","","G4String");
190  ID = "PostVPath";
191  (*store)[ID] = G4AttDef(ID,"Post-step Volume Path",
192  "Physics","","G4String");
193  ID = "PreW";
194  (*store)[ID] = G4AttDef(ID,"Pre-step-point weight",
195  "Physics","","G4double");
196  ID = "PostW";
197  (*store)[ID] = G4AttDef(ID,"Post-step-point weight",
198  "Physics","","G4double");
199  }
200  return store;
201 }
bool G4bool
Definition: G4Types.hh:79
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const

Here is the call graph for this function:

Here is the caller graph for this function:

const std::vector<G4ThreeVector>* G4RichTrajectoryPoint::GetAuxiliaryPoints ( ) const
inlinevirtual

Reimplemented from G4VTrajectoryPoint.

Definition at line 90 of file G4RichTrajectoryPoint.hh.

91  { return fpAuxiliaryPointVector; }
void G4RichTrajectoryPoint::operator delete ( void aRichTrajectoryPoint)
inline

Definition at line 130 of file G4RichTrajectoryPoint.hh.

131 {
133  ((G4RichTrajectoryPoint *) aRichTrajectoryPoint);
134 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4RichTrajectoryPoint > * aRichTrajectoryPointAllocator
void * G4RichTrajectoryPoint::operator new ( size_t  )
inline

Definition at line 123 of file G4RichTrajectoryPoint.hh.

124 {
127  return (void *) aRichTrajectoryPointAllocator->MallocSingle();
128 }
G4TRACKING_DLL G4ThreadLocal G4Allocator< G4RichTrajectoryPoint > * aRichTrajectoryPointAllocator
int G4RichTrajectoryPoint::operator== ( const G4RichTrajectoryPoint right) const
inline

Definition at line 96 of file G4RichTrajectoryPoint.hh.

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

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