Geant4 11.1.1
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#include "G4ParticleTable.hh"
38#include "G4AttDefStore.hh"
39#include "G4AttDef.hh"
40#include "G4AttValue.hh"
41#include "G4UIcommand.hh"
42#include "G4UnitsTable.hh"
43
44//#define G4ATTDEBUG
45#ifdef G4ATTDEBUG
46#include "G4AttCheck.hh"
47#endif
48
50{
52 return _instance;
53}
54
56 : initialMomentum( G4ThreeVector() )
57{
58}
59
61{
62 G4ParticleDefinition* fpParticleDefinition = aTrack->GetDefinition();
63 ParticleName = fpParticleDefinition->GetParticleName();
64 PDGCharge = fpParticleDefinition->GetPDGCharge();
65 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
66 fTrackID = aTrack->GetTrackID();
67 fParentID = aTrack->GetParentID();
68 initialKineticEnergy = aTrack->GetKineticEnergy();
69 initialMomentum = aTrack->GetMomentum();
70 positionRecord = new G4TrajectoryPointContainer();
71
72 // Following is for the first trajectory point
73 //
74 positionRecord
75 ->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition()));
76
77 // The first point has no auxiliary points, so set the auxiliary
78 // points vector to NULL
79 //
80 positionRecord
81 ->push_back(new G4SmoothTrajectoryPoint(aTrack->GetPosition(), nullptr));
82}
83
86{
87 ParticleName = right.ParticleName;
88 PDGCharge = right.PDGCharge;
89 PDGEncoding = right.PDGEncoding;
90 fTrackID = right.fTrackID;
91 fParentID = right.fParentID;
92 initialKineticEnergy = right.initialKineticEnergy;
93 initialMomentum = right.initialMomentum;
94 positionRecord = new G4TrajectoryPointContainer();
95
96 for(std::size_t i=0; i<right.positionRecord->size(); ++i)
97 {
98 G4SmoothTrajectoryPoint* rightPoint
99 = (G4SmoothTrajectoryPoint*)((*(right.positionRecord))[i]);
100 positionRecord->push_back(new G4SmoothTrajectoryPoint(*rightPoint));
101 }
102}
103
105{
106 if (positionRecord)
107 {
108 for(std::size_t i=0; i<positionRecord->size(); ++i)
109 {
110 delete (*positionRecord)[i];
111 }
112 positionRecord->clear();
113 delete positionRecord;
114 }
115}
116
117void G4SmoothTrajectory::ShowTrajectory(std::ostream& os) const
118{
119 // Invoke the default implementation in G4VTrajectory...
120 //
122
123 // ... or override with your own code here.
124}
125
127{
128 // Invoke the default implementation in G4VTrajectory...
129 //
131
132 // ... or override with your own code here.
133}
134
135const std::map<G4String,G4AttDef>* G4SmoothTrajectory::GetAttDefs() const
136{
137 G4bool isNew;
138 std::map<G4String,G4AttDef>* store
139 = G4AttDefStore::GetInstance("G4SmoothTrajectory",isNew);
140 if (isNew)
141 {
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",
159 "Physics","G4BestUnit","G4double");
160
161 G4String IMom("IMom");
162 (*store)[IMom] = G4AttDef(IMom, "Initial momentum",
163 "Physics","G4BestUnit","G4ThreeVector");
164
165 G4String IMag("IMag");
166 (*store)[IMag] = G4AttDef(IMag, "Initial momentum magnitude",
167 "Physics","G4BestUnit","G4double");
168
169 G4String NTP("NTP");
170 (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
171 }
172 return store;
173}
174
175
176std::vector<G4AttValue>* G4SmoothTrajectory::CreateAttValues() const
177{
178 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
179
180 values->push_back
181 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
182
183 values->push_back
184 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
185
186 values->push_back(G4AttValue("PN",ParticleName,""));
187
188 values->push_back
189 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
190
191 values->push_back
192 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
193
194 values->push_back
195 (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
196
197 values->push_back
198 (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
199
200 values->push_back
201 (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
202
203 values->push_back
205
206#ifdef G4ATTDEBUG
207 G4cout << G4AttCheck(values,GetAttDefs());
208#endif
209
210 return values;
211}
212
214{
215 positionRecord->push_back(
218}
219
221{
222 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
223}
224
226{
227 if(secondTrajectory == nullptr) return;
228
229 G4SmoothTrajectory* seco = (G4SmoothTrajectory*)secondTrajectory;
230 G4int ent = seco->GetPointEntries();
231
232 // initial point of the second trajectory should not be merged
233 //
234 for(G4int i=1; i<ent; ++i)
235 {
236 positionRecord->push_back((*(seco->positionRecord))[i]);
237 }
238 delete (*seco->positionRecord)[0];
239 seco->positionRecord->clear();
240}
G4Allocator< G4SmoothTrajectory > *& aSmoothTrajectoryAllocator()
#define G4BestUnit(a, b)
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
double mag() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4ParticleTable * GetParticleTable()
virtual void DrawTrajectory() const
virtual G4int GetPointEntries() const
virtual void AppendStep(const G4Step *aStep)
G4ParticleDefinition * GetParticleDefinition()
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:62
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)
Definition: G4UIcommand.cc:446
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