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

#include <G4ClonedRichTrajectoryPoint.hh>

+ Inheritance diagram for G4ClonedRichTrajectoryPoint:

Public Member Functions

 G4ClonedRichTrajectoryPoint ()=default
 
 G4ClonedRichTrajectoryPoint (const G4Track *)
 
 G4ClonedRichTrajectoryPoint (const G4Step *)
 
 G4ClonedRichTrajectoryPoint (const G4RichTrajectoryPoint &right)
 
 ~G4ClonedRichTrajectoryPoint () override
 
G4ClonedRichTrajectoryPointoperator= (const G4ClonedRichTrajectoryPoint &)=delete
 
G4bool operator== (const G4ClonedRichTrajectoryPoint &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
 

Detailed Description

Definition at line 57 of file G4ClonedRichTrajectoryPoint.hh.

Constructor & Destructor Documentation

◆ G4ClonedRichTrajectoryPoint() [1/4]

G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint ( )
default

◆ G4ClonedRichTrajectoryPoint() [2/4]

G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint ( const G4Track * aTrack)

Definition at line 56 of file G4ClonedRichTrajectoryPoint.cc.

57 : fPosition(aTrack->GetPosition()),
58 fPreStepPointGlobalTime(aTrack->GetGlobalTime()),
59 fPostStepPointGlobalTime(aTrack->GetGlobalTime()),
60 fpPreStepPointVolume(aTrack->GetTouchableHandle()),
61 fpPostStepPointVolume(aTrack->GetNextTouchableHandle()),
62 fPreStepPointWeight(aTrack->GetWeight()),
63 fPostStepPointWeight(aTrack->GetWeight())
64{}
const G4TouchableHandle & GetNextTouchableHandle() const
G4double GetWeight() const
const G4ThreeVector & GetPosition() const
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const

◆ G4ClonedRichTrajectoryPoint() [3/4]

G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint ( const G4Step * aStep)

Definition at line 66 of file G4ClonedRichTrajectoryPoint.cc.

67 : fPosition(aStep->GetPostStepPoint()->GetPosition()),
68 fpAuxiliaryPointVector(aStep->GetPointerToVectorOfAuxiliaryPoints()),
69 fTotEDep(aStep->GetTotalEnergyDeposit())
70{
71 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
72 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
73 if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) // First step
74 {
75 fRemainingEnergy = aStep->GetTrack()->GetKineticEnergy();
76 }
77 else {
78 fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep;
79 }
80 fpProcess = postStepPoint->GetProcessDefinedStep();
81 fPreStepPointStatus = preStepPoint->GetStepStatus();
82 fPostStepPointStatus = postStepPoint->GetStepStatus();
83 fPreStepPointGlobalTime = preStepPoint->GetGlobalTime();
84 fPostStepPointGlobalTime = postStepPoint->GetGlobalTime();
85 fpPreStepPointVolume = preStepPoint->GetTouchableHandle();
86 fpPostStepPointVolume = postStepPoint->GetTouchableHandle();
87 fPreStepPointWeight = preStepPoint->GetWeight();
88 fPostStepPointWeight = postStepPoint->GetWeight();
89}
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

◆ G4ClonedRichTrajectoryPoint() [4/4]

G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint ( const G4RichTrajectoryPoint & right)

Definition at line 91 of file G4ClonedRichTrajectoryPoint.cc.

92{
93 fPosition = right.fPosition;
94 fTotEDep = right.fTotEDep;
95 fRemainingEnergy = right.fRemainingEnergy;
96 fpProcess = right.fpProcess;
97 fPreStepPointStatus = right.fPreStepPointStatus;
98 fPostStepPointStatus = right.fPostStepPointStatus;
99 fPreStepPointGlobalTime = right.fPreStepPointGlobalTime;
100 fPostStepPointGlobalTime = right.fPostStepPointGlobalTime;
101 fpPreStepPointVolume = right.fpPreStepPointVolume;
102 fpPostStepPointVolume = right.fpPostStepPointVolume;
103 fPreStepPointWeight = right.fPreStepPointWeight;
104 fPostStepPointWeight = right.fPostStepPointWeight;
105 if(right.fpAuxiliaryPointVector!=nullptr && right.fpAuxiliaryPointVector->size()>0)
106 {
107 fpAuxiliaryPointVector = new std::vector<G4ThreeVector>;
108 for(auto& i : *(right.fpAuxiliaryPointVector))
109 { fpAuxiliaryPointVector->push_back(i); }
110 }
111}

◆ ~G4ClonedRichTrajectoryPoint()

G4ClonedRichTrajectoryPoint::~G4ClonedRichTrajectoryPoint ( )
override

Definition at line 113 of file G4ClonedRichTrajectoryPoint.cc.

113{ delete fpAuxiliaryPointVector; }

Member Function Documentation

◆ CreateAttValues()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 201 of file G4ClonedRichTrajectoryPoint.cc.

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 115 of file G4ClonedRichTrajectoryPoint.cc.

116{
117 G4bool isNew;
118 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectoryPoint", isNew);
119 if (isNew) {
120 G4String ID;
121
122 ID = "Pos";
123 (*store)[ID] = G4AttDef(ID, "Position", "Physics", "G4BestUnit", "G4ThreeVector");
124 ID = "Aux";
125 (*store)[ID] =
126 G4AttDef(ID, "Auxiliary Point Position", "Physics", "G4BestUnit", "G4ThreeVector");
127 ID = "TED";
128 (*store)[ID] = G4AttDef(ID, "Total Energy Deposit", "Physics", "G4BestUnit", "G4double");
129 ID = "RE";
130 (*store)[ID] = G4AttDef(ID, "Remaining Energy", "Physics", "G4BestUnit", "G4double");
131 ID = "PDS";
132 (*store)[ID] = G4AttDef(ID, "Process Defined Step", "Physics", "", "G4String");
133 ID = "PTDS";
134 (*store)[ID] = G4AttDef(ID, "Process Type Defined Step", "Physics", "", "G4String");
135 ID = "PreStatus";
136 (*store)[ID] = G4AttDef(ID, "Pre-step-point status", "Physics", "", "G4String");
137 ID = "PostStatus";
138 (*store)[ID] = G4AttDef(ID, "Post-step-point status", "Physics", "", "G4String");
139 ID = "PreT";
140 (*store)[ID] = G4AttDef(ID, "Pre-step-point global time", "Physics", "G4BestUnit", "G4double");
141 ID = "PostT";
142 (*store)[ID] = G4AttDef(ID, "Post-step-point global time", "Physics", "G4BestUnit", "G4double");
143 ID = "PreVPath";
144 (*store)[ID] = G4AttDef(ID, "Pre-step Volume Path", "Physics", "", "G4String");
145 ID = "PostVPath";
146 (*store)[ID] = G4AttDef(ID, "Post-step Volume Path", "Physics", "", "G4String");
147 ID = "PreW";
148 (*store)[ID] = G4AttDef(ID, "Pre-step-point weight", "Physics", "", "G4double");
149 ID = "PostW";
150 (*store)[ID] = G4AttDef(ID, "Post-step-point weight", "Physics", "", "G4double");
151 }
152 return store;
153}
bool G4bool
Definition G4Types.hh:86
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)

Referenced by CreateAttValues().

◆ GetAuxiliaryPoints()

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

Reimplemented from G4VTrajectoryPoint.

Definition at line 123 of file G4ClonedRichTrajectoryPoint.hh.

124{
125 return fpAuxiliaryPointVector;
126}

◆ GetPosition()

const G4ThreeVector G4ClonedRichTrajectoryPoint::GetPosition ( ) const
inlineoverridevirtual

Implements G4VTrajectoryPoint.

Definition at line 78 of file G4ClonedRichTrajectoryPoint.hh.

78{ return fPosition; }

Referenced by G4ClonedRichTrajectoryPoint(), and G4ClonedRichTrajectoryPoint().

◆ operator delete()

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

Definition at line 113 of file G4ClonedRichTrajectoryPoint.hh.

114{
115 aClonedRichTrajectoryPointAllocator()->FreeSingle((G4ClonedRichTrajectoryPoint*)aRichTrajectoryPoint);
116}
G4TRACKING_DLL G4Allocator< G4ClonedRichTrajectoryPoint > *& aClonedRichTrajectoryPointAllocator()

◆ operator new()

void * G4ClonedRichTrajectoryPoint::operator new ( size_t )
inline

Definition at line 105 of file G4ClonedRichTrajectoryPoint.hh.

106{
107 if (aClonedRichTrajectoryPointAllocator() == nullptr) {
108 aClonedRichTrajectoryPointAllocator() = new G4Allocator<G4ClonedRichTrajectoryPoint>;
109 }
110 return (void*)aClonedRichTrajectoryPointAllocator()->MallocSingle();
111}

◆ operator=()

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

◆ operator==()

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

Definition at line 118 of file G4ClonedRichTrajectoryPoint.hh.

119{
120 return (this == &r);
121}

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