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

Applications | User Support | Publications | Collaboration