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
G4RichTrajectory.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// G4RichTrajectory class implementation
27//
28// Contact:
29// Questions and comments on G4Trajectory, 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
39#include "G4RichTrajectory.hh"
41#include "G4AttDefStore.hh"
42#include "G4AttDef.hh"
43#include "G4AttValue.hh"
44#include "G4UIcommand.hh"
45#include "G4UnitsTable.hh"
46#include "G4VProcess.hh"
48
49//#define G4ATTDEBUG
50#ifdef G4ATTDEBUG
51#include "G4AttCheck.hh"
52#endif
53
54#include <sstream>
55
57{
59 return _instance;
60}
61
63{
64}
65
67 : G4Trajectory(aTrack) // Note: this initialises the base class data
68 // members and, unfortunately but never mind,
69 // creates a G4TrajectoryPoint in
70 // TrajectoryPointContainer that we cannot
71 // access because it's private. We store the
72 // same information (plus more) in a
73 // G4RichTrajectoryPoint in the
74 // RichTrajectoryPointsContainer
75{
76 fpInitialVolume = aTrack->GetTouchableHandle();
77 fpInitialNextVolume = aTrack->GetNextTouchableHandle();
78 fpCreatorProcess = aTrack->GetCreatorProcess();
79 fCreatorModelID = aTrack->GetCreatorModelID();
80
81 // On construction, set final values to initial values.
82 // Final values are updated at the addition of every step - see AppendStep.
83 //
84 fpFinalVolume = aTrack->GetTouchableHandle();
85 fpFinalNextVolume = aTrack->GetNextTouchableHandle();
86 fpEndingProcess = aTrack->GetCreatorProcess();
87 fFinalKineticEnergy = aTrack->GetKineticEnergy();
88
89 // Insert the first rich trajectory point (see note above)...
90 //
91 fpRichPointsContainer = new RichTrajectoryPointsContainer;
92 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
93}
94
96 : G4Trajectory(right)
97{
98 fpInitialVolume = right.fpInitialVolume;
99 fpInitialNextVolume = right.fpInitialNextVolume;
100 fpCreatorProcess = right.fpCreatorProcess;
101 fCreatorModelID = right.fCreatorModelID;
102 fpFinalVolume = right.fpFinalVolume;
103 fpFinalNextVolume = right.fpFinalNextVolume;
104 fpEndingProcess = right.fpEndingProcess;
105 fFinalKineticEnergy = right.fFinalKineticEnergy;
106 fpRichPointsContainer = new RichTrajectoryPointsContainer;
107 for(std::size_t i=0; i<right.fpRichPointsContainer->size(); ++i)
108 {
109 G4RichTrajectoryPoint* rightPoint =
110 (G4RichTrajectoryPoint*)((*(right.fpRichPointsContainer))[i]);
111 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
112 }
113}
114
116{
117 if (fpRichPointsContainer)
118 {
119 for(std::size_t i=0; i<fpRichPointsContainer->size(); ++i)
120 {
121 delete (*fpRichPointsContainer)[i];
122 }
123 fpRichPointsContainer->clear();
124 delete fpRichPointsContainer;
125 }
126}
127
129{
130 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
131
132 // Except for first step, which is a sort of virtual step to start
133 // the track, compute the final values...
134 //
135 const G4Track* track = aStep->GetTrack();
136 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
137 if (track->GetCurrentStepNumber() > 0)
138 {
139 fpFinalVolume = track->GetTouchableHandle();
140 fpFinalNextVolume = track->GetNextTouchableHandle();
141 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
142 fFinalKineticEnergy = aStep->GetPreStepPoint()->GetKineticEnergy()
143 - aStep->GetTotalEnergyDeposit();
144 }
145}
146
148{
149 if(secondTrajectory == nullptr) return;
150
151 G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
152 G4int ent = seco->GetPointEntries();
153 for(G4int i=1; i<ent; ++i)
154 {
155 // initial point of the second trajectory should not be merged
156 //
157 fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
158 }
159 delete (*seco->fpRichPointsContainer)[0];
160 seco->fpRichPointsContainer->clear();
161}
162
163void G4RichTrajectory::ShowTrajectory(std::ostream& os) const
164{
165 // Invoke the default implementation in G4VTrajectory...
166 //
168
169 // ... or override with your own code here.
170}
171
173{
174 // Invoke the default implementation in G4VTrajectory...
175 //
177
178 // ... or override with your own code here.
179}
180
181const std::map<G4String,G4AttDef>* G4RichTrajectory::GetAttDefs() const
182{
183 G4bool isNew;
184 std::map<G4String,G4AttDef>* store
185 = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
186 if (isNew)
187 {
188 // Get att defs from base class...
189 //
190 *store = *(G4Trajectory::GetAttDefs());
191
192 G4String ID;
193
194 ID = "IVPath";
195 (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
196 "Physics","","G4String");
197
198 ID = "INVPath";
199 (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
200 "Physics","","G4String");
201
202 ID = "CPN";
203 (*store)[ID] = G4AttDef(ID,"Creator Process Name",
204 "Physics","","G4String");
205
206 ID = "CPTN";
207 (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
208 "Physics","","G4String");
209
210 ID = "CMID";
211 (*store)[ID] = G4AttDef(ID,"Creator Model ID",
212 "Physics","","G4int");
213
214 ID = "CMN";
215 (*store)[ID] = G4AttDef(ID,"Creator Model Name",
216 "Physics","","G4String");
217
218 ID = "FVPath";
219 (*store)[ID] = G4AttDef(ID,"Final Volume Path",
220 "Physics","","G4String");
221
222 ID = "FNVPath";
223 (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
224 "Physics","","G4String");
225
226 ID = "EPN";
227 (*store)[ID] = G4AttDef(ID,"Ending Process Name",
228 "Physics","","G4String");
229
230 ID = "EPTN";
231 (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
232 "Physics","","G4String");
233
234 ID = "FKE";
235 (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
236 "Physics","G4BestUnit","G4double");
237 }
238
239 return store;
240}
241
242static G4String Path(const G4TouchableHandle& th)
243{
244 std::ostringstream oss;
245 G4int depth = th->GetHistoryDepth();
246 for (G4int i = depth; i >= 0; --i)
247 {
248 oss << th->GetVolume(i)->GetName()
249 << ':' << th->GetCopyNumber(i);
250 if (i != 0) oss << '/';
251 }
252 return oss.str();
253}
254
255std::vector<G4AttValue>* G4RichTrajectory::CreateAttValues() const
256{
257 // Create base class att values...
258 std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
259
260 if (fpInitialVolume && fpInitialVolume->GetVolume())
261 {
262 values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
263 }
264 else
265 {
266 values->push_back(G4AttValue("IVPath","None",""));
267 }
268
269 if (fpInitialNextVolume && fpInitialNextVolume->GetVolume())
270 {
271 values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
272 }
273 else
274 {
275 values->push_back(G4AttValue("INVPath","None",""));
276 }
277
278 if (fpCreatorProcess != nullptr)
279 {
280 values->push_back
281 (G4AttValue("CPN",fpCreatorProcess->GetProcessName(),""));
282 G4ProcessType type = fpCreatorProcess->GetProcessType();
283 values->push_back
285 values->push_back
286 (G4AttValue("CMID",G4UIcommand::ConvertToString(fCreatorModelID),""));
287 const G4String& creatorModelName =
289 values->push_back(G4AttValue("CMN",creatorModelName,""));
290 }
291 else
292 {
293 values->push_back(G4AttValue("CPN","None",""));
294 values->push_back(G4AttValue("CPTN","None",""));
295 values->push_back(G4AttValue("CMID","None",""));
296 values->push_back(G4AttValue("CMN","None",""));
297 }
298
299 if (fpFinalVolume && fpFinalVolume->GetVolume())
300 {
301 values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
302 }
303 else
304 {
305 values->push_back(G4AttValue("FVPath","None",""));
306 }
307
308 if (fpFinalNextVolume && fpFinalNextVolume->GetVolume())
309 {
310 values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
311 }
312 else
313 {
314 values->push_back(G4AttValue("FNVPath","None",""));
315 }
316
317 if (fpEndingProcess != nullptr)
318 {
319 values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
320 G4ProcessType type = fpEndingProcess->GetProcessType();
321 values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
322 }
323 else
324 {
325 values->push_back(G4AttValue("EPN","None",""));
326 values->push_back(G4AttValue("EPTN","None",""));
327 }
328
329 values->push_back
330 (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
331
332#ifdef G4ATTDEBUG
333 G4cout << G4AttCheck(values,GetAttDefs());
334#endif
335
336 return values;
337}
G4ProcessType
G4Allocator< G4RichTrajectory > *& aRichTrajectoryAllocator()
#define G4BestUnit(a, b)
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
static const G4String GetModelNameFromID(const G4int modelID)
virtual ~G4RichTrajectory()
void ShowTrajectory(std::ostream &os=G4cout) const
void AppendStep(const G4Step *aStep)
virtual std::vector< G4AttValue > * CreateAttValues() const
void DrawTrajectory() const
void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4int GetPointEntries() const
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:62
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
G4int GetCreatorModelID() const
G4int GetCurrentStepNumber() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:446
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
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual void DrawTrajectory() const
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
#define G4ThreadLocalStatic
Definition: tls.hh:76