Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ClonedRichTrajectory.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// G4ClonedRichTrajectory class implementation
27//
28// Makoto Asai (JLab) - Oct.2024
29// --------------------------------------------------------------------
30
32#include "G4RichTrajectory.hh"
33
34#include "G4ParticleTable.hh"
35#include "G4AttDef.hh"
36#include "G4AttDefStore.hh"
37#include "G4AttValue.hh"
40#include "G4UIcommand.hh"
41#include "G4UnitsTable.hh"
42#include "G4VProcess.hh"
43
44// #define G4ATTDEBUG
45#ifdef G4ATTDEBUG
46# include "G4AttCheck.hh"
47#endif
48
49#include <sstream>
50
56
58{
59 ParticleName = right.ParticleName;
60 PDGCharge = right.PDGCharge;
61 PDGEncoding = right.PDGEncoding;
62 fTrackID = right.fTrackID;
63 fParentID = right.fParentID;
64 initialKineticEnergy = right.initialKineticEnergy;
65 initialMomentum = right.initialMomentum;
66 //positionRecord = new G4ClonedTrajectoryPointContainer();
67 //for (auto& i : *right.positionRecord) {
68 // auto rightPoint = (G4ClonedTrajectoryPoint*)i;
69 // positionRecord->push_back(new G4ClonedTrajectoryPoint(*rightPoint));
70 //}
71 fpInitialVolume = right.fpInitialVolume;
72 fpInitialNextVolume = right.fpInitialNextVolume;
73 fpCreatorProcess = right.fpCreatorProcess;
74 fCreatorModelID = right.fCreatorModelID;
75 fpFinalVolume = right.fpFinalVolume;
76 fpFinalNextVolume = right.fpFinalNextVolume;
77 fpEndingProcess = right.fpEndingProcess;
78 fFinalKineticEnergy = right.fFinalKineticEnergy;
79 if(right.fpRichPointContainer!=nullptr && right.fpRichPointContainer->size()>0)
80 {
81 fpRichPointContainer = new G4TrajectoryPointContainer;
82 for (auto& i : *right.fpRichPointContainer) {
83 auto rightPoint = (G4RichTrajectoryPoint*)i;
84 fpRichPointContainer->push_back(new G4ClonedRichTrajectoryPoint(*rightPoint));
85 }
86 }
87}
88
90{
91 if (fpRichPointContainer != nullptr) {
92 for (auto& i : *fpRichPointContainer) {
93 delete i;
94 }
95 fpRichPointContainer->clear();
96 delete fpRichPointContainer;
97 }
98}
99
101{
102 fpRichPointContainer->push_back(new G4ClonedRichTrajectoryPoint(aStep));
103
104 // Except for first step, which is a sort of virtual step to start
105 // the track, compute the final values...
106 //
107 const G4Track* track = aStep->GetTrack();
108 const G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
109 if (track->GetCurrentStepNumber() > 0) {
110 fpFinalVolume = track->GetTouchableHandle();
111 fpFinalNextVolume = track->GetNextTouchableHandle();
112 fpEndingProcess = postStepPoint->GetProcessDefinedStep();
113 fFinalKineticEnergy =
115 }
116}
117
119{
120 if (secondTrajectory == nullptr) return;
121
122 auto seco = (G4ClonedRichTrajectory*)secondTrajectory;
123 G4int ent = seco->GetPointEntries();
124 for (G4int i = 1; i < ent; ++i) {
125 // initial point of the second trajectory should not be merged
126 //
127 fpRichPointContainer->push_back((*(seco->fpRichPointContainer))[i]);
128 }
129 delete (*seco->fpRichPointContainer)[0];
130 seco->fpRichPointContainer->clear();
131}
132
133void G4ClonedRichTrajectory::ShowTrajectory(std::ostream& os) const
134{
135 // Invoke the default implementation in G4VTrajectory...
136 //
138
139 // ... or override with your own code here.
140}
141
143{
144 // Invoke the default implementation in G4VTrajectory...
145 //
147
148 // ... or override with your own code here.
149}
150
151const std::map<G4String, G4AttDef>* G4ClonedRichTrajectory::GetAttDefs() const
152{
153 G4bool isNew;
154 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4ClonedRichTrajectory", isNew);
155 if (isNew) {
156 G4String ID;
157
158 ID = "ID";
159 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
160
161 ID = "PID";
162 (*store)[ID] = G4AttDef(ID, "Parent ID", "Physics", "", "G4int");
163
164 ID = "PN";
165 (*store)[ID] = G4AttDef(ID, "Particle Name", "Physics", "", "G4String");
166
167 ID = "Ch";
168 (*store)[ID] = G4AttDef(ID, "Charge", "Physics", "e+", "G4double");
169
170 ID = "PDG";
171 (*store)[ID] = G4AttDef(ID, "PDG Encoding", "Physics", "", "G4int");
172
173 ID = "IKE";
174 (*store)[ID] = G4AttDef(ID, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
175
176 ID = "IMom";
177 (*store)[ID] = G4AttDef(ID, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
178
179 ID = "IMag";
180 (*store)[ID] = G4AttDef(ID, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
181
182 ID = "NTP";
183 (*store)[ID] = G4AttDef(ID, "No. of points", "Physics", "", "G4int");
184
185 ID = "IVPath";
186 (*store)[ID] = G4AttDef(ID, "Initial Volume Path", "Physics", "", "G4String");
187
188 ID = "INVPath";
189 (*store)[ID] = G4AttDef(ID, "Initial Next Volume Path", "Physics", "", "G4String");
190
191 ID = "CPN";
192 (*store)[ID] = G4AttDef(ID, "Creator Process Name", "Physics", "", "G4String");
193
194 ID = "CPTN";
195 (*store)[ID] = G4AttDef(ID, "Creator Process Type Name", "Physics", "", "G4String");
196
197 ID = "CMID";
198 (*store)[ID] = G4AttDef(ID, "Creator Model ID", "Physics", "", "G4int");
199
200 ID = "CMN";
201 (*store)[ID] = G4AttDef(ID, "Creator Model Name", "Physics", "", "G4String");
202
203 ID = "FVPath";
204 (*store)[ID] = G4AttDef(ID, "Final Volume Path", "Physics", "", "G4String");
205
206 ID = "FNVPath";
207 (*store)[ID] = G4AttDef(ID, "Final Next Volume Path", "Physics", "", "G4String");
208
209 ID = "EPN";
210 (*store)[ID] = G4AttDef(ID, "Ending Process Name", "Physics", "", "G4String");
211
212 ID = "EPTN";
213 (*store)[ID] = G4AttDef(ID, "Ending Process Type Name", "Physics", "", "G4String");
214
215 ID = "FKE";
216 (*store)[ID] = G4AttDef(ID, "Final kinetic energy", "Physics", "G4BestUnit", "G4double");
217 }
218
219 return store;
220}
221
222static G4String Path(const G4TouchableHandle& th)
223{
224 std::ostringstream oss;
225 G4int depth = th->GetHistoryDepth();
226 for (G4int i = depth; i >= 0; --i) {
227 oss << th->GetVolume(i)->GetName() << ':' << th->GetCopyNumber(i);
228 if (i != 0) oss << '/';
229 }
230 return oss.str();
231}
232
233std::vector<G4AttValue>* G4ClonedRichTrajectory::CreateAttValues() const
234{
235 // Create base class att values...
236 //std::vector<G4AttValue>* values = G4VTrajectory::CreateAttValues();
237 auto values = new std::vector<G4AttValue>;
238 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
239 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
240 values->push_back(G4AttValue("PN", ParticleName, ""));
241 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
242 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
243 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
244 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
245 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
246 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
247
248 if (fpInitialVolume && (fpInitialVolume->GetVolume() != nullptr)) {
249 values->push_back(G4AttValue("IVPath", Path(fpInitialVolume), ""));
250 }
251 else {
252 values->push_back(G4AttValue("IVPath", "None", ""));
253 }
254
255 if (fpInitialNextVolume && (fpInitialNextVolume->GetVolume() != nullptr)) {
256 values->push_back(G4AttValue("INVPath", Path(fpInitialNextVolume), ""));
257 }
258 else {
259 values->push_back(G4AttValue("INVPath", "None", ""));
260 }
261
262 if (fpCreatorProcess != nullptr) {
263 values->push_back(G4AttValue("CPN", fpCreatorProcess->GetProcessName(), ""));
264 G4ProcessType type = fpCreatorProcess->GetProcessType();
265 values->push_back(G4AttValue("CPTN", G4VProcess::GetProcessTypeName(type), ""));
266 values->push_back(G4AttValue("CMID", G4UIcommand::ConvertToString(fCreatorModelID), ""));
267 const G4String& creatorModelName = G4PhysicsModelCatalog::GetModelNameFromID(fCreatorModelID);
268 values->push_back(G4AttValue("CMN", creatorModelName, ""));
269 }
270 else {
271 values->push_back(G4AttValue("CPN", "None", ""));
272 values->push_back(G4AttValue("CPTN", "None", ""));
273 values->push_back(G4AttValue("CMID", "None", ""));
274 values->push_back(G4AttValue("CMN", "None", ""));
275 }
276
277 if (fpFinalVolume && (fpFinalVolume->GetVolume() != nullptr)) {
278 values->push_back(G4AttValue("FVPath", Path(fpFinalVolume), ""));
279 }
280 else {
281 values->push_back(G4AttValue("FVPath", "None", ""));
282 }
283
284 if (fpFinalNextVolume && (fpFinalNextVolume->GetVolume() != nullptr)) {
285 values->push_back(G4AttValue("FNVPath", Path(fpFinalNextVolume), ""));
286 }
287 else {
288 values->push_back(G4AttValue("FNVPath", "None", ""));
289 }
290
291 if (fpEndingProcess != nullptr) {
292 values->push_back(G4AttValue("EPN", fpEndingProcess->GetProcessName(), ""));
293 G4ProcessType type = fpEndingProcess->GetProcessType();
294 values->push_back(G4AttValue("EPTN", G4VProcess::GetProcessTypeName(type), ""));
295 }
296 else {
297 values->push_back(G4AttValue("EPN", "None", ""));
298 values->push_back(G4AttValue("EPTN", "None", ""));
299 }
300
301 values->push_back(G4AttValue("FKE", G4BestUnit(fFinalKineticEnergy, "Energy"), ""));
302
303#ifdef G4ATTDEBUG
304 G4cout << G4AttCheck(values, GetAttDefs());
305#endif
306
307 return values;
308}
309
314
315
G4Allocator< G4ClonedRichTrajectory > *& aClonedRichTrajectoryAllocator()
G4ProcessType
#define G4BestUnit(a, b)
G4ReferenceCountedHandle< G4VTouchable > G4TouchableHandle
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
void AppendStep(const G4Step *aStep) override
void DrawTrajectory() const override
std::vector< G4AttValue > * CreateAttValues() const override
G4ClonedRichTrajectory()=default
void ShowTrajectory(std::ostream &os=G4cout) const override
G4ParticleDefinition * GetParticleDefinition()
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4int GetPointEntries() const override
static G4ParticleTable * GetParticleTable()
static const G4String GetModelNameFromID(const G4int modelID)
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
G4int GetCurrentStepNumber() const
const G4TouchableHandle & GetTouchableHandle() 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)