Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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: [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
39#include "G4RichTrajectory.hh"
40
41#include "G4AttDef.hh"
42#include "G4AttDefStore.hh"
43#include "G4AttValue.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 : 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}
91
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}
108
110{
111 if (fpRichPointsContainer != nullptr) {
112 for (auto& i : *fpRichPointsContainer) {
113 delete i;
114 }
115 fpRichPointsContainer->clear();
116 delete fpRichPointsContainer;
117 }
118}
119
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}
137
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}
152
153void G4RichTrajectory::ShowTrajectory(std::ostream& os) const
154{
155 // Invoke the default implementation in G4VTrajectory...
156 //
158
159 // ... or override with your own code here.
160}
161
163{
164 // Invoke the default implementation in G4VTrajectory...
165 //
167
168 // ... or override with your own code here.
169}
170
171const std::map<G4String, G4AttDef>* G4RichTrajectory::GetAttDefs() const
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}
218
219static G4String Path(const G4TouchableHandle& th)
220{
221 std::ostringstream oss;
222 G4int depth = th->GetHistoryDepth();
223 for (G4int i = depth; i >= 0; --i) {
224 oss << th->GetVolume(i)->GetName() << ':' << th->GetCopyNumber(i);
225 if (i != 0) oss << '/';
226 }
227 return oss.str();
228}
229
230std::vector<G4AttValue>* G4RichTrajectory::CreateAttValues() const
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
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)
void AppendStep(const G4Step *aStep) override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
void DrawTrajectory() const override
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
~G4RichTrajectory() override
G4RichTrajectory()=default
void ShowTrajectory(std::ostream &os=G4cout) const override
std::vector< G4AttValue > * CreateAttValues() const override
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
G4int GetCopyNumber(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
virtual G4int GetHistoryDepth() const
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
G4int GetCreatorModelID() const
G4int GetCurrentStepNumber() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
const std::map< G4String, G4AttDef > * GetAttDefs() const override
std::vector< G4AttValue > * CreateAttValues() const override
static G4String ConvertToString(G4bool boolVal)
const G4String & GetName() const
static const G4String & GetProcessTypeName(G4ProcessType)
G4ProcessType GetProcessType() const
const G4String & GetProcessName() const
virtual G4int GetPointEntries() const =0
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