Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RichTrajectoryPoint.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4RichTrajectoryPoint class implementation
27//
28// Contact:
29// Questions and comments on G4TrajectoryPoint, on which this is based,
30// should be sent to
31// Katsuya Amako (e-mail: [email protected])
32// Makoto Asai (e-mail: [email protected])
33// Takashi Sasaki (e-mail: [email protected])
34// and on the extended code to:
35// John Allison (e-mail: [email protected])
36// Joseph Perl (e-mail: [email protected])
37// --------------------------------------------------------------------
38
40
41#include "G4AttDef.hh"
42#include "G4AttDefStore.hh"
43#include "G4AttValue.hh"
44#include "G4Step.hh"
45#include "G4Track.hh"
46#include "G4UIcommand.hh"
47#include "G4UnitsTable.hh"
48#include "G4VProcess.hh"
49
50// #define G4ATTDEBUG
51#ifdef G4ATTDEBUG
52# include "G4AttCheck.hh"
53#endif
54
55#include <sstream>
56
62
64
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{}
74
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}
99
100// This is certainly wrong. the fpAuxiliaryPointVector is owned by this class,
101// and must be deep copied. Note also the copy-assignment operator is deleted.
102// However, G4RichTrajectory requires a copy constructor....
104
105G4RichTrajectoryPoint::~G4RichTrajectoryPoint() { delete fpAuxiliaryPointVector; }
106
107const std::map<G4String, G4AttDef>* G4RichTrajectoryPoint::GetAttDefs() const
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}
146
147static G4String Status(G4StepStatus stps)
148{
149 G4String status;
150 switch (stps) {
151 case fWorldBoundary:
152 status = "fWorldBoundary";
153 break;
154 case fGeomBoundary:
155 status = "fGeomBoundary";
156 break;
157 case fAtRestDoItProc:
158 status = "fAtRestDoItProc";
159 break;
161 status = "fAlongStepDoItProc";
162 break;
164 status = "fPostStepDoItProc";
165 break;
167 status = "fUserDefinedLimit";
168 break;
170 status = "fExclusivelyForcedProc";
171 break;
172 case fUndefined:
173 status = "fUndefined";
174 break;
175 default:
176 status = "Not recognised";
177 break;
178 }
179 return status;
180}
181
182static G4String Path(const G4TouchableHandle& th)
183{
184 std::ostringstream oss;
185 G4int depth = th->GetHistoryDepth();
186 for (G4int i = depth; i >= 0; --i) {
187 oss << th->GetVolume(i)->GetName() << ':' << th->GetCopyNumber(i);
188 if (i != 0) oss << '/';
189 }
190 return oss.str();
191}
192
193std::vector<G4AttValue>* G4RichTrajectoryPoint::CreateAttValues() const
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}
G4Allocator< G4RichTrajectoryPoint > *& aRichTrajectoryPointAllocator()
G4StepStatus
@ fGeomBoundary
@ fWorldBoundary
@ fUserDefinedLimit
@ fUndefined
@ fPostStepDoItProc
@ fAtRestDoItProc
@ fAlongStepDoItProc
@ fExclusivelyForcedProc
#define G4BestUnit(a, b)
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
const std::map< G4String, G4AttDef > * GetAttDefs() const override
const G4ThreeVector GetPosition() const override
std::vector< G4AttValue > * CreateAttValues() const override
G4StepStatus GetStepStatus() const
G4double GetGlobalTime() const
const G4VProcess * GetProcessDefinedStep() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetWeight() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetCopyNumber(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
virtual G4int GetHistoryDepth() const
G4int GetCurrentStepNumber() const
const G4String & GetName() const
static const G4String & GetProcessTypeName(G4ProcessType)
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
#define G4ThreadLocalStatic
Definition tls.hh:76