Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ClonedSmoothTrajectory.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// G4ClonedSmoothTrajectory class implementation
27//
28// Makoto Asai (JLab) - Oct.2024
29// --------------------------------------------------------------------
30
32#include "G4SmoothTrajectory.hh"
33
34#include "G4AttDef.hh"
35#include "G4AttDefStore.hh"
36#include "G4AttValue.hh"
37#include "G4ParticleTable.hh"
39#include "G4UIcommand.hh"
40#include "G4UnitsTable.hh"
41#include "G4AutoLock.hh"
42
43// #define G4ATTDEBUG
44#ifdef G4ATTDEBUG
45# include "G4AttCheck.hh"
46#endif
47
53
55{
56 ParticleName = right.ParticleName;
57 PDGCharge = right.PDGCharge;
58 PDGEncoding = right.PDGEncoding;
59 fTrackID = right.fTrackID;
60 fParentID = right.fParentID;
61 initialKineticEnergy = right.initialKineticEnergy;
62 initialMomentum = right.initialMomentum;
63 positionRecord = new G4TrajectoryPointContainer();
64
65 for (auto& i : *right.positionRecord) {
66 auto rightPoint = (G4ClonedSmoothTrajectoryPoint*)i;
67 positionRecord->push_back(new G4ClonedSmoothTrajectoryPoint(*rightPoint));
68 }
69}
70
72{
73 if (positionRecord != nullptr) {
74 for (auto& i : *positionRecord) {
75 delete i;
76 }
77 positionRecord->clear();
78 delete positionRecord;
79 }
80}
81
82void G4ClonedSmoothTrajectory::ShowTrajectory(std::ostream& os) const
83{
84 // Invoke the default implementation in G4VTrajectory...
85 //
87
88 // ... or override with your own code here.
89}
90
92{
93 // Invoke the default implementation in G4VTrajectory...
94 //
96
97 // ... or override with your own code here.
98}
99
100const std::map<G4String, G4AttDef>* G4ClonedSmoothTrajectory::GetAttDefs() const
101{
102 G4bool isNew;
103 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4ClonedSmoothTrajectory", isNew);
104 if (isNew) {
105 G4String ID("ID");
106 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
107
108 G4String PID("PID");
109 (*store)[PID] = G4AttDef(PID, "Parent ID", "Physics", "", "G4int");
110
111 G4String PN("PN");
112 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
113
114 G4String Ch("Ch");
115 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
116
117 G4String PDG("PDG");
118 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
119
120 G4String IKE("IKE");
121 (*store)[IKE] = G4AttDef(IKE, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
122
123 G4String IMom("IMom");
124 (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
125
126 G4String IMag("IMag");
127 (*store)[IMag] =
128 G4AttDef(IMag, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
129
130 G4String NTP("NTP");
131 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Physics", "", "G4int");
132 }
133 return store;
134}
135
136std::vector<G4AttValue>* G4ClonedSmoothTrajectory::CreateAttValues() const
137{
138 auto values = new std::vector<G4AttValue>;
139
140 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
141
142 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
143
144 values->push_back(G4AttValue("PN", ParticleName, ""));
145
146 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
147
148 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
149
150 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
151
152 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
153
154 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
155
156 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
157
158#ifdef G4ATTDEBUG
159 G4cout << G4AttCheck(values, GetAttDefs());
160#endif
161
162 return values;
163}
164
166{
167 positionRecord->push_back(new G4ClonedSmoothTrajectoryPoint(
169}
170
175
177{
178 if (secondTrajectory == nullptr) return;
179
180 auto seco = (G4ClonedSmoothTrajectory*)secondTrajectory;
181 G4int ent = seco->GetPointEntries();
182
183 // initial point of the second trajectory should not be merged
184 //
185 for (G4int i = 1; i < ent; ++i) {
186 positionRecord->push_back((*(seco->positionRecord))[i]);
187 }
188 delete (*seco->positionRecord)[0];
189 seco->positionRecord->clear();
190}
G4Allocator< G4ClonedSmoothTrajectory > *& aClonedSmoothTrajectoryAllocator()
#define G4BestUnit(a, b)
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
void DrawTrajectory() const override
G4int GetPointEntries() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
void ShowTrajectory(std::ostream &os=G4cout) const override
std::vector< G4AttValue > * CreateAttValues() const override
G4ParticleDefinition * GetParticleDefinition()
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
void AppendStep(const G4Step *aStep) override
G4ClonedSmoothTrajectory()=default
static G4ParticleTable * GetParticleTable()
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4StepPoint * GetPostStepPoint() const
static G4String ConvertToString(G4bool boolVal)
G4VTrajectory()=default
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual void DrawTrajectory() const
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)