Geant4 11.3.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// 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"
41
42#include "G4ParticleTable.hh"
43#include "G4AttDef.hh"
44#include "G4AttDefStore.hh"
45#include "G4AttValue.hh"
48#include "G4UIcommand.hh"
49#include "G4UnitsTable.hh"
50#include "G4VProcess.hh"
51
52namespace {
53 G4Mutex CloneRichTrajectoryMutex = G4MUTEX_INITIALIZER;
54}
55
56// #define G4ATTDEBUG
57#ifdef G4ATTDEBUG
58# include "G4AttCheck.hh"
59#endif
60
61#include <sstream>
62
68
70{
71 G4ParticleDefinition* fpParticleDefinition = aTrack->GetDefinition();
72 ParticleName = fpParticleDefinition->GetParticleName();
73 PDGCharge = fpParticleDefinition->GetPDGCharge();
74 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
75 fTrackID = aTrack->GetTrackID();
76 fParentID = aTrack->GetParentID();
77 initialKineticEnergy = aTrack->GetKineticEnergy();
78 initialMomentum = aTrack->GetMomentum();
79 positionRecord = new G4TrajectoryPointContainer();
80
81 // Following is for the first trajectory point
82 positionRecord->push_back(new G4RichTrajectoryPoint(aTrack));
83
84 fpInitialVolume = aTrack->GetTouchableHandle();
85 fpInitialNextVolume = aTrack->GetNextTouchableHandle();
86 fpCreatorProcess = aTrack->GetCreatorProcess();
87 fCreatorModelID = aTrack->GetCreatorModelID();
88
89 // On construction, set final values to initial values.
90 // Final values are updated at the addition of every step - see AppendStep.
91 //
92 fpFinalVolume = aTrack->GetTouchableHandle();
93 fpFinalNextVolume = aTrack->GetNextTouchableHandle();
94 fpEndingProcess = aTrack->GetCreatorProcess();
95 fFinalKineticEnergy = aTrack->GetKineticEnergy();
96
97 // Insert the first rich trajectory point (see note above)...
98 //
99 fpRichPointContainer = new G4TrajectoryPointContainer;
100 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(aTrack));
101}
102
104{
105 ParticleName = right.ParticleName;
106 PDGCharge = right.PDGCharge;
107 PDGEncoding = right.PDGEncoding;
108 fTrackID = right.fTrackID;
109 fParentID = right.fParentID;
110 initialKineticEnergy = right.initialKineticEnergy;
111 initialMomentum = right.initialMomentum;
112 positionRecord = new G4TrajectoryPointContainer();
113
114 for (auto& i : *right.positionRecord) {
115 auto rightPoint = (G4RichTrajectoryPoint*)i;
116 positionRecord->push_back(new G4RichTrajectoryPoint(*rightPoint));
117 }
118
119 fpInitialVolume = right.fpInitialVolume;
120 fpInitialNextVolume = right.fpInitialNextVolume;
121 fpCreatorProcess = right.fpCreatorProcess;
122 fCreatorModelID = right.fCreatorModelID;
123 fpFinalVolume = right.fpFinalVolume;
124 fpFinalNextVolume = right.fpFinalNextVolume;
125 fpEndingProcess = right.fpEndingProcess;
126 fFinalKineticEnergy = right.fFinalKineticEnergy;
127 fpRichPointContainer = new G4TrajectoryPointContainer;
128 for (auto& i : *right.fpRichPointContainer) {
129 auto rightPoint = (G4RichTrajectoryPoint*)i;
130 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(*rightPoint));
131 }
132}
133
135{
136 if (fpRichPointContainer != nullptr) {
137 for (auto& i : *fpRichPointContainer) {
138 delete i;
139 }
140 fpRichPointContainer->clear();
141 delete fpRichPointContainer;
142 }
143}
144
146{
147 fpRichPointContainer->push_back(new G4RichTrajectoryPoint(aStep));
148
149 // Except for first step, which is a sort of virtual step to start
150 // the track, compute the final values...
151 //
152 const G4Track* track = aStep->GetTrack();
153 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
154 if (track->GetCurrentStepNumber() > 0) {
155 fpFinalVolume = track->GetTouchableHandle();
156 fpFinalNextVolume = track->GetNextTouchableHandle();
157 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
158 fFinalKineticEnergy =
160 }
161}
162
164{
165 if (secondTrajectory == nullptr) return;
166
167 auto seco = (G4RichTrajectory*)secondTrajectory;
168 G4int ent = seco->GetPointEntries();
169 for (G4int i = 1; i < ent; ++i) {
170 // initial point of the second trajectory should not be merged
171 //
172 fpRichPointContainer->push_back((*(seco->fpRichPointContainer))[i]);
173 }
174 delete (*seco->fpRichPointContainer)[0];
175 seco->fpRichPointContainer->clear();
176}
177
178void G4RichTrajectory::ShowTrajectory(std::ostream& os) const
179{
180 // Invoke the default implementation in G4VTrajectory...
181 //
183
184 // ... or override with your own code here.
185}
186
188{
189 // Invoke the default implementation in G4VTrajectory...
190 //
192
193 // ... or override with your own code here.
194}
195
196const std::map<G4String, G4AttDef>* G4RichTrajectory::GetAttDefs() const
197{
198 G4bool isNew;
199 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectory", isNew);
200 if (isNew) {
201 G4String ID;
202
203 ID = "ID";
204 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
205
206 ID = "PID";
207 (*store)[ID] = G4AttDef(ID, "Parent ID", "Physics", "", "G4int");
208
209 ID = "PN";
210 (*store)[ID] = G4AttDef(ID, "Particle Name", "Physics", "", "G4String");
211
212 ID = "Ch";
213 (*store)[ID] = G4AttDef(ID, "Charge", "Physics", "e+", "G4double");
214
215 ID = "PDG";
216 (*store)[ID] = G4AttDef(ID, "PDG Encoding", "Physics", "", "G4int");
217
218 ID = "IKE";
219 (*store)[ID] = G4AttDef(ID, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
220
221 ID = "IMom";
222 (*store)[ID] = G4AttDef(ID, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
223
224 ID = "IMag";
225 (*store)[ID] = G4AttDef(ID, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
226
227 ID = "NTP";
228 (*store)[ID] = G4AttDef(ID, "No. of points", "Physics", "", "G4int");
229
230 ID = "IVPath";
231 (*store)[ID] = G4AttDef(ID, "Initial Volume Path", "Physics", "", "G4String");
232
233 ID = "INVPath";
234 (*store)[ID] = G4AttDef(ID, "Initial Next Volume Path", "Physics", "", "G4String");
235
236 ID = "CPN";
237 (*store)[ID] = G4AttDef(ID, "Creator Process Name", "Physics", "", "G4String");
238
239 ID = "CPTN";
240 (*store)[ID] = G4AttDef(ID, "Creator Process Type Name", "Physics", "", "G4String");
241
242 ID = "CMID";
243 (*store)[ID] = G4AttDef(ID, "Creator Model ID", "Physics", "", "G4int");
244
245 ID = "CMN";
246 (*store)[ID] = G4AttDef(ID, "Creator Model Name", "Physics", "", "G4String");
247
248 ID = "FVPath";
249 (*store)[ID] = G4AttDef(ID, "Final Volume Path", "Physics", "", "G4String");
250
251 ID = "FNVPath";
252 (*store)[ID] = G4AttDef(ID, "Final Next Volume Path", "Physics", "", "G4String");
253
254 ID = "EPN";
255 (*store)[ID] = G4AttDef(ID, "Ending Process Name", "Physics", "", "G4String");
256
257 ID = "EPTN";
258 (*store)[ID] = G4AttDef(ID, "Ending Process Type Name", "Physics", "", "G4String");
259
260 ID = "FKE";
261 (*store)[ID] = G4AttDef(ID, "Final kinetic energy", "Physics", "G4BestUnit", "G4double");
262 }
263
264 return store;
265}
266
267static G4String Path(const G4TouchableHandle& th)
268{
269 std::ostringstream oss;
270 G4int depth = th->GetHistoryDepth();
271 for (G4int i = depth; i >= 0; --i) {
272 oss << th->GetVolume(i)->GetName() << ':' << th->GetCopyNumber(i);
273 if (i != 0) oss << '/';
274 }
275 return oss.str();
276}
277
278std::vector<G4AttValue>* G4RichTrajectory::CreateAttValues() const
279{
280 // Create base class att values...
281 //std::vector<G4AttValue>* values = G4VTrajectory::CreateAttValues();
282 auto values = new std::vector<G4AttValue>;
283 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
284 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
285 values->push_back(G4AttValue("PN", ParticleName, ""));
286 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
287 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
288 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
289 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
290 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
291 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
292
293 if (fpInitialVolume && (fpInitialVolume->GetVolume() != nullptr)) {
294 values->push_back(G4AttValue("IVPath", Path(fpInitialVolume), ""));
295 }
296 else {
297 values->push_back(G4AttValue("IVPath", "None", ""));
298 }
299
300 if (fpInitialNextVolume && (fpInitialNextVolume->GetVolume() != nullptr)) {
301 values->push_back(G4AttValue("INVPath", Path(fpInitialNextVolume), ""));
302 }
303 else {
304 values->push_back(G4AttValue("INVPath", "None", ""));
305 }
306
307 if (fpCreatorProcess != nullptr) {
308 values->push_back(G4AttValue("CPN", fpCreatorProcess->GetProcessName(), ""));
309 G4ProcessType type = fpCreatorProcess->GetProcessType();
310 values->push_back(G4AttValue("CPTN", G4VProcess::GetProcessTypeName(type), ""));
311 values->push_back(G4AttValue("CMID", G4UIcommand::ConvertToString(fCreatorModelID), ""));
312 const G4String& creatorModelName = G4PhysicsModelCatalog::GetModelNameFromID(fCreatorModelID);
313 values->push_back(G4AttValue("CMN", creatorModelName, ""));
314 }
315 else {
316 values->push_back(G4AttValue("CPN", "None", ""));
317 values->push_back(G4AttValue("CPTN", "None", ""));
318 values->push_back(G4AttValue("CMID", "None", ""));
319 values->push_back(G4AttValue("CMN", "None", ""));
320 }
321
322 if (fpFinalVolume && (fpFinalVolume->GetVolume() != nullptr)) {
323 values->push_back(G4AttValue("FVPath", Path(fpFinalVolume), ""));
324 }
325 else {
326 values->push_back(G4AttValue("FVPath", "None", ""));
327 }
328
329 if (fpFinalNextVolume && (fpFinalNextVolume->GetVolume() != nullptr)) {
330 values->push_back(G4AttValue("FNVPath", Path(fpFinalNextVolume), ""));
331 }
332 else {
333 values->push_back(G4AttValue("FNVPath", "None", ""));
334 }
335
336 if (fpEndingProcess != nullptr) {
337 values->push_back(G4AttValue("EPN", fpEndingProcess->GetProcessName(), ""));
338 G4ProcessType type = fpEndingProcess->GetProcessType();
339 values->push_back(G4AttValue("EPTN", G4VProcess::GetProcessTypeName(type), ""));
340 }
341 else {
342 values->push_back(G4AttValue("EPN", "None", ""));
343 values->push_back(G4AttValue("EPTN", "None", ""));
344 }
345
346 values->push_back(G4AttValue("FKE", G4BestUnit(fFinalKineticEnergy, "Energy"), ""));
347
348#ifdef G4ATTDEBUG
349 G4cout << G4AttCheck(values, GetAttDefs());
350#endif
351
352 return values;
353}
354
359
361{
362 G4AutoLock lock(&CloneRichTrajectoryMutex);
363 auto* cloned = new G4ClonedRichTrajectory(*this);
364 return cloned;
365}
366
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4ProcessType
G4Allocator< G4RichTrajectory > *& aRichTrajectoryAllocator()
#define G4BestUnit(a, b)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleName() const
static G4ParticleTable * GetParticleTable()
static const G4String GetModelNameFromID(const G4int modelID)
void AppendStep(const G4Step *aStep) override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4int GetPointEntries() const override
void DrawTrajectory() const override
G4VTrajectory * CloneForMaster() const override
friend class G4ClonedRichTrajectory
G4ParticleDefinition * GetParticleDefinition()
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
G4int GetTrackID() const
const G4TouchableHandle & GetNextTouchableHandle() const
const G4VProcess * GetCreatorProcess() const
G4int GetCreatorModelID() const
G4int GetCurrentStepNumber() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
G4int GetParentID() const
static G4String ConvertToString(G4bool boolVal)
const G4String & GetName() const
static const G4String & GetProcessTypeName(G4ProcessType)
G4VTrajectory()=default
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