source: trunk/examples/extended/field/field04/src/F04Trajectory.cc@ 1337

Last change on this file since 1337 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 6.9 KB
Line 
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//
28
29#include "G4AttDef.hh"
30#include "G4AttValue.hh"
31#include "G4AttDefStore.hh"
32
33#include "G4UIcommand.hh"
34#include "G4UnitsTable.hh"
35
36#include "F04Trajectory.hh"
37#include "F04TrajectoryPoint.hh"
38#include "G4ParticleTable.hh"
39
40//#define G4ATTDEBUG
41#ifdef G4ATTDEBUG
42#include "G4AttCheck.hh"
43#endif
44
45G4Allocator<F04Trajectory> aTrajectoryAllocator;
46
47F04Trajectory::F04Trajectory()
48 : fpPointsContainer(0), fTrackID(0), fParentID(0),
49 PDGCharge(0.0), PDGEncoding(0), ParticleName(""),
50 initialMomentum(G4ThreeVector()) {;}
51
52F04Trajectory::F04Trajectory(const G4Track* aTrack)
53{
54 G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
55 ParticleName = fpParticleDefinition->GetParticleName();
56 PDGCharge = fpParticleDefinition->GetPDGCharge();
57 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
58 fTrackID = aTrack->GetTrackID();
59 fParentID = aTrack->GetParentID();
60 initialMomentum = aTrack->GetMomentum();
61 fpPointsContainer = new TrajectoryPointContainer();
62 // Following is for the first trajectory point
63 fpPointsContainer->push_back(new F04TrajectoryPoint(aTrack));
64}
65
66F04Trajectory::F04Trajectory(F04Trajectory & right) : G4VTrajectory()
67{
68 ParticleName = right.ParticleName;
69 PDGCharge = right.PDGCharge;
70 PDGEncoding = right.PDGEncoding;
71 fTrackID = right.fTrackID;
72 fParentID = right.fParentID;
73 initialMomentum = right.initialMomentum;
74 fpPointsContainer = new TrajectoryPointContainer();
75
76 for(size_t i=0;i<right.fpPointsContainer->size();++i) {
77 F04TrajectoryPoint* rightPoint
78 = (F04TrajectoryPoint*)((*(right.fpPointsContainer))[i]);
79 fpPointsContainer->push_back(new F04TrajectoryPoint(*rightPoint));
80 }
81}
82
83F04Trajectory::~F04Trajectory()
84{
85 for(size_t i=0;i<fpPointsContainer->size();++i){
86 delete (*fpPointsContainer)[i];
87 }
88 fpPointsContainer->clear();
89
90 delete fpPointsContainer;
91}
92
93void F04Trajectory::ShowTrajectory(std::ostream& os) const
94{
95 // Invoke the default implementation in G4VTrajectory...
96 G4VTrajectory::ShowTrajectory(os);
97 // ... or override with your own code here.
98}
99
100void F04Trajectory::DrawTrajectory() const
101{
102 // Invoke the default implementation in G4VTrajectory...
103 G4VTrajectory::DrawTrajectory();
104 // ... or override with your own code here.
105}
106
107void F04Trajectory::DrawTrajectory(G4int i_mode) const
108{
109 // Invoke the default implementation in G4VTrajectory...
110 G4VTrajectory::DrawTrajectory(i_mode);
111 // ... or override with your own code here.
112}
113
114void F04Trajectory::AppendStep(const G4Step* aStep)
115{
116 fpPointsContainer->push_back(new F04TrajectoryPoint(aStep));
117}
118
119G4ParticleDefinition* F04Trajectory::GetParticleDefinition()
120{
121 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
122}
123
124void F04Trajectory::MergeTrajectory(G4VTrajectory* secondTrajectory)
125{
126 if(!secondTrajectory) return;
127
128 F04Trajectory* second = (F04Trajectory*)secondTrajectory;
129 G4int ent = second->GetPointEntries();
130 // initial point of the second trajectory should not be merged
131 for(G4int i=1; i<ent; ++i) {
132 fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
133 }
134 delete (*second->fpPointsContainer)[0];
135 second->fpPointsContainer->clear();
136}
137
138const std::map<G4String,G4AttDef>* F04Trajectory::GetAttDefs() const
139{
140 G4bool isNew;
141 std::map<G4String,G4AttDef>* store
142 = G4AttDefStore::GetInstance("Trajectory",isNew);
143
144 if (isNew) {
145
146 G4String ID("ID");
147 (*store)[ID] = G4AttDef(ID,"Track ID","Bookkeeping","","G4int");
148
149 G4String PID("PID");
150 (*store)[PID] = G4AttDef(PID,"Parent ID","Bookkeeping","","G4int");
151
152 G4String PN("PN");
153 (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
154
155 G4String Ch("Ch");
156 (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
157
158 G4String PDG("PDG");
159 (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
160
161 G4String IMom("IMom");
162 (*store)[IMom] = G4AttDef(IMom,
163 "Momentum of track at start of trajectory",
164 "Physics","G4BestUnit","G4ThreeVector");
165
166 G4String IMag("IMag");
167 (*store)[IMag] = G4AttDef(IMag,
168 "Magnitude of momentum of track at start of trajectory",
169 "Physics","G4BestUnit","G4double");
170
171 G4String NTP("NTP");
172 (*store)[NTP] = G4AttDef(NTP,"No. of points","Bookkeeping","","G4int");
173
174 }
175 return store;
176}
177
178std::vector<G4AttValue>* F04Trajectory::CreateAttValues() const
179{
180 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
181
182 values->push_back
183 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
184
185 values->push_back
186 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
187
188 values->push_back(G4AttValue("PN",ParticleName,""));
189
190 values->push_back
191 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
192
193 values->push_back
194 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
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
203 (G4AttValue("NTP",G4UIcommand::ConvertToString(GetPointEntries()),""));
204
205#ifdef G4ATTDEBUG
206 G4cout << G4AttCheck(values,GetAttDefs());
207#endif
208 return values;
209}
Note: See TracBrowser for help on using the repository browser.