Geant4 11.3.0
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 G4ThreeVector GetPosition () const override
 
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 G4VTrajectoryPoint
 G4VTrajectoryPoint ()=default
 
virtual ~G4VTrajectoryPoint ()=default
 
G4bool operator== (const G4VTrajectoryPoint &right) const
 

Friends

class G4ClonedRichTrajectoryPoint
 

Detailed Description

Definition at line 71 of file G4RichTrajectoryPoint.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectoryPoint() [1/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( )
default

◆ G4RichTrajectoryPoint() [2/4]

G4RichTrajectoryPoint::G4RichTrajectoryPoint ( const G4Track * aTrack)

Definition at line 65 of file G4RichTrajectoryPoint.cc.

66 : fPosition(aTrack->GetPosition()),
67 fPreStepPointGlobalTime(aTrack->GetGlobalTime()),
68 fPostStepPointGlobalTime(aTrack->GetGlobalTime()),
69 fpPreStepPointVolume(aTrack->GetTouchableHandle()),
70 fpPostStepPointVolume(aTrack->GetNextTouchableHandle()),
71 fPreStepPointWeight(aTrack->GetWeight()),
72 fPostStepPointWeight(aTrack->GetWeight())
73{}
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 75 of file G4RichTrajectoryPoint.cc.

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

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

◆ GetAttDefs()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 107 of file G4RichTrajectoryPoint.cc.

108{
109 G4bool isNew;
110 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectoryPoint", isNew);
111 if (isNew) {
112 G4String ID;
113
114 ID = "Pos";
115 (*store)[ID] = G4AttDef(ID, "Position", "Physics", "G4BestUnit", "G4ThreeVector");
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
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 140 of file G4RichTrajectoryPoint.hh.

141{
142 return fpAuxiliaryPointVector;
143}

◆ GetPosition()

const G4ThreeVector G4RichTrajectoryPoint::GetPosition ( ) const
inlineoverridevirtual

Implements G4VTrajectoryPoint.

Definition at line 95 of file G4RichTrajectoryPoint.hh.

95{ return fPosition; }

Referenced by G4RichTrajectoryPoint(), and G4RichTrajectoryPoint().

◆ operator delete()

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

Definition at line 130 of file G4RichTrajectoryPoint.hh.

131{
132 aRichTrajectoryPointAllocator()->FreeSingle((G4RichTrajectoryPoint*)aRichTrajectoryPoint);
133}
G4TRACKING_DLL G4Allocator< G4RichTrajectoryPoint > *& aRichTrajectoryPointAllocator()

◆ operator new()

void * G4RichTrajectoryPoint::operator new ( size_t )
inline

Definition at line 122 of file G4RichTrajectoryPoint.hh.

123{
124 if (aRichTrajectoryPointAllocator() == nullptr) {
125 aRichTrajectoryPointAllocator() = new G4Allocator<G4RichTrajectoryPoint>;
126 }
127 return (void*)aRichTrajectoryPointAllocator()->MallocSingle();
128}

◆ operator=()

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

◆ operator==()

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

Definition at line 135 of file G4RichTrajectoryPoint.hh.

136{
137 return (this == &r);
138}

Friends And Related Symbol Documentation

◆ G4ClonedRichTrajectoryPoint

friend class G4ClonedRichTrajectoryPoint
friend

Definition at line 74 of file G4RichTrajectoryPoint.hh.

Referenced by G4ClonedRichTrajectoryPoint.


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