Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Trajectory.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//
27// $Id$
28//
29// ---------------------------------------------------------------
30//
31// G4Trajectory.cc
32//
33// Contact:
34// Questions and comments to this code should be sent to
35// Katsuya Amako (e-mail: [email protected])
36// Makoto Asai (e-mail: [email protected])
37// Takashi Sasaki (e-mail: [email protected])
38//
39// ---------------------------------------------------------------
40
41#include "G4Trajectory.hh"
42#include "G4TrajectoryPoint.hh"
43#include "G4ParticleTable.hh"
44#include "G4AttDefStore.hh"
45#include "G4AttDef.hh"
46#include "G4AttValue.hh"
47#include "G4UIcommand.hh"
48#include "G4UnitsTable.hh"
49
50//#define G4ATTDEBUG
51#ifdef G4ATTDEBUG
52#include "G4AttCheck.hh"
53#endif
54
56
58: positionRecord(0), fTrackID(0), fParentID(0),
59 PDGEncoding( 0 ), PDGCharge(0.0), ParticleName(""),
60 initialKineticEnergy( 0. ), initialMomentum( G4ThreeVector() )
61{;}
62
64{
65 G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
66 ParticleName = fpParticleDefinition->GetParticleName();
67 PDGCharge = fpParticleDefinition->GetPDGCharge();
68 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
69 fTrackID = aTrack->GetTrackID();
70 fParentID = aTrack->GetParentID();
71 initialKineticEnergy = aTrack->GetKineticEnergy();
72 initialMomentum = aTrack->GetMomentum();
73 positionRecord = new TrajectoryPointContainer();
74 // Following is for the first trajectory point
75 positionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
76}
77
79{
80 ParticleName = right.ParticleName;
81 PDGCharge = right.PDGCharge;
82 PDGEncoding = right.PDGEncoding;
83 fTrackID = right.fTrackID;
84 fParentID = right.fParentID;
85 initialKineticEnergy = right.initialKineticEnergy;
86 initialMomentum = right.initialMomentum;
87 positionRecord = new TrajectoryPointContainer();
88
89 for(size_t i=0;i<right.positionRecord->size();i++)
90 {
91 G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right.positionRecord))[i]);
92 positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
93 }
94}
95
97{
98 if (positionRecord) {
99 // positionRecord->clearAndDestroy();
100 size_t i;
101 for(i=0;i<positionRecord->size();i++){
102 delete (*positionRecord)[i];
103 }
104 positionRecord->clear();
105 delete positionRecord;
106 }
107}
108
109void G4Trajectory::ShowTrajectory(std::ostream& os) const
110{
111 // Invoke the default implementation in G4VTrajectory...
113 // ... or override with your own code here.
114}
115
116/***
117void G4Trajectory::DrawTrajectory() const
118{
119 // Invoke the default implementation in G4VTrajectory...
120 G4VTrajectory::DrawTrajectory();
121 // ... or override with your own code here.
122}
123***/
124
126{
127 // Invoke the default implementation in G4VTrajectory...
129 // ... or override with your own code here.
130}
131
132const std::map<G4String,G4AttDef>* G4Trajectory::GetAttDefs() const
133{
134 G4bool isNew;
135 std::map<G4String,G4AttDef>* store
136 = G4AttDefStore::GetInstance("G4Trajectory",isNew);
137 if (isNew) {
138
139 G4String ID("ID");
140 (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int");
141
142 G4String PID("PID");
143 (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int");
144
145 G4String PN("PN");
146 (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
147
148 G4String Ch("Ch");
149 (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
150
151 G4String PDG("PDG");
152 (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
153
154 G4String IKE("IKE");
155 (*store)[IKE] =
156 G4AttDef(IKE, "Initial kinetic energy",
157 "Physics","G4BestUnit","G4double");
158
159 G4String IMom("IMom");
160 (*store)[IMom] = G4AttDef(IMom, "Initial momentum",
161 "Physics","G4BestUnit","G4ThreeVector");
162
163 G4String IMag("IMag");
164 (*store)[IMag] =
165 G4AttDef(IMag, "Initial momentum magnitude",
166 "Physics","G4BestUnit","G4double");
167
168 G4String NTP("NTP");
169 (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
170
171 }
172 return store;
173}
174
175std::vector<G4AttValue>* G4Trajectory::CreateAttValues() const
176{
177 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
178
179 values->push_back
180 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
181
182 values->push_back
183 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
184
185 values->push_back(G4AttValue("PN",ParticleName,""));
186
187 values->push_back
188 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
189
190 values->push_back
191 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
192
193 values->push_back
194 (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
195
196 values->push_back
197 (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
198
199 values->push_back
200 (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
201
202 values->push_back
204
205#ifdef G4ATTDEBUG
206 G4cout << G4AttCheck(values,GetAttDefs());
207#endif
208
209 return values;
210}
211
213{
214 positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
215 GetPosition() ));
216}
217
219{
220 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
221}
222
224{
225 if(!secondTrajectory) return;
226
227 G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
228 G4int ent = seco->GetPointEntries();
229 for(G4int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
230 {
231 positionRecord->push_back((*(seco->positionRecord))[i]);
232 // positionRecord->push_back(seco->positionRecord->removeAt(1));
233 }
234 delete (*seco->positionRecord)[0];
235 seco->positionRecord->clear();
236}
237
238
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4Allocator< G4Trajectory > aTrajectoryAllocator
Definition: G4Trajectory.cc:55
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
Definition: G4Trajectory.hh:67
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4DLLIMPORT std::ostream G4cout
double mag() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4ParticleTable * GetParticleTable()
Definition: G4Step.hh:78
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4int GetParentID() const
G4ParticleDefinition * GetParticleDefinition()
virtual void DrawTrajectory(G4int i_mode=0) const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual void AppendStep(const G4Step *aStep)
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual int GetPointEntries() const
virtual ~G4Trajectory()
Definition: G4Trajectory.cc:96
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:349
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual void DrawTrajectory(G4int i_mode=0) const
std::map< G4String, G4AttDef > * GetInstance(G4String storeKey, G4bool &isNew)