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

#include <G4RichTrajectory.hh>

+ Inheritance diagram for G4RichTrajectory:

Public Member Functions

 G4RichTrajectory ()=default
 
 G4RichTrajectory (const G4Track *aTrack)
 
 ~G4RichTrajectory () override
 
 G4RichTrajectory (G4RichTrajectory &)
 
G4RichTrajectoryoperator= (const G4RichTrajectory &)=delete
 
G4bool operator== (const G4RichTrajectory &r) const
 
void * operator new (size_t)
 
void operator delete (void *)
 
void ShowTrajectory (std::ostream &os=G4cout) const override
 
void DrawTrajectory () const override
 
void AppendStep (const G4Step *aStep) override
 
void MergeTrajectory (G4VTrajectory *secondTrajectory) override
 
G4int GetPointEntries () const override
 
G4VTrajectoryPointGetPoint (G4int i) const override
 
const std::map< G4String, G4AttDef > * GetAttDefs () const override
 
std::vector< G4AttValue > * CreateAttValues () const override
 
- Public Member Functions inherited from G4Trajectory
 G4Trajectory ()=default
 
 G4Trajectory (const G4Track *aTrack)
 
 G4Trajectory (G4Trajectory &)
 
 ~G4Trajectory () override
 
void * operator new (size_t)
 
void operator delete (void *)
 
G4bool operator== (const G4Trajectory &r) const
 
G4int GetTrackID () const override
 
G4int GetParentID () const override
 
G4String GetParticleName () const override
 
G4double GetCharge () const override
 
G4int GetPDGEncoding () const override
 
G4double GetInitialKineticEnergy () const
 
G4ThreeVector GetInitialMomentum () const override
 
void ShowTrajectory (std::ostream &os=G4cout) const override
 
void DrawTrajectory () const override
 
void AppendStep (const G4Step *aStep) override
 
G4int GetPointEntries () const override
 
G4VTrajectoryPointGetPoint (G4int i) const override
 
void MergeTrajectory (G4VTrajectory *secondTrajectory) override
 
G4ParticleDefinitionGetParticleDefinition ()
 
const std::map< G4String, G4AttDef > * GetAttDefs () const override
 
std::vector< G4AttValue > * CreateAttValues () const override
 
- Public Member Functions inherited from G4VTrajectory
 G4VTrajectory ()=default
 
virtual ~G4VTrajectory ()=default
 
G4bool operator== (const G4VTrajectory &right) const
 

Detailed Description

Definition at line 58 of file G4RichTrajectory.hh.

Constructor & Destructor Documentation

◆ G4RichTrajectory() [1/3]

G4RichTrajectory::G4RichTrajectory ( )
default

◆ G4RichTrajectory() [2/3]

G4RichTrajectory::G4RichTrajectory ( const G4Track * aTrack)

Definition at line 63 of file G4RichTrajectory.cc.

64 : G4Trajectory(aTrack) // Note: this initialises the base class data
65 // members and, unfortunately but never mind,
66 // creates a G4TrajectoryPoint in
67 // TrajectoryPointContainer that we cannot
68 // access because it's private. We store the
69 // same information (plus more) in a
70 // G4RichTrajectoryPoint in the
71 // RichTrajectoryPointsContainer
72{
73 fpInitialVolume = aTrack->GetTouchableHandle();
74 fpInitialNextVolume = aTrack->GetNextTouchableHandle();
75 fpCreatorProcess = aTrack->GetCreatorProcess();
76 fCreatorModelID = aTrack->GetCreatorModelID();
77
78 // On construction, set final values to initial values.
79 // Final values are updated at the addition of every step - see AppendStep.
80 //
81 fpFinalVolume = aTrack->GetTouchableHandle();
82 fpFinalNextVolume = aTrack->GetNextTouchableHandle();
83 fpEndingProcess = aTrack->GetCreatorProcess();
84 fFinalKineticEnergy = aTrack->GetKineticEnergy();
85
86 // Insert the first rich trajectory point (see note above)...
87 //
88 fpRichPointsContainer = new RichTrajectoryPointsContainer;
89 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
90}
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
G4int GetCreatorModelID() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4Trajectory()=default

◆ ~G4RichTrajectory()

G4RichTrajectory::~G4RichTrajectory ( )
override

Definition at line 109 of file G4RichTrajectory.cc.

110{
111 if (fpRichPointsContainer != nullptr) {
112 for (auto& i : *fpRichPointsContainer) {
113 delete i;
114 }
115 fpRichPointsContainer->clear();
116 delete fpRichPointsContainer;
117 }
118}

◆ G4RichTrajectory() [3/3]

G4RichTrajectory::G4RichTrajectory ( G4RichTrajectory & right)

Definition at line 92 of file G4RichTrajectory.cc.

92 : G4Trajectory(right)
93{
94 fpInitialVolume = right.fpInitialVolume;
95 fpInitialNextVolume = right.fpInitialNextVolume;
96 fpCreatorProcess = right.fpCreatorProcess;
97 fCreatorModelID = right.fCreatorModelID;
98 fpFinalVolume = right.fpFinalVolume;
99 fpFinalNextVolume = right.fpFinalNextVolume;
100 fpEndingProcess = right.fpEndingProcess;
101 fFinalKineticEnergy = right.fFinalKineticEnergy;
102 fpRichPointsContainer = new RichTrajectoryPointsContainer;
103 for (auto& i : *right.fpRichPointsContainer) {
104 auto rightPoint = (G4RichTrajectoryPoint*)i;
105 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
106 }
107}

Member Function Documentation

◆ AppendStep()

void G4RichTrajectory::AppendStep ( const G4Step * aStep)
overridevirtual

Implements G4VTrajectory.

Definition at line 120 of file G4RichTrajectory.cc.

121{
122 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
123
124 // Except for first step, which is a sort of virtual step to start
125 // the track, compute the final values...
126 //
127 const G4Track* track = aStep->GetTrack();
128 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
129 if (track->GetCurrentStepNumber() > 0) {
130 fpFinalVolume = track->GetTouchableHandle();
131 fpFinalNextVolume = track->GetNextTouchableHandle();
132 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
133 fFinalKineticEnergy =
135 }
136}
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4int GetCurrentStepNumber() const

◆ CreateAttValues()

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

Reimplemented from G4VTrajectory.

Definition at line 230 of file G4RichTrajectory.cc.

231{
232 // Create base class att values...
233 std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
234
235 if (fpInitialVolume && (fpInitialVolume->GetVolume() != nullptr)) {
236 values->push_back(G4AttValue("IVPath", Path(fpInitialVolume), ""));
237 }
238 else {
239 values->push_back(G4AttValue("IVPath", "None", ""));
240 }
241
242 if (fpInitialNextVolume && (fpInitialNextVolume->GetVolume() != nullptr)) {
243 values->push_back(G4AttValue("INVPath", Path(fpInitialNextVolume), ""));
244 }
245 else {
246 values->push_back(G4AttValue("INVPath", "None", ""));
247 }
248
249 if (fpCreatorProcess != nullptr) {
250 values->push_back(G4AttValue("CPN", fpCreatorProcess->GetProcessName(), ""));
251 G4ProcessType type = fpCreatorProcess->GetProcessType();
252 values->push_back(G4AttValue("CPTN", G4VProcess::GetProcessTypeName(type), ""));
253 values->push_back(G4AttValue("CMID", G4UIcommand::ConvertToString(fCreatorModelID), ""));
254 const G4String& creatorModelName = G4PhysicsModelCatalog::GetModelNameFromID(fCreatorModelID);
255 values->push_back(G4AttValue("CMN", creatorModelName, ""));
256 }
257 else {
258 values->push_back(G4AttValue("CPN", "None", ""));
259 values->push_back(G4AttValue("CPTN", "None", ""));
260 values->push_back(G4AttValue("CMID", "None", ""));
261 values->push_back(G4AttValue("CMN", "None", ""));
262 }
263
264 if (fpFinalVolume && (fpFinalVolume->GetVolume() != nullptr)) {
265 values->push_back(G4AttValue("FVPath", Path(fpFinalVolume), ""));
266 }
267 else {
268 values->push_back(G4AttValue("FVPath", "None", ""));
269 }
270
271 if (fpFinalNextVolume && (fpFinalNextVolume->GetVolume() != nullptr)) {
272 values->push_back(G4AttValue("FNVPath", Path(fpFinalNextVolume), ""));
273 }
274 else {
275 values->push_back(G4AttValue("FNVPath", "None", ""));
276 }
277
278 if (fpEndingProcess != nullptr) {
279 values->push_back(G4AttValue("EPN", fpEndingProcess->GetProcessName(), ""));
280 G4ProcessType type = fpEndingProcess->GetProcessType();
281 values->push_back(G4AttValue("EPTN", G4VProcess::GetProcessTypeName(type), ""));
282 }
283 else {
284 values->push_back(G4AttValue("EPN", "None", ""));
285 values->push_back(G4AttValue("EPTN", "None", ""));
286 }
287
288 values->push_back(G4AttValue("FKE", G4BestUnit(fFinalKineticEnergy, "Energy"), ""));
289
290#ifdef G4ATTDEBUG
291 G4cout << G4AttCheck(values, GetAttDefs());
292#endif
293
294 return values;
295}
G4ProcessType
#define G4BestUnit(a, b)
G4GLOB_DLL std::ostream G4cout
static const G4String GetModelNameFromID(const G4int modelID)
const std::map< G4String, G4AttDef > * GetAttDefs() const override
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
std::vector< G4AttValue > * CreateAttValues() const override
static G4String ConvertToString(G4bool boolVal)
static const G4String & GetProcessTypeName(G4ProcessType)
G4ProcessType GetProcessType() const
const G4String & GetProcessName() const

◆ DrawTrajectory()

void G4RichTrajectory::DrawTrajectory ( ) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 162 of file G4RichTrajectory.cc.

163{
164 // Invoke the default implementation in G4VTrajectory...
165 //
167
168 // ... or override with your own code here.
169}
virtual void DrawTrajectory() const

◆ GetAttDefs()

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

Reimplemented from G4VTrajectory.

Definition at line 171 of file G4RichTrajectory.cc.

172{
173 G4bool isNew;
174 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectory", isNew);
175 if (isNew) {
176 // Get att defs from base class...
177 //
178 *store = *(G4Trajectory::GetAttDefs());
179
180 G4String ID;
181
182 ID = "IVPath";
183 (*store)[ID] = G4AttDef(ID, "Initial Volume Path", "Physics", "", "G4String");
184
185 ID = "INVPath";
186 (*store)[ID] = G4AttDef(ID, "Initial Next Volume Path", "Physics", "", "G4String");
187
188 ID = "CPN";
189 (*store)[ID] = G4AttDef(ID, "Creator Process Name", "Physics", "", "G4String");
190
191 ID = "CPTN";
192 (*store)[ID] = G4AttDef(ID, "Creator Process Type Name", "Physics", "", "G4String");
193
194 ID = "CMID";
195 (*store)[ID] = G4AttDef(ID, "Creator Model ID", "Physics", "", "G4int");
196
197 ID = "CMN";
198 (*store)[ID] = G4AttDef(ID, "Creator Model Name", "Physics", "", "G4String");
199
200 ID = "FVPath";
201 (*store)[ID] = G4AttDef(ID, "Final Volume Path", "Physics", "", "G4String");
202
203 ID = "FNVPath";
204 (*store)[ID] = G4AttDef(ID, "Final Next Volume Path", "Physics", "", "G4String");
205
206 ID = "EPN";
207 (*store)[ID] = G4AttDef(ID, "Ending Process Name", "Physics", "", "G4String");
208
209 ID = "EPTN";
210 (*store)[ID] = G4AttDef(ID, "Ending Process Type Name", "Physics", "", "G4String");
211
212 ID = "FKE";
213 (*store)[ID] = G4AttDef(ID, "Final kinetic energy", "Physics", "G4BestUnit", "G4double");
214 }
215
216 return store;
217}
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().

◆ GetPoint()

G4VTrajectoryPoint * G4RichTrajectory::GetPoint ( G4int i) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 126 of file G4RichTrajectory.hh.

127{
128 return (*fpRichPointsContainer)[i];
129}

Referenced by G4TrajectoryDrawByEncounteredVolume::Draw(), and G4TrajectoryEncounteredVolumeFilter::Evaluate().

◆ GetPointEntries()

G4int G4RichTrajectory::GetPointEntries ( ) const
inlineoverridevirtual

Implements G4VTrajectory.

Definition at line 121 of file G4RichTrajectory.hh.

122{
123 return G4int(fpRichPointsContainer->size());
124}
int G4int
Definition G4Types.hh:85

Referenced by G4TrajectoryDrawByEncounteredVolume::Draw(), and G4TrajectoryEncounteredVolumeFilter::Evaluate().

◆ MergeTrajectory()

void G4RichTrajectory::MergeTrajectory ( G4VTrajectory * secondTrajectory)
overridevirtual

Implements G4VTrajectory.

Definition at line 138 of file G4RichTrajectory.cc.

139{
140 if (secondTrajectory == nullptr) return;
141
142 auto seco = (G4RichTrajectory*)secondTrajectory;
143 G4int ent = seco->GetPointEntries();
144 for (G4int i = 1; i < ent; ++i) {
145 // initial point of the second trajectory should not be merged
146 //
147 fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
148 }
149 delete (*seco->fpRichPointsContainer)[0];
150 seco->fpRichPointsContainer->clear();
151}

◆ operator delete()

void G4RichTrajectory::operator delete ( void * aRichTrajectory)
inline

Definition at line 116 of file G4RichTrajectory.hh.

117{
118 aRichTrajectoryAllocator()->FreeSingle((G4RichTrajectory*)aRichTrajectory);
119}
G4TRACKING_DLL G4Allocator< G4RichTrajectory > *& aRichTrajectoryAllocator()

◆ operator new()

void * G4RichTrajectory::operator new ( size_t )
inline

Definition at line 108 of file G4RichTrajectory.hh.

109{
110 if (aRichTrajectoryAllocator() == nullptr) {
112 }
113 return (void*)aRichTrajectoryAllocator()->MallocSingle();
114}

◆ operator=()

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

◆ operator==()

G4bool G4RichTrajectory::operator== ( const G4RichTrajectory & r) const
inline

Definition at line 73 of file G4RichTrajectory.hh.

73{ return (this == &r); }

◆ ShowTrajectory()

void G4RichTrajectory::ShowTrajectory ( std::ostream & os = G4cout) const
overridevirtual

Reimplemented from G4VTrajectory.

Definition at line 153 of file G4RichTrajectory.cc.

154{
155 // Invoke the default implementation in G4VTrajectory...
156 //
158
159 // ... or override with your own code here.
160}
virtual void ShowTrajectory(std::ostream &os=G4cout) const

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