Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RichTrajectoryPoint Class Reference

#include <G4RichTrajectoryPoint.hh>

+ Inheritance diagram for G4RichTrajectoryPoint:

Public Member Functions

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

Detailed Description

Definition at line 67 of file G4RichTrajectoryPoint.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectoryPoint() [1/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( )
default

◆ G4RichTrajectoryPoint() [2/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Track * aTrack)

Definition at line 64 of file G4RichTrajectoryPoint.cc.

65 : G4TrajectoryPoint(aTrack->GetPosition()),
66 fPreStepPointGlobalTime(aTrack->GetGlobalTime()),
67 fPostStepPointGlobalTime(aTrack->GetGlobalTime()),
68 fpPreStepPointVolume(aTrack->GetTouchableHandle()),
69 fpPostStepPointVolume(aTrack->GetNextTouchableHandle()),
70 fPreStepPointWeight(aTrack->GetWeight()),
71 fPostStepPointWeight(aTrack->GetWeight())
72{}
const G4TouchableHandle & GetNextTouchableHandle() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
G4TrajectoryPoint()=default

◆ G4RichTrajectoryPoint() [3/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Step * aStep)

Definition at line 74 of file G4RichTrajectoryPoint.cc.

76 fpAuxiliaryPointVector(aStep->GetPointerToVectorOfAuxiliaryPoints()),
77 fTotEDep(aStep->GetTotalEnergyDeposit())
78{
79 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
80 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
81 if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) // First step
82 {
83 fRemainingEnergy = aStep->GetTrack()->GetKineticEnergy();
84 }
85 else {
86 fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep;
87 }
88 fpProcess = postStepPoint->GetProcessDefinedStep();
89 fPreStepPointStatus = preStepPoint->GetStepStatus();
90 fPostStepPointStatus = postStepPoint->GetStepStatus();
91 fPreStepPointGlobalTime = preStepPoint->GetGlobalTime();
92 fPostStepPointGlobalTime = postStepPoint->GetGlobalTime();
93 fpPreStepPointVolume = preStepPoint->GetTouchableHandle();
94 fpPostStepPointVolume = postStepPoint->GetTouchableHandle();
95 fPreStepPointWeight = preStepPoint->GetWeight();
96 fPostStepPointWeight = postStepPoint->GetWeight();
97}
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

◆ G4RichTrajectoryPoint() [4/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4RichTrajectoryPoint & right)
default

◆ ~G4RichTrajectoryPoint()

G4RichTrajectoryPoint::~G4RichTrajectoryPoint ( )
override

Definition at line 104 of file G4RichTrajectoryPoint.cc.

104{ delete fpAuxiliaryPointVector; }

Member Function Documentation

◆ CreateAttValues()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 193 of file G4RichTrajectoryPoint.cc.

194{
195 // Create base class att values...
196
197 std::vector<G4AttValue>* values = G4TrajectoryPoint::CreateAttValues();
198
199 if (fpAuxiliaryPointVector != nullptr) {
200 for (const auto& iAux : *fpAuxiliaryPointVector) {
201 values->push_back(G4AttValue("Aux", G4BestUnit(iAux, "Length"), ""));
202 }
203 }
204
205 values->push_back(G4AttValue("TED", G4BestUnit(fTotEDep, "Energy"), ""));
206 values->push_back(G4AttValue("RE", G4BestUnit(fRemainingEnergy, "Energy"), ""));
207
208 if (fpProcess != nullptr) {
209 values->push_back(G4AttValue("PDS", fpProcess->GetProcessName(), ""));
210 values->push_back(
211 G4AttValue("PTDS", G4VProcess::GetProcessTypeName(fpProcess->GetProcessType()), ""));
212 }
213 else {
214 values->push_back(G4AttValue("PDS", "None", ""));
215 values->push_back(G4AttValue("PTDS", "None", ""));
216 }
217
218 values->push_back(G4AttValue("PreStatus", Status(fPreStepPointStatus), ""));
219
220 values->push_back(G4AttValue("PostStatus", Status(fPostStepPointStatus), ""));
221
222 values->push_back(G4AttValue("PreT", G4BestUnit(fPreStepPointGlobalTime, "Time"), ""));
223
224 values->push_back(G4AttValue("PostT", G4BestUnit(fPostStepPointGlobalTime, "Time"), ""));
225
226 if (fpPreStepPointVolume && (fpPreStepPointVolume->GetVolume() != nullptr)) {
227 values->push_back(G4AttValue("PreVPath", Path(fpPreStepPointVolume), ""));
228 }
229 else {
230 values->push_back(G4AttValue("PreVPath", "None", ""));
231 }
232
233 if (fpPostStepPointVolume && (fpPostStepPointVolume->GetVolume() != nullptr)) {
234 values->push_back(G4AttValue("PostVPath", Path(fpPostStepPointVolume), ""));
235 }
236 else {
237 values->push_back(G4AttValue("PostVPath", "None", ""));
238 }
239
240 std::ostringstream oss1;
241 oss1 << fPreStepPointWeight;
242 values->push_back(G4AttValue("PreW", oss1.str(), ""));
243
244 std::ostringstream oss2;
245 oss2 << fPostStepPointWeight;
246 values->push_back(G4AttValue("PostW", oss2.str(), ""));
247
248#ifdef G4ATTDEBUG
249 G4cout << G4AttCheck(values, GetAttDefs());
250#endif
251
252 return values;
253}
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
const std::map< G4String, G4AttDef > * GetAttDefs() const override
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
std::vector< G4AttValue > * CreateAttValues() const override
static const G4String & GetProcessTypeName(G4ProcessType)
G4ProcessType GetProcessType() const
const G4String & GetProcessName() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 106 of file G4RichTrajectoryPoint.cc.

107{
108 G4bool isNew;
109 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectoryPoint", isNew);
110 if (isNew) {
111 // Copy base class att defs...
112 *store = *(G4TrajectoryPoint::GetAttDefs());
113
114 G4String ID;
115
116 ID = "Aux";
117 (*store)[ID] =
118 G4AttDef(ID, "Auxiliary Point Position", "Physics", "G4BestUnit", "G4ThreeVector");
119 ID = "TED";
120 (*store)[ID] = G4AttDef(ID, "Total Energy Deposit", "Physics", "G4BestUnit", "G4double");
121 ID = "RE";
122 (*store)[ID] = G4AttDef(ID, "Remaining Energy", "Physics", "G4BestUnit", "G4double");
123 ID = "PDS";
124 (*store)[ID] = G4AttDef(ID, "Process Defined Step", "Physics", "", "G4String");
125 ID = "PTDS";
126 (*store)[ID] = G4AttDef(ID, "Process Type Defined Step", "Physics", "", "G4String");
127 ID = "PreStatus";
128 (*store)[ID] = G4AttDef(ID, "Pre-step-point status", "Physics", "", "G4String");
129 ID = "PostStatus";
130 (*store)[ID] = G4AttDef(ID, "Post-step-point status", "Physics", "", "G4String");
131 ID = "PreT";
132 (*store)[ID] = G4AttDef(ID, "Pre-step-point global time", "Physics", "G4BestUnit", "G4double");
133 ID = "PostT";
134 (*store)[ID] = G4AttDef(ID, "Post-step-point global time", "Physics", "G4BestUnit", "G4double");
135 ID = "PreVPath";
136 (*store)[ID] = G4AttDef(ID, "Pre-step Volume Path", "Physics", "", "G4String");
137 ID = "PostVPath";
138 (*store)[ID] = G4AttDef(ID, "Post-step Volume Path", "Physics", "", "G4String");
139 ID = "PreW";
140 (*store)[ID] = G4AttDef(ID, "Pre-step-point weight", "Physics", "", "G4double");
141 ID = "PostW";
142 (*store)[ID] = G4AttDef(ID, "Post-step-point weight", "Physics", "", "G4double");
143 }
144 return store;
145}
bool G4bool
Definition G4Types.hh:86
const std::map< G4String, G4AttDef > * GetAttDefs() const override
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

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

◆ GetAuxiliaryPoints()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 128 of file G4RichTrajectoryPoint.hh.

129{
130 return fpAuxiliaryPointVector;
131}

◆ operator delete()

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

Definition at line 118 of file G4RichTrajectoryPoint.hh.

119{
120 aRichTrajectoryPointAllocator()->FreeSingle((G4RichTrajectoryPoint*)aRichTrajectoryPoint);
121}
G4TRACKING_DLL G4Allocator< G4RichTrajectoryPoint > *& aRichTrajectoryPointAllocator()

◆ operator new()

void * G4RichTrajectoryPoint::operator new ( size_t )
inline

Definition at line 110 of file G4RichTrajectoryPoint.hh.

111{
112 if (aRichTrajectoryPointAllocator() == nullptr) {
114 }
115 return (void*)aRichTrajectoryPointAllocator()->MallocSingle();
116}

◆ operator=()

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

◆ operator==()

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

Definition at line 123 of file G4RichTrajectoryPoint.hh.

124{
125 return (this == &r);
126}

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