Geant4  10.00.p01
F04Trajectory.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 // $Id: F04Trajectory.cc 76690 2013-11-14 08:45:07Z gcosmo $
27 //
30 //
31 
32 #include "G4AttDef.hh"
33 #include "G4AttValue.hh"
34 #include "G4AttDefStore.hh"
35 
36 #include "G4UIcommand.hh"
37 #include "G4UnitsTable.hh"
38 
39 #include "F04Trajectory.hh"
40 #include "F04TrajectoryPoint.hh"
41 #include "G4ParticleTable.hh"
42 
43 //#define G4ATTDEBUG
44 #ifdef G4ATTDEBUG
45 #include "G4AttCheck.hh"
46 #endif
47 
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
53  : fpPointsContainer(0), fTrackID(0), fParentID(0),
54  fPDGCharge(0.0), fPDGEncoding(0), fParticleName(""),
55  fInitialMomentum(G4ThreeVector()) {;}
56 
57 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58 
60 {
61  G4ParticleDefinition * particleDefinition = aTrack->GetDefinition();
62  fParticleName = particleDefinition->GetParticleName();
63  fPDGCharge = particleDefinition->GetPDGCharge();
64  fPDGEncoding = particleDefinition->GetPDGEncoding();
65  fTrackID = aTrack->GetTrackID();
66  fParentID = aTrack->GetParentID();
67  fInitialMomentum = aTrack->GetMomentum();
69  // Following is for the first trajectory point
70  fpPointsContainer->push_back(new F04TrajectoryPoint(aTrack));
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74 
76 {
78  fPDGCharge = right.fPDGCharge;
79  fPDGEncoding = right.fPDGEncoding;
80  fTrackID = right.fTrackID;
81  fParentID = right.fParentID;
84 
85  for(size_t i=0;i<right.fpPointsContainer->size();++i) {
86  F04TrajectoryPoint* rightPoint
87  = (F04TrajectoryPoint*)((*(right.fpPointsContainer))[i]);
88  fpPointsContainer->push_back(new F04TrajectoryPoint(*rightPoint));
89  }
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
95 {
96  for(size_t i=0;i<fpPointsContainer->size();++i){
97  delete (*fpPointsContainer)[i];
98  }
99  fpPointsContainer->clear();
100 
101  delete fpPointsContainer;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
106 void F04Trajectory::ShowTrajectory(std::ostream& os) const
107 {
108  // Invoke the default implementation in G4VTrajectory...
110  // ... or override with your own code here.
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
114 
116 {
117  // Invoke the default implementation in G4VTrajectory...
119  // ... or override with your own code here.
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123 
125 {
126  fpPointsContainer->push_back(new F04TrajectoryPoint(aStep));
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
132 {
133  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137 
139 {
140  if(!secondTrajectory) return;
141 
142  F04Trajectory* second = (F04Trajectory*)secondTrajectory;
143  G4int ent = second->GetPointEntries();
144  // initial point of the second trajectory should not be merged
145  for(G4int i=1; i<ent; ++i) {
146  fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
147  }
148  delete (*second->fpPointsContainer)[0];
149  second->fpPointsContainer->clear();
150 }
151 
152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
153 
154 const std::map<G4String,G4AttDef>* F04Trajectory::GetAttDefs() const
155 {
156  G4bool isNew;
157  std::map<G4String,G4AttDef>* store
158  = G4AttDefStore::GetInstance("Trajectory",isNew);
159 
160  if (isNew) {
161 
162  G4String ID("ID");
163  (*store)[ID] = G4AttDef(ID,"Track ID","Bookkeeping","","G4int");
164 
165  G4String PID("PID");
166  (*store)[PID] = G4AttDef(PID,"Parent ID","Bookkeeping","","G4int");
167 
168  G4String PN("PN");
169  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
170 
171  G4String Ch("Ch");
172  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
173 
174  G4String PDG("PDG");
175  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
176 
177  G4String IMom("IMom");
178  (*store)[IMom] = G4AttDef(IMom,
179  "Momentum of track at start of trajectory",
180  "Physics","G4BestUnit","G4ThreeVector");
181 
182  G4String IMag("IMag");
183  (*store)[IMag] = G4AttDef(IMag,
184  "Magnitude of momentum of track at start of trajectory",
185  "Physics","G4BestUnit","G4double");
186 
187  G4String NTP("NTP");
188  (*store)[NTP] = G4AttDef(NTP,"No. of points","Bookkeeping","","G4int");
189 
190  }
191  return store;
192 }
193 
194 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
195 
196 std::vector<G4AttValue>* F04Trajectory::CreateAttValues() const
197 {
198  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
199 
200  values->push_back
202 
203  values->push_back
205 
206  values->push_back(G4AttValue("PN",fParticleName,""));
207 
208  values->push_back
210 
211  values->push_back
213 
214  values->push_back
215  (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
216 
217  values->push_back
218  (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
219 
220  values->push_back
222 
223 #ifdef G4ATTDEBUG
224  G4cout << G4AttCheck(values,GetAttDefs());
225 #endif
226  return values;
227 }
TrajectoryPointContainer * fpPointsContainer
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
virtual ~F04Trajectory()
Definition of the F04Trajectory class.
CLHEP::Hep3Vector G4ThreeVector
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
virtual int GetPointEntries() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
G4ThreeVector fInitialMomentum
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
Definition of the F04TrajectoryPoint class.
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4ThreadLocal G4Allocator< F04Trajectory > * F04TrajectoryAllocator
G4String fParticleName
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual void DrawTrajectory() const
G4double fPDGCharge
G4GLOB_DLL std::ostream G4cout
virtual void ShowTrajectory(std::ostream &os=G4cout) const
bool G4bool
Definition: G4Types.hh:79
static const double second
Definition: G4SIunits.hh:138
Definition: G4Step.hh:76
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4int GetTrackID() const
G4ThreeVector GetMomentum() const
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * GetParticleDefinition()
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4double GetPDGCharge() const
virtual void DrawTrajectory() const
virtual void AppendStep(const G4Step *aStep)