Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4SmoothTrajectory.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// G4SmoothTrajectory class implementation
27//
28// Contact:
29// Questions and comments to this code should be sent to
30// Katsuya Amako (e-mail: [email protected])
31// Makoto Asai (e-mail: [email protected])
32// Takashi Sasaki (e-mail: [email protected])
33// --------------------------------------------------------------------
34
35#include "G4SmoothTrajectory.hh"
37
38#include "G4AttDef.hh"
39#include "G4AttDefStore.hh"
40#include "G4AttValue.hh"
41#include "G4ParticleTable.hh"
43#include "G4UIcommand.hh"
44#include "G4UnitsTable.hh"
45#include "G4AutoLock.hh"
46
47namespace {
48 G4Mutex CloneSmoothTrajectoryMutex = G4MUTEX_INITIALIZER;
49}
50
51// #define G4ATTDEBUG
52#ifdef G4ATTDEBUG
53# include "G4AttCheck.hh"
54#endif
55
61
63{
64 G4ParticleDefinition* fpParticleDefinition = aTrack->GetDefinition();
65 ParticleName = fpParticleDefinition->GetParticleName();
66 PDGCharge = fpParticleDefinition->GetPDGCharge();
67 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
68 fTrackID = aTrack->GetTrackID();
69 fParentID = aTrack->GetParentID();
70 initialKineticEnergy = aTrack->GetKineticEnergy();
71 initialMomentum = aTrack->GetMomentum();
72 positionRecord = new G4TrajectoryPointContainer();
73
74 // Following is for the first trajectory point
75 //
76 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition()));
77
78 // The first point has no auxiliary points, so set the auxiliary
79 // points vector to NULL
80 //
81 positionRecord->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), nullptr));
82}
83
85{
86 ParticleName = right.ParticleName;
87 PDGCharge = right.PDGCharge;
88 PDGEncoding = right.PDGEncoding;
89 fTrackID = right.fTrackID;
90 fParentID = right.fParentID;
91 initialKineticEnergy = right.initialKineticEnergy;
92 initialMomentum = right.initialMomentum;
93 positionRecord = new G4TrajectoryPointContainer();
94
95 for (auto& i : *right.positionRecord) {
96 auto rightPoint = (G4SmoothTrajectoryPoint*)i;
97 positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint));
98 }
99}
100
102{
103 G4AutoLock lock(&CloneSmoothTrajectoryMutex);
104 auto* cloned = new G4ClonedSmoothTrajectory(*this);
105 return cloned;
106}
107
109{
110 if (positionRecord != nullptr) {
111 for (auto& i : *positionRecord) {
112 delete i;
113 }
114 positionRecord->clear();
115 delete positionRecord;
116 }
117}
118
119void G4SmoothTrajectory::ShowTrajectory(std::ostream& os) const
120{
121 // Invoke the default implementation in G4VTrajectory...
122 //
124
125 // ... or override with your own code here.
126}
127
129{
130 // Invoke the default implementation in G4VTrajectory...
131 //
133
134 // ... or override with your own code here.
135}
136
137const std::map<G4String, G4AttDef>* G4SmoothTrajectory::GetAttDefs() const
138{
139 G4bool isNew;
140 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4SmoothTrajectory", isNew);
141 if (isNew) {
142 G4String ID("ID");
143 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int");
144
145 G4String PID("PID");
146 (*store)[PID] = G4AttDef(PID, "Parent ID", "Physics", "", "G4int");
147
148 G4String PN("PN");
149 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
150
151 G4String Ch("Ch");
152 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
153
154 G4String PDG("PDG");
155 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
156
157 G4String IKE("IKE");
158 (*store)[IKE] = G4AttDef(IKE, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double");
159
160 G4String IMom("IMom");
161 (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector");
162
163 G4String IMag("IMag");
164 (*store)[IMag] =
165 G4AttDef(IMag, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double");
166
167 G4String NTP("NTP");
168 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Physics", "", "G4int");
169 }
170 return store;
171}
172
173std::vector<G4AttValue>* G4SmoothTrajectory::CreateAttValues() const
174{
175 auto values = new std::vector<G4AttValue>;
176
177 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
178
179 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
180
181 values->push_back(G4AttValue("PN", ParticleName, ""));
182
183 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), ""));
184
185 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), ""));
186
187 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), ""));
188
189 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), ""));
190
191 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), ""));
192
193 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
194
195#ifdef G4ATTDEBUG
196 G4cout << G4AttCheck(values, GetAttDefs());
197#endif
198
199 return values;
200}
201
203{
204 positionRecord->push_back(new G4SmoothTrajectoryPoint(
206}
207
212
214{
215 if (secondTrajectory == nullptr) return;
216
217 auto seco = (G4SmoothTrajectory*)secondTrajectory;
218 G4int ent = seco->GetPointEntries();
219
220 // initial point of the second trajectory should not be merged
221 //
222 for (G4int i = 1; i < ent; ++i) {
223 positionRecord->push_back((*(seco->positionRecord))[i]);
224 }
225 delete (*seco->positionRecord)[0];
226 seco->positionRecord->clear();
227}
G4TemplateAutoLock< G4Mutex > G4AutoLock
G4Allocator< G4SmoothTrajectory > *& aSmoothTrajectoryAllocator()
#define G4BestUnit(a, b)
#define G4MUTEX_INITIALIZER
std::mutex G4Mutex
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
const G4String & GetParticleName() const
static G4ParticleTable * GetParticleTable()
friend class G4ClonedSmoothTrajectory
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
G4SmoothTrajectory()=default
void AppendStep(const G4Step *aStep) override
G4int GetPointEntries() const override
std::vector< G4AttValue > * CreateAttValues() const override
G4VTrajectory * CloneForMaster() const override
void DrawTrajectory() const override
G4ParticleDefinition * GetParticleDefinition()
void ShowTrajectory(std::ostream &os=G4cout) const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
const G4ThreeVector & GetPosition() const
std::vector< G4ThreeVector > * GetPointerToVectorOfAuxiliaryPoints() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetParentID() 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)
#define G4ThreadLocalStatic
Definition tls.hh:76