Geant4 9.6.0
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//
27// $Id$
28//
29// ---------------------------------------------------------------
30//
31// G4RichTrajectory.cc
32//
33// Contact:
34// Questions and comments on G4Trajectory, on which this is based,
35// should be sent to
36// Katsuya Amako (e-mail: [email protected])
37// Makoto Asai (e-mail: [email protected])
38// Takashi Sasaki (e-mail: [email protected])
39// and on the extended code to:
40// John Allison (e-mail: [email protected])
41// Joseph Perl (e-mail: [email protected])
42//
43// ---------------------------------------------------------------
44
45#include "G4RichTrajectory.hh"
47#include "G4AttDefStore.hh"
48#include "G4AttDef.hh"
49#include "G4AttValue.hh"
50#include "G4UnitsTable.hh"
51#include "G4VProcess.hh"
52
53//#define G4ATTDEBUG
54#ifdef G4ATTDEBUG
55#include "G4AttCheck.hh"
56#endif
57
58#include <sstream>
59
61
63 fpRichPointsContainer(0),
64 fpCreatorProcess(0),
65 fpEndingProcess(0),
66 fFinalKineticEnergy(0.)
67{}
68
70 G4Trajectory(aTrack) // Note: this initialises the base class data
71 // members and, unfortunately but never mind,
72 // creates a G4TrajectoryPoint in
73 // TrajectoryPointContainer that we cannot
74 // access because it's private. We store the
75 // same information (plus more) in a
76 // G4RichTrajectoryPoint in the
77 // RichTrajectoryPointsContainer
78{
79 fpInitialVolume = aTrack->GetTouchableHandle();
80 fpInitialNextVolume = aTrack->GetNextTouchableHandle();
81 fpCreatorProcess = aTrack->GetCreatorProcess();
82 // On construction, set final values to initial values.
83 // Final values are updated at the addition of every step - see AppendStep.
84 fpFinalVolume = aTrack->GetTouchableHandle();
85 fpFinalNextVolume = aTrack->GetNextTouchableHandle();
86 fpEndingProcess = aTrack->GetCreatorProcess();
87 fFinalKineticEnergy = aTrack->GetKineticEnergy();
88 // Insert the first rich trajectory point (see note above)...
89 fpRichPointsContainer = new RichTrajectoryPointsContainer;
90 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aTrack));
91}
92
94 G4Trajectory(right)
95{
96 fpInitialVolume = right.fpInitialVolume;
97 fpInitialNextVolume = right.fpInitialNextVolume;
98 fpCreatorProcess = right.fpCreatorProcess;
99 fpFinalVolume = right.fpFinalVolume;
100 fpFinalNextVolume = right.fpFinalNextVolume;
101 fpEndingProcess = right.fpEndingProcess;
102 fFinalKineticEnergy = right.fFinalKineticEnergy;
103 fpRichPointsContainer = new RichTrajectoryPointsContainer;
104 for(size_t i=0;i<right.fpRichPointsContainer->size();i++)
105 {
106 G4RichTrajectoryPoint* rightPoint =
107 (G4RichTrajectoryPoint*)((*(right.fpRichPointsContainer))[i]);
108 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
109 }
110}
111
113{
114 if (fpRichPointsContainer) {
115 // fpRichPointsContainer->clearAndDestroy();
116 size_t i;
117 for(i=0;i<fpRichPointsContainer->size();i++){
118 delete (*fpRichPointsContainer)[i];
119 }
120 fpRichPointsContainer->clear();
121 delete fpRichPointsContainer;
122 }
123}
124
126{
127 fpRichPointsContainer->push_back(new G4RichTrajectoryPoint(aStep));
128 // Except for first step, which is a sort of virtual step to start
129 // the track, compute the final values...
130 const G4Track* track = aStep->GetTrack();
131 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
132 if (track->GetCurrentStepNumber() > 0) {
133 fpFinalVolume = track->GetTouchableHandle();
134 fpFinalNextVolume = track->GetNextTouchableHandle();
135 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
136 fFinalKineticEnergy =
138 aStep->GetTotalEnergyDeposit();
139 }
140}
141
143{
144 if(!secondTrajectory) return;
145
146 G4RichTrajectory* seco = (G4RichTrajectory*)secondTrajectory;
147 G4int ent = seco->GetPointEntries();
148 for(G4int i=1;i<ent;i++) {
149 // initial point of the second trajectory should not be merged
150 fpRichPointsContainer->push_back((*(seco->fpRichPointsContainer))[i]);
151 // fpRichPointsContainer->push_back(seco->fpRichPointsContainer->removeAt(1));
152 }
153 delete (*seco->fpRichPointsContainer)[0];
154 seco->fpRichPointsContainer->clear();
155}
156
157const std::map<G4String,G4AttDef>* G4RichTrajectory::GetAttDefs() const
158{
159 G4bool isNew;
160 std::map<G4String,G4AttDef>* store
161 = G4AttDefStore::GetInstance("G4RichTrajectory",isNew);
162 if (isNew) {
163
164 // Get att defs from base class...
165 *store = *(G4Trajectory::GetAttDefs());
166
167 G4String ID;
168
169 ID = "IVPath";
170 (*store)[ID] = G4AttDef(ID,"Initial Volume Path",
171 "Physics","","G4String");
172
173 ID = "INVPath";
174 (*store)[ID] = G4AttDef(ID,"Initial Next Volume Path",
175 "Physics","","G4String");
176
177 ID = "CPN";
178 (*store)[ID] = G4AttDef(ID,"Creator Process Name",
179 "Physics","","G4String");
180
181 ID = "CPTN";
182 (*store)[ID] = G4AttDef(ID,"Creator Process Type Name",
183 "Physics","","G4String");
184
185 ID = "FVPath";
186 (*store)[ID] = G4AttDef(ID,"Final Volume Path",
187 "Physics","","G4String");
188
189 ID = "FNVPath";
190 (*store)[ID] = G4AttDef(ID,"Final Next Volume Path",
191 "Physics","","G4String");
192
193 ID = "EPN";
194 (*store)[ID] = G4AttDef(ID,"Ending Process Name",
195 "Physics","","G4String");
196
197 ID = "EPTN";
198 (*store)[ID] = G4AttDef(ID,"Ending Process Type Name",
199 "Physics","","G4String");
200
201 ID = "FKE";
202 (*store)[ID] = G4AttDef(ID,"Final kinetic energy",
203 "Physics","G4BestUnit","G4double");
204
205 }
206
207 return store;
208}
209
210static G4String Path(const G4TouchableHandle& th)
211{
212 std::ostringstream oss;
213 G4int depth = th->GetHistoryDepth();
214 for (G4int i = depth; i >= 0; --i) {
215 oss << th->GetVolume(i)->GetName()
216 << ':' << th->GetCopyNumber(i);
217 if (i != 0) oss << '/';
218 }
219 return oss.str();
220}
221
222std::vector<G4AttValue>* G4RichTrajectory::CreateAttValues() const
223{
224 // Create base class att values...
225 std::vector<G4AttValue>* values = G4Trajectory::CreateAttValues();
226
227 if (fpInitialVolume && fpInitialVolume->GetVolume()) {
228 values->push_back(G4AttValue("IVPath",Path(fpInitialVolume),""));
229 } else {
230 values->push_back(G4AttValue("IVPath","None",""));
231 }
232
233 if (fpInitialNextVolume && fpInitialNextVolume->GetVolume()) {
234 values->push_back(G4AttValue("INVPath",Path(fpInitialNextVolume),""));
235 } else {
236 values->push_back(G4AttValue("INVPath","None",""));
237 }
238
239 if (fpCreatorProcess) {
240 values->push_back(G4AttValue("CPN",fpCreatorProcess->GetProcessName(),""));
241 G4ProcessType type = fpCreatorProcess->GetProcessType();
242 values->push_back(G4AttValue("CPTN",G4VProcess::GetProcessTypeName(type),""));
243 } else {
244 values->push_back(G4AttValue("CPN","None",""));
245 values->push_back(G4AttValue("CPTN","None",""));
246 }
247
248 if (fpFinalVolume && fpFinalVolume->GetVolume()) {
249 values->push_back(G4AttValue("FVPath",Path(fpFinalVolume),""));
250 } else {
251 values->push_back(G4AttValue("FVPath","None",""));
252 }
253
254 if (fpFinalNextVolume && fpFinalNextVolume->GetVolume()) {
255 values->push_back(G4AttValue("FNVPath",Path(fpFinalNextVolume),""));
256 } else {
257 values->push_back(G4AttValue("FNVPath","None",""));
258 }
259
260 if (fpEndingProcess) {
261 values->push_back(G4AttValue("EPN",fpEndingProcess->GetProcessName(),""));
262 G4ProcessType type = fpEndingProcess->GetProcessType();
263 values->push_back(G4AttValue("EPTN",G4VProcess::GetProcessTypeName(type),""));
264 } else {
265 values->push_back(G4AttValue("EPN","None",""));
266 values->push_back(G4AttValue("EPTN","None",""));
267 }
268
269 values->push_back
270 (G4AttValue("FKE",G4BestUnit(fFinalKineticEnergy,"Energy"),""));
271
272#ifdef G4ATTDEBUG
273 G4cout << G4AttCheck(values,GetAttDefs());
274#endif
275
276 return values;
277}
G4ProcessType
G4Allocator< G4RichTrajectory > aRichTrajectoryAllocator
std::vector< G4VTrajectoryPoint * > RichTrajectoryPointsContainer
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4DLLIMPORT std::ostream G4cout
virtual ~G4RichTrajectory()
void AppendStep(const G4Step *aStep)
virtual std::vector< G4AttValue > * CreateAttValues() const
int GetPointEntries() const
void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
const G4VProcess * GetProcessDefinedStep() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:78
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetTotalEnergyDeposit() const
G4StepPoint * GetPostStepPoint() const
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() 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
const G4String & GetName() const
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:150
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:385
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
G4int GetCopyNumber(G4int depth=0) const
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:44
virtual G4int GetHistoryDepth() const
Definition: G4VTouchable.cc:79
std::map< G4String, G4AttDef > * GetInstance(G4String storeKey, G4bool &isNew)