Loading...
Searching...
No Matches
RE04Trajectory.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/// \file runAndEvent/RE04/src/RE04Trajectory.cc
27/// \brief Implementation of the RE04Trajectory class
28//
29//
30#include "RE04Trajectory.hh"
32#include "G4ParticleTable.hh"
33#include "G4AttDefStore.hh"
34#include "G4AttDef.hh"
35#include "G4AttValue.hh"
36#include "G4UIcommand.hh"
37#include "G4UnitsTable.hh"
38
39//#define G4ATTDEBUG
40#ifdef G4ATTDEBUG
41#include "G4AttCheck.hh"
42#endif
43
45
46//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 fPositionRecord(0), fTrackID(0), fParentID(0),
50 fPDGEncoding( 0 ), fPDGCharge(0.0), fParticleName(""),
51 fInitialKineticEnergy( 0. ), fInitialMomentum( G4ThreeVector() )
52{;}
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56{
57 G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
58 fParticleName = fpParticleDefinition->GetParticleName();
59 fPDGCharge = fpParticleDefinition->GetPDGCharge();
60 fPDGEncoding = fpParticleDefinition->GetPDGEncoding();
61 fTrackID = aTrack->GetTrackID();
62 fParentID = aTrack->GetParentID();
63 fInitialKineticEnergy = aTrack->GetKineticEnergy();
64 fInitialMomentum = aTrack->GetMomentum();
66 // Following is for the first trajectory point
68 aTrack->GetPosition(),aTrack->GetMaterial()));
69}
70
71//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
73{
74 fParticleName = right.fParticleName;
75 fPDGCharge = right.fPDGCharge;
76 fPDGEncoding = right.fPDGEncoding;
77 fTrackID = right.fTrackID;
78 fParentID = right.fParentID;
79 fInitialKineticEnergy = right.fInitialKineticEnergy;
80 fInitialMomentum = right.fInitialMomentum;
82
83 for(size_t i=0;i<right.fPositionRecord->size();i++)
84 {
85 RE04TrajectoryPoint* rightPoint
86 = (RE04TrajectoryPoint*)((*(right.fPositionRecord))[i]);
87 fPositionRecord->push_back(new RE04TrajectoryPoint(*rightPoint));
88 }
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93{
94 if (fPositionRecord) {
95 // fPositionRecord->clearAndDestroy();
96 size_t i;
97 for(i=0;i<fPositionRecord->size();i++){
98 delete (*fPositionRecord)[i];
99 }
100 fPositionRecord->clear();
101 delete fPositionRecord;
102 }
103}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106void RE04Trajectory::ShowTrajectory(std::ostream& os) const
107{
108 // Invoke the default implementation in G4VTrajectory...
109 G4VTrajectory::ShowTrajectory(os);
110 // ... or override with your own code here.
111}
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115{
116 // Invoke the default implementation in G4VTrajectory...
117 G4VTrajectory::DrawTrajectory();
118 // ... or override with your own code here.
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122const std::map<G4String,G4AttDef>* RE04Trajectory::GetAttDefs() const
123{
124 G4bool isNew;
125 std::map<G4String,G4AttDef>* store
126 = G4AttDefStore::GetInstance("RE04Trajectory",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//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
166std::vector<G4AttValue>* RE04Trajectory::CreateAttValues() const
167{
168 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
169
170 values->push_back
171 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
172
173 values->push_back
174 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
175
176 values->push_back(G4AttValue("PN",fParticleName,""));
177
178 values->push_back
179 (G4AttValue("Ch",G4UIcommand::ConvertToString(fPDGCharge),""));
180
181 values->push_back
182 (G4AttValue("PDG",G4UIcommand::ConvertToString(fPDGEncoding),""));
183
184 values->push_back
185 (G4AttValue("IKE",G4BestUnit(fInitialKineticEnergy,"Energy"),""));
186
187 values->push_back
188 (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
189
190 values->push_back
191 (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
192
193 values->push_back
194 (G4AttValue("NTP",G4UIcommand::ConvertToString(GetPointEntries()),""));
195
196#ifdef G4ATTDEBUG
197 G4cout << G4AttCheck(values,GetAttDefs());
198#endif
199
200 return values;
201}
202
203//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
205{
206 fPositionRecord->push_back( new RE04TrajectoryPoint(
207 aStep->GetPostStepPoint()->GetPosition(),
208 aStep->GetPreStepPoint()->GetMaterial() ));
209}
210
211//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
213{
214 return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
219{
220 if(!secondTrajectory) return;
221
222 RE04Trajectory* seco = (RE04Trajectory*)secondTrajectory;
223 G4int ent = seco->GetPointEntries();
224 for(G4int i=1;i<ent;i++) // initial point of the second trajectory
225 // should not be merged
226 {
227 fPositionRecord->push_back((*(seco->fPositionRecord))[i]);
228 // fPositionRecord->push_back(seco->fPositionRecord->removeAt(1));
229 }
230 delete (*seco->fPositionRecord)[0];
231 seco->fPositionRecord->clear();
232}
233
234
Definition of the RE04TrajectoryPoint class.
G4ThreadLocal G4Allocator< RE04Trajectory > * faTrajAllocator
Definition of the RE04Trajectory class.
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
Trajectory point class.
User trajectory class.
virtual void AppendStep(const G4Step *aStep)
G4double fInitialKineticEnergy
TrajectoryPointContainer * fPositionRecord
virtual int GetPointEntries() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
G4ParticleDefinition * GetParticleDefinition()
virtual ~RE04Trajectory()
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual void DrawTrajectory() const
G4ThreeVector fInitialMomentum
virtual void ShowTrajectory(std::ostream &os=G4cout) const

Applications | User Support | Publications | Collaboration