Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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: Katsuya.Amako@kek.jp)
32// Makoto Asai (e-mail: asai@slac.stanford.edu)
33// Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
34// and on the extended code to:
35// John Allison (e-mail: John.Allison@manchester.ac.uk)
36// Joseph Perl (e-mail: perl@slac.stanford.edu)
37// --------------------------------------------------------------------
38
40
41#include "G4Track.hh"
42#include "G4Step.hh"
43#include "G4VProcess.hh"
44
45#include "G4AttDefStore.hh"
46#include "G4AttDef.hh"
47#include "G4AttValue.hh"
48#include "G4UnitsTable.hh"
49
50//#define G4ATTDEBUG
51#ifdef G4ATTDEBUG
52#include "G4AttCheck.hh"
53#endif
54
55#include <sstream>
56
58{
60 return _instance;
61}
62
64{
65}
66
68 : G4TrajectoryPoint(aTrack->GetPosition()),
69 fPreStepPointGlobalTime(aTrack->GetGlobalTime()),
70 fPostStepPointGlobalTime(aTrack->GetGlobalTime()),
71 fpPreStepPointVolume(aTrack->GetTouchableHandle()),
72 fpPostStepPointVolume(aTrack->GetNextTouchableHandle()),
73 fPreStepPointWeight(aTrack->GetWeight()),
74 fPostStepPointWeight(aTrack->GetWeight())
75{
76}
77
79 : G4TrajectoryPoint(aStep->GetPostStepPoint()->GetPosition()),
80 fpAuxiliaryPointVector(aStep->GetPointerToVectorOfAuxiliaryPoints()),
81 fTotEDep(aStep->GetTotalEnergyDeposit())
82{
83 G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
84 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
85 if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) // First step
86 {
87 fRemainingEnergy = aStep->GetTrack()->GetKineticEnergy();
88 }
89 else
90 {
91 fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep;
92 }
93 fpProcess = postStepPoint->GetProcessDefinedStep();
94 fPreStepPointStatus = preStepPoint->GetStepStatus();
95 fPostStepPointStatus = postStepPoint->GetStepStatus();
96 fPreStepPointGlobalTime = preStepPoint->GetGlobalTime();
97 fPostStepPointGlobalTime = postStepPoint->GetGlobalTime();
98 fpPreStepPointVolume = preStepPoint->GetTouchableHandle();
99 fpPostStepPointVolume = postStepPoint->GetTouchableHandle();
100 fPreStepPointWeight = preStepPoint->GetWeight();
101 fPostStepPointWeight = postStepPoint->GetWeight();
102}
103
106 fpAuxiliaryPointVector(r.fpAuxiliaryPointVector),
107 fTotEDep(r.fTotEDep),
108 fRemainingEnergy(r.fRemainingEnergy),
109 fpProcess(r.fpProcess),
110 fPreStepPointStatus(r.fPreStepPointStatus),
111 fPostStepPointStatus(r.fPostStepPointStatus),
112 fPreStepPointGlobalTime(r.fPreStepPointGlobalTime),
113 fPostStepPointGlobalTime(r.fPostStepPointGlobalTime),
114 fpPreStepPointVolume(r.fpPreStepPointVolume),
115 fpPostStepPointVolume(r.fpPostStepPointVolume),
116 fPreStepPointWeight(r.fPreStepPointWeight),
117 fPostStepPointWeight(r.fPostStepPointWeight)
118{
119}
120
122{
123 delete fpAuxiliaryPointVector;
124}
125
126const std::map<G4String,G4AttDef>* G4RichTrajectoryPoint::GetAttDefs() const
127{
128 G4bool isNew;
129 std::map<G4String,G4AttDef>* store
130 = G4AttDefStore::GetInstance("G4RichTrajectoryPoint",isNew);
131 if (isNew)
132 {
133 // Copy base class att defs...
134 *store = *(G4TrajectoryPoint::GetAttDefs());
135
136 G4String ID;
137
138 ID = "Aux";
139 (*store)[ID] = G4AttDef(ID, "Auxiliary Point Position",
140 "Physics","G4BestUnit","G4ThreeVector");
141 ID = "TED";
142 (*store)[ID] = G4AttDef(ID,"Total Energy Deposit",
143 "Physics","G4BestUnit","G4double");
144 ID = "RE";
145 (*store)[ID] = G4AttDef(ID,"Remaining Energy",
146 "Physics","G4BestUnit","G4double");
147 ID = "PDS";
148 (*store)[ID] = G4AttDef(ID,"Process Defined Step",
149 "Physics","","G4String");
150 ID = "PTDS";
151 (*store)[ID] = G4AttDef(ID,"Process Type Defined Step",
152 "Physics","","G4String");
153 ID = "PreStatus";
154 (*store)[ID] = G4AttDef(ID,"Pre-step-point status",
155 "Physics","","G4String");
156 ID = "PostStatus";
157 (*store)[ID] = G4AttDef(ID,"Post-step-point status",
158 "Physics","","G4String");
159 ID = "PreT";
160 (*store)[ID] = G4AttDef(ID,"Pre-step-point global time",
161 "Physics","G4BestUnit","G4double");
162 ID = "PostT";
163 (*store)[ID] = G4AttDef(ID,"Post-step-point global time",
164 "Physics","G4BestUnit","G4double");
165 ID = "PreVPath";
166 (*store)[ID] = G4AttDef(ID,"Pre-step Volume Path",
167 "Physics","","G4String");
168 ID = "PostVPath";
169 (*store)[ID] = G4AttDef(ID,"Post-step Volume Path",
170 "Physics","","G4String");
171 ID = "PreW";
172 (*store)[ID] = G4AttDef(ID,"Pre-step-point weight",
173 "Physics","","G4double");
174 ID = "PostW";
175 (*store)[ID] = G4AttDef(ID,"Post-step-point weight",
176 "Physics","","G4double");
177 }
178 return store;
179}
180
181static G4String Status(G4StepStatus stps)
182{
183 G4String status;
184 switch (stps)
185 {
186 case fWorldBoundary: status = "fWorldBoundary"; break;
187 case fGeomBoundary: status = "fGeomBoundary"; break;
188 case fAtRestDoItProc: status = "fAtRestDoItProc"; break;
189 case fAlongStepDoItProc: status = "fAlongStepDoItProc"; break;
190 case fPostStepDoItProc: status = "fPostStepDoItProc"; break;
191 case fUserDefinedLimit: status = "fUserDefinedLimit"; break;
192 case fExclusivelyForcedProc: status = "fExclusivelyForcedProc"; break;
193 case fUndefined: status = "fUndefined"; break;
194 default: status = "Not recognised"; break;
195 }
196 return status;
197}
198
199static G4String Path(const G4TouchableHandle& th)
200{
201 std::ostringstream oss;
202 G4int depth = th->GetHistoryDepth();
203 for (G4int i = depth; i >= 0; --i)
204 {
205 oss << th->GetVolume(i)->GetName()
206 << ':' << th->GetCopyNumber(i);
207 if (i != 0) oss << '/';
208 }
209 return oss.str();
210}
211
212std::vector<G4AttValue>* G4RichTrajectoryPoint::CreateAttValues() const
213{
214 // Create base class att values...
215
216 std::vector<G4AttValue>* values = G4TrajectoryPoint::CreateAttValues();
217
218 if (fpAuxiliaryPointVector != nullptr)
219 {
220 for (auto iAux = fpAuxiliaryPointVector->cbegin();
221 iAux != fpAuxiliaryPointVector->cend(); ++iAux)
222 {
223 values->push_back(G4AttValue("Aux",G4BestUnit(*iAux,"Length"),""));
224 }
225 }
226
227 values->push_back(G4AttValue("TED",G4BestUnit(fTotEDep,"Energy"),""));
228 values->push_back(G4AttValue("RE",G4BestUnit(fRemainingEnergy,"Energy"),""));
229
230 if (fpProcess != nullptr)
231 {
232 values->push_back
233 (G4AttValue("PDS",fpProcess->GetProcessName(),""));
234 values->push_back
236 ("PTDS",G4VProcess::GetProcessTypeName(fpProcess->GetProcessType()),""));
237 }
238 else
239 {
240 values->push_back(G4AttValue("PDS","None",""));
241 values->push_back(G4AttValue("PTDS","None",""));
242 }
243
244 values->push_back
245 (G4AttValue("PreStatus",Status(fPreStepPointStatus),""));
246
247 values->push_back
248 (G4AttValue("PostStatus",Status(fPostStepPointStatus),""));
249
250 values->push_back
251 (G4AttValue("PreT",G4BestUnit(fPreStepPointGlobalTime,"Time"),""));
252
253 values->push_back
254 (G4AttValue("PostT",G4BestUnit(fPostStepPointGlobalTime,"Time"),""));
255
256 if (fpPreStepPointVolume && fpPreStepPointVolume->GetVolume())
257 {
258 values->push_back(G4AttValue("PreVPath",Path(fpPreStepPointVolume),""));
259 }
260 else
261 {
262 values->push_back(G4AttValue("PreVPath","None",""));
263 }
264
265 if (fpPostStepPointVolume && fpPostStepPointVolume->GetVolume())
266 {
267 values->push_back(G4AttValue("PostVPath",Path(fpPostStepPointVolume),""));
268 }
269 else
270 {
271 values->push_back(G4AttValue("PostVPath","None",""));
272 }
273
274 std::ostringstream oss1;
275 oss1 << fPreStepPointWeight;
276 values->push_back(G4AttValue("PreW",oss1.str(),""));
277
278 std::ostringstream oss2;
279 oss2 << fPostStepPointWeight;
280 values->push_back(G4AttValue("PostW",oss2.str(),""));
281
282#ifdef G4ATTDEBUG
283 G4cout << G4AttCheck(values,GetAttDefs());
284#endif
285
286 return values;
287}
G4Allocator< G4RichTrajectoryPoint > *& aRichTrajectoryPointAllocator()
G4StepStatus
Definition: G4StepStatus.hh:40
@ fGeomBoundary
Definition: G4StepStatus.hh:43
@ fWorldBoundary
Definition: G4StepStatus.hh:41
@ fUserDefinedLimit
Definition: G4StepStatus.hh:51
@ fUndefined
Definition: G4StepStatus.hh:55
@ fPostStepDoItProc
Definition: G4StepStatus.hh:49
@ fAtRestDoItProc
Definition: G4StepStatus.hh:45
@ fAlongStepDoItProc
Definition: G4StepStatus.hh:47
@ fExclusivelyForcedProc
Definition: G4StepStatus.hh:53
#define G4BestUnit(a, b)
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
G4StepStatus GetStepStatus() const
G4double GetGlobalTime() const
const G4VProcess * GetProcessDefinedStep() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4double GetWeight() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4int GetCurrentStepNumber() const
G4double GetKineticEnergy() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
const G4String & GetName() const
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:134
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:392
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386
G4int GetCopyNumber(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
virtual G4int GetHistoryDepth() const
Definition: G4VTouchable.cc:74
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
#define G4ThreadLocalStatic
Definition: tls.hh:76