Geant4  10.01.p01
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: G4Trajectory.cc 69003 2013-04-15 09:25:23Z gcosmo $
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: Katsuya.Amako@kek.jp)
36 // Makoto Asai (e-mail: asai@kekvax.kek.jp)
37 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
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();
72  initialMomentum = aTrack->GetMomentum();
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;
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  size_t i;
100  for(i=0;i<positionRecord->size();i++){
101  delete (*positionRecord)[i];
102  }
103  positionRecord->clear();
104  delete positionRecord;
105  }
106 }
107 
108 void G4Trajectory::ShowTrajectory(std::ostream& os) const
109 {
110  // Invoke the default implementation in G4VTrajectory...
112  // ... or override with your own code here.
113 }
114 
116 {
117  // Invoke the default implementation in G4VTrajectory...
119  // ... or override with your own code here.
120 }
121 
122 const std::map<G4String,G4AttDef>* G4Trajectory::GetAttDefs() const
123 {
124  G4bool isNew;
125  std::map<G4String,G4AttDef>* store
126  = G4AttDefStore::GetInstance("G4Trajectory",isNew);
127  if (isNew) {
128 
129  G4String ID("ID");
130  (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int");
131 
132  G4String PID("PID");
133  (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int");
134 
135  G4String PN("PN");
136  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
137 
138  G4String Ch("Ch");
139  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
140 
141  G4String PDG("PDG");
142  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
143 
144  G4String IKE("IKE");
145  (*store)[IKE] =
146  G4AttDef(IKE, "Initial kinetic energy",
147  "Physics","G4BestUnit","G4double");
148 
149  G4String IMom("IMom");
150  (*store)[IMom] = G4AttDef(IMom, "Initial momentum",
151  "Physics","G4BestUnit","G4ThreeVector");
152 
153  G4String IMag("IMag");
154  (*store)[IMag] =
155  G4AttDef(IMag, "Initial momentum magnitude",
156  "Physics","G4BestUnit","G4double");
157 
158  G4String NTP("NTP");
159  (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
160 
161  }
162  return store;
163 }
164 
165 std::vector<G4AttValue>* G4Trajectory::CreateAttValues() const
166 {
167  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
168 
169  values->push_back
171 
172  values->push_back
174 
175  values->push_back(G4AttValue("PN",ParticleName,""));
176 
177  values->push_back
179 
180  values->push_back
182 
183  values->push_back
184  (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
185 
186  values->push_back
187  (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
188 
189  values->push_back
190  (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
191 
192  values->push_back
194 
195 #ifdef G4ATTDEBUG
196  G4cout << G4AttCheck(values,GetAttDefs());
197 #endif
198 
199  return values;
200 }
201 
203 {
204  positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
205  GetPosition() ));
206 }
207 
209 {
210  return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
211 }
212 
214 {
215  if(!secondTrajectory) return;
216 
217  G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
218  G4int ent = seco->GetPointEntries();
219  for(G4int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
220  {
221  positionRecord->push_back((*(seco->positionRecord))[i]);
222  }
223  delete (*seco->positionRecord)[0];
224  seco->positionRecord->clear();
225 }
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
TrajectoryPointContainer * positionRecord
CLHEP::Hep3Vector G4ThreeVector
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
const G4ThreeVector & GetPosition() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:371
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4double initialKineticEnergy
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
const G4String & GetParticleName() const
virtual void DrawTrajectory() const
virtual int GetPointEntries() const
G4double GetKineticEnergy() const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4GLOB_DLL std::ostream G4cout
G4String ParticleName
virtual ~G4Trajectory()
Definition: G4Trajectory.cc:96
G4ThreeVector initialMomentum
bool G4bool
Definition: G4Types.hh:79
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4Step.hh:76
G4int GetTrackID() const
G4double PDGCharge
G4ThreeVector GetMomentum() const
static G4ParticleTable * GetParticleTable()
virtual void DrawTrajectory() const
G4ThreadLocal G4Allocator< G4Trajectory > * aTrajectoryAllocator
Definition: G4Trajectory.cc:55
G4StepPoint * GetPostStepPoint() const
virtual void AppendStep(const G4Step *aStep)
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
G4ParticleDefinition * GetParticleDefinition()
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4double GetPDGCharge() const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)