Geant4-11
Public Member Functions | Private Attributes
G4RichTrajectoryPoint Class Reference

#include <G4RichTrajectoryPoint.hh>

Inheritance diagram for G4RichTrajectoryPoint:
G4TrajectoryPoint G4VTrajectoryPoint

Public Member Functions

virtual std::vector< G4AttValue > * CreateAttValues () const
 
 G4RichTrajectoryPoint ()
 
 G4RichTrajectoryPoint (const G4RichTrajectoryPoint &right)
 
 G4RichTrajectoryPoint (const G4Step *)
 
 G4RichTrajectoryPoint (const G4Track *)
 
virtual const std::map< G4String, G4AttDef > * GetAttDefs () const
 
const std::vector< G4ThreeVector > * GetAuxiliaryPoints () const
 
const G4ThreeVector GetPosition () const
 
void operator delete (void *aRichTrajectoryPoint)
 
void * operator new (size_t)
 
G4RichTrajectoryPointoperator= (const G4RichTrajectoryPoint &)=delete
 
G4bool operator== (const G4RichTrajectoryPoint &right) const
 
G4bool operator== (const G4TrajectoryPoint &right) const
 
G4bool operator== (const G4VTrajectoryPoint &right) const
 
virtual ~G4RichTrajectoryPoint ()
 

Private Attributes

std::vector< G4ThreeVector > * fpAuxiliaryPointVector = nullptr
 
G4ThreeVector fPosition
 
G4double fPostStepPointGlobalTime = 0.0
 
G4StepStatus fPostStepPointStatus = fUndefined
 
G4double fPostStepPointWeight = 0.0
 
G4TouchableHandle fpPostStepPointVolume
 
G4TouchableHandle fpPreStepPointVolume
 
const G4VProcessfpProcess = nullptr
 
G4double fPreStepPointGlobalTime = 0.0
 
G4StepStatus fPreStepPointStatus = fUndefined
 
G4double fPreStepPointWeight = 0.0
 
G4double fRemainingEnergy = 0.0
 
G4double fTotEDep = 0.0
 

Detailed Description

Definition at line 66 of file G4RichTrajectoryPoint.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectoryPoint() [1/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( )

Definition at line 63 of file G4RichTrajectoryPoint.cc.

64{
65}

◆ G4RichTrajectoryPoint() [2/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Track aTrack)

Definition at line 67 of file G4RichTrajectoryPoint.cc.

68 : G4TrajectoryPoint(aTrack->GetPosition()),
75{
76}
G4TouchableHandle fpPreStepPointVolume
G4TouchableHandle fpPostStepPointVolume
const G4TouchableHandle & GetNextTouchableHandle() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const

◆ G4RichTrajectoryPoint() [3/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Step aStep)

Definition at line 78 of file G4RichTrajectoryPoint.cc.

82{
83 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
84 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
85 if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) // First step
86 {
88 }
89 else
90 {
91 fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep;
92 }
93 fpProcess = postStepPoint->GetProcessDefinedStep();
94 fPreStepPointStatus = preStepPoint->GetStepStatus();
95 fPostStepPointStatus = postStepPoint->GetStepStatus();
96 fPreStepPointGlobalTime = preStepPoint->GetGlobalTime();
97 fPostStepPointGlobalTime = postStepPoint->GetGlobalTime();
100 fPreStepPointWeight = preStepPoint->GetWeight();
101 fPostStepPointWeight = postStepPoint->GetWeight();
102}
const G4VProcess * fpProcess
std::vector< G4ThreeVector > * fpAuxiliaryPointVector
G4StepStatus GetStepStatus() const
G4double GetGlobalTime() const
const G4VProcess * GetProcessDefinedStep() const
const G4ThreeVector & GetPosition() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4double GetWeight() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4int GetCurrentStepNumber() const
G4double GetKineticEnergy() const

References fPostStepPointGlobalTime, fPostStepPointStatus, fPostStepPointWeight, fpPostStepPointVolume, fpPreStepPointVolume, fpProcess, fPreStepPointGlobalTime, fPreStepPointStatus, fPreStepPointWeight, fRemainingEnergy, fTotEDep, G4Track::GetCurrentStepNumber(), G4StepPoint::GetGlobalTime(), G4StepPoint::GetKineticEnergy(), G4Track::GetKineticEnergy(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4StepPoint::GetProcessDefinedStep(), G4StepPoint::GetStepStatus(), G4StepPoint::GetTouchableHandle(), G4Step::GetTrack(), and G4StepPoint::GetWeight().

◆ G4RichTrajectoryPoint() [4/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4RichTrajectoryPoint right)

Definition at line 104 of file G4RichTrajectoryPoint.cc.

106 fpAuxiliaryPointVector(r.fpAuxiliaryPointVector),
107 fTotEDep(r.fTotEDep),
108 fRemainingEnergy(r.fRemainingEnergy),
109 fpProcess(r.fpProcess),
110 fPreStepPointStatus(r.fPreStepPointStatus),
111 fPostStepPointStatus(r.fPostStepPointStatus),
112 fPreStepPointGlobalTime(r.fPreStepPointGlobalTime),
113 fPostStepPointGlobalTime(r.fPostStepPointGlobalTime),
114 fpPreStepPointVolume(r.fpPreStepPointVolume),
115 fpPostStepPointVolume(r.fpPostStepPointVolume),
116 fPreStepPointWeight(r.fPreStepPointWeight),
117 fPostStepPointWeight(r.fPostStepPointWeight)
118{
119}

◆ ~G4RichTrajectoryPoint()

G4RichTrajectoryPoint::~G4RichTrajectoryPoint ( )
virtual

Definition at line 121 of file G4RichTrajectoryPoint.cc.

122{
124}

References fpAuxiliaryPointVector.

Member Function Documentation

◆ CreateAttValues()

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

Reimplemented from G4TrajectoryPoint.

Definition at line 212 of file G4RichTrajectoryPoint.cc.

213{
214 // Create base class att values...
215
216 std::vector<G4AttValue>* values = G4TrajectoryPoint::CreateAttValues();
217
218 if (fpAuxiliaryPointVector != nullptr)
219 {
220 for (auto iAux = fpAuxiliaryPointVector->cbegin();
221 iAux != fpAuxiliaryPointVector->cend(); ++iAux)
222 {
223 values->push_back(G4AttValue("Aux",G4BestUnit(*iAux,"Length"),""));
224 }
225 }
226
227 values->push_back(G4AttValue("TED",G4BestUnit(fTotEDep,"Energy"),""));
228 values->push_back(G4AttValue("RE",G4BestUnit(fRemainingEnergy,"Energy"),""));
229
230 if (fpProcess != nullptr)
231 {
232 values->push_back
233 (G4AttValue("PDS",fpProcess->GetProcessName(),""));
234 values->push_back
237 }
238 else
239 {
240 values->push_back(G4AttValue("PDS","None",""));
241 values->push_back(G4AttValue("PTDS","None",""));
242 }
243
244 values->push_back
245 (G4AttValue("PreStatus",Status(fPreStepPointStatus),""));
246
247 values->push_back
248 (G4AttValue("PostStatus",Status(fPostStepPointStatus),""));
249
250 values->push_back
251 (G4AttValue("PreT",G4BestUnit(fPreStepPointGlobalTime,"Time"),""));
252
253 values->push_back
254 (G4AttValue("PostT",G4BestUnit(fPostStepPointGlobalTime,"Time"),""));
255
257 {
258 values->push_back(G4AttValue("PreVPath",Path(fpPreStepPointVolume),""));
259 }
260 else
261 {
262 values->push_back(G4AttValue("PreVPath","None",""));
263 }
264
266 {
267 values->push_back(G4AttValue("PostVPath",Path(fpPostStepPointVolume),""));
268 }
269 else
270 {
271 values->push_back(G4AttValue("PostVPath","None",""));
272 }
273
274 std::ostringstream oss1;
275 oss1 << fPreStepPointWeight;
276 values->push_back(G4AttValue("PreW",oss1.str(),""));
277
278 std::ostringstream oss2;
279 oss2 << fPostStepPointWeight;
280 values->push_back(G4AttValue("PostW",oss2.str(),""));
281
282#ifdef G4ATTDEBUG
283 G4cout << G4AttCheck(values,GetAttDefs());
284#endif
285
286 return values;
287}
static G4String Path(const G4TouchableHandle &th)
static G4String Status(G4StepStatus stps)
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:134
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:388
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34

References G4TrajectoryPoint::CreateAttValues(), fpAuxiliaryPointVector, fPostStepPointGlobalTime, fPostStepPointStatus, fPostStepPointWeight, fpPostStepPointVolume, fpPreStepPointVolume, fpProcess, fPreStepPointGlobalTime, fPreStepPointStatus, fPreStepPointWeight, fRemainingEnergy, fTotEDep, G4BestUnit, G4cout, GetAttDefs(), G4VProcess::GetProcessName(), G4VProcess::GetProcessType(), G4VProcess::GetProcessTypeName(), G4VTouchable::GetVolume(), Path(), and Status().

◆ GetAttDefs()

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

Reimplemented from G4TrajectoryPoint.

Definition at line 126 of file G4RichTrajectoryPoint.cc.

127{
128 G4bool isNew;
129 std::map<G4String,G4AttDef>* store
130 = G4AttDefStore::GetInstance("G4RichTrajectoryPoint",isNew);
131 if (isNew)
132 {
133 // Copy base class att defs...
134 *store = *(G4TrajectoryPoint::GetAttDefs());
135
136 G4String ID;
137
138 ID = "Aux";
139 (*store)[ID] = G4AttDef(ID, "Auxiliary Point Position",
140 "Physics","G4BestUnit","G4ThreeVector");
141 ID = "TED";
142 (*store)[ID] = G4AttDef(ID,"Total Energy Deposit",
143 "Physics","G4BestUnit","G4double");
144 ID = "RE";
145 (*store)[ID] = G4AttDef(ID,"Remaining Energy",
146 "Physics","G4BestUnit","G4double");
147 ID = "PDS";
148 (*store)[ID] = G4AttDef(ID,"Process Defined Step",
149 "Physics","","G4String");
150 ID = "PTDS";
151 (*store)[ID] = G4AttDef(ID,"Process Type Defined Step",
152 "Physics","","G4String");
153 ID = "PreStatus";
154 (*store)[ID] = G4AttDef(ID,"Pre-step-point status",
155 "Physics","","G4String");
156 ID = "PostStatus";
157 (*store)[ID] = G4AttDef(ID,"Post-step-point status",
158 "Physics","","G4String");
159 ID = "PreT";
160 (*store)[ID] = G4AttDef(ID,"Pre-step-point global time",
161 "Physics","G4BestUnit","G4double");
162 ID = "PostT";
163 (*store)[ID] = G4AttDef(ID,"Post-step-point global time",
164 "Physics","G4BestUnit","G4double");
165 ID = "PreVPath";
166 (*store)[ID] = G4AttDef(ID,"Pre-step Volume Path",
167 "Physics","","G4String");
168 ID = "PostVPath";
169 (*store)[ID] = G4AttDef(ID,"Post-step Volume Path",
170 "Physics","","G4String");
171 ID = "PreW";
172 (*store)[ID] = G4AttDef(ID,"Pre-step-point weight",
173 "Physics","","G4double");
174 ID = "PostW";
175 (*store)[ID] = G4AttDef(ID,"Post-step-point weight",
176 "Physics","","G4double");
177 }
178 return store;
179}
bool G4bool
Definition: G4Types.hh:86
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

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

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

◆ GetAuxiliaryPoints()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 137 of file G4RichTrajectoryPoint.hh.

138{
140}

References fpAuxiliaryPointVector.

◆ GetPosition()

const G4ThreeVector G4TrajectoryPoint::GetPosition ( void  ) const
inlinevirtualinherited

Implements G4VTrajectoryPoint.

Definition at line 65 of file G4TrajectoryPoint.hh.

66 { return fPosition; }
G4ThreeVector fPosition

References G4TrajectoryPoint::fPosition.

◆ operator delete()

void G4RichTrajectoryPoint::operator delete ( void *  aRichTrajectoryPoint)
inline

Definition at line 124 of file G4RichTrajectoryPoint.hh.

125{
127 ((G4RichTrajectoryPoint *) aRichTrajectoryPoint);
128}
G4TRACKING_DLL G4Allocator< G4RichTrajectoryPoint > *& aRichTrajectoryPointAllocator()

References aRichTrajectoryPointAllocator().

◆ operator new()

void * G4RichTrajectoryPoint::operator new ( size_t  )
inline

Definition at line 114 of file G4RichTrajectoryPoint.hh.

115{
116 if (aRichTrajectoryPointAllocator() == nullptr)
117 {
119 }
120 return (void *) aRichTrajectoryPointAllocator()->MallocSingle();
121}

References aRichTrajectoryPointAllocator().

◆ operator=()

G4RichTrajectoryPoint & G4RichTrajectoryPoint::operator= ( const G4RichTrajectoryPoint )
delete

◆ operator==() [1/3]

G4bool G4RichTrajectoryPoint::operator== ( const G4RichTrajectoryPoint right) const
inline

Definition at line 131 of file G4RichTrajectoryPoint.hh.

132{
133 return (this==&r);
134}

◆ operator==() [2/3]

G4bool G4TrajectoryPoint::operator== ( const G4TrajectoryPoint right) const
inlineinherited

Definition at line 98 of file G4TrajectoryPoint.hh.

99{
100 return (this==&right);
101}

◆ operator==() [3/3]

G4bool G4VTrajectoryPoint::operator== ( const G4VTrajectoryPoint right) const
inherited

Definition at line 44 of file G4VTrajectoryPoint.cc.

45{
46 return (this==&right);
47}

Field Documentation

◆ fpAuxiliaryPointVector

std::vector<G4ThreeVector>* G4RichTrajectoryPoint::fpAuxiliaryPointVector = nullptr
private

◆ fPosition

G4ThreeVector G4TrajectoryPoint::fPosition
privateinherited

◆ fPostStepPointGlobalTime

G4double G4RichTrajectoryPoint::fPostStepPointGlobalTime = 0.0
private

Definition at line 103 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fPostStepPointStatus

G4StepStatus G4RichTrajectoryPoint::fPostStepPointStatus = fUndefined
private

Definition at line 101 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fPostStepPointWeight

G4double G4RichTrajectoryPoint::fPostStepPointWeight = 0.0
private

Definition at line 107 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fpPostStepPointVolume

G4TouchableHandle G4RichTrajectoryPoint::fpPostStepPointVolume
private

Definition at line 105 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fpPreStepPointVolume

G4TouchableHandle G4RichTrajectoryPoint::fpPreStepPointVolume
private

Definition at line 104 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fpProcess

const G4VProcess* G4RichTrajectoryPoint::fpProcess = nullptr
private

Definition at line 99 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fPreStepPointGlobalTime

G4double G4RichTrajectoryPoint::fPreStepPointGlobalTime = 0.0
private

Definition at line 102 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fPreStepPointStatus

G4StepStatus G4RichTrajectoryPoint::fPreStepPointStatus = fUndefined
private

Definition at line 100 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fPreStepPointWeight

G4double G4RichTrajectoryPoint::fPreStepPointWeight = 0.0
private

Definition at line 106 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fRemainingEnergy

G4double G4RichTrajectoryPoint::fRemainingEnergy = 0.0
private

Definition at line 98 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().

◆ fTotEDep

G4double G4RichTrajectoryPoint::fTotEDep = 0.0
private

Definition at line 97 of file G4RichTrajectoryPoint.hh.

Referenced by CreateAttValues(), and G4RichTrajectoryPoint().


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