Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: $
30 //
31 #include "RE04Trajectory.hh"
32 #include "RE04TrajectoryPoint.hh"
33 #include "G4ParticleTable.hh"
34 #include "G4AttDefStore.hh"
35 #include "G4AttDef.hh"
36 #include "G4AttValue.hh"
37 #include "G4UIcommand.hh"
38 #include "G4UnitsTable.hh"
39 
40 //#define G4ATTDEBUG
41 #ifdef G4ATTDEBUG
42 #include "G4AttCheck.hh"
43 #endif
44 
46 
47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49 : G4VTrajectory(),
50  fPositionRecord(0), fTrackID(0), fParentID(0),
51  fPDGEncoding( 0 ), fPDGCharge(0.0), fParticleName(""),
52  fInitialKineticEnergy( 0. ), fInitialMomentum( G4ThreeVector() )
53 {;}
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
57 {
58  G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
59  fParticleName = fpParticleDefinition->GetParticleName();
60  fPDGCharge = fpParticleDefinition->GetPDGCharge();
61  fPDGEncoding = fpParticleDefinition->GetPDGEncoding();
62  fTrackID = aTrack->GetTrackID();
63  fParentID = aTrack->GetParentID();
64  fInitialKineticEnergy = aTrack->GetKineticEnergy();
65  fInitialMomentum = aTrack->GetMomentum();
66  fPositionRecord = new TrajectoryPointContainer();
67  // Following is for the first trajectory point
68  fPositionRecord->push_back(new RE04TrajectoryPoint(
69  aTrack->GetPosition(),aTrack->GetMaterial()));
70 }
71 
72 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
74 {
75  fParticleName = right.fParticleName;
76  fPDGCharge = right.fPDGCharge;
77  fPDGEncoding = right.fPDGEncoding;
78  fTrackID = right.fTrackID;
79  fParentID = right.fParentID;
80  fInitialKineticEnergy = right.fInitialKineticEnergy;
81  fInitialMomentum = right.fInitialMomentum;
82  fPositionRecord = new TrajectoryPointContainer();
83 
84  for(size_t i=0;i<right.fPositionRecord->size();i++)
85  {
86  RE04TrajectoryPoint* rightPoint
87  = (RE04TrajectoryPoint*)((*(right.fPositionRecord))[i]);
88  fPositionRecord->push_back(new RE04TrajectoryPoint(*rightPoint));
89  }
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 {
95  if (fPositionRecord) {
96  // fPositionRecord->clearAndDestroy();
97  size_t i;
98  for(i=0;i<fPositionRecord->size();i++){
99  delete (*fPositionRecord)[i];
100  }
101  fPositionRecord->clear();
102  delete fPositionRecord;
103  }
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 void RE04Trajectory::ShowTrajectory(std::ostream& os) const
108 {
109  // Invoke the default implementation in G4VTrajectory...
111  // ... or override with your own code here.
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 {
117  // Invoke the default implementation in G4VTrajectory...
119  // ... or override with your own code here.
120 }
121 
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123 const std::map<G4String,G4AttDef>* RE04Trajectory::GetAttDefs() const
124 {
125  G4bool isNew;
126  std::map<G4String,G4AttDef>* store
127  = G4AttDefStore::GetInstance("RE04Trajectory",isNew);
128  if (isNew) {
129 
130  G4String id("ID");
131  (*store)[id] = G4AttDef(id,"Track ID","Physics","","G4int");
132 
133  G4String pid("PID");
134  (*store)[pid] = G4AttDef(pid,"Parent ID","Physics","","G4int");
135 
136  G4String pn("PN");
137  (*store)[pn] = G4AttDef(pn,"Particle Name","Physics","","G4String");
138 
139  G4String ch("Ch");
140  (*store)[ch] = G4AttDef(ch,"Charge","Physics","e+","G4double");
141 
142  G4String pdg("PDG");
143  (*store)[pdg] = G4AttDef(pdg,"PDG Encoding","Physics","","G4int");
144 
145  G4String ike("IKE");
146  (*store)[ike] =
147  G4AttDef(ike, "Initial kinetic energy",
148  "Physics","G4BestUnit","G4double");
149 
150  G4String iMom("IMom");
151  (*store)[iMom] = G4AttDef(iMom, "Initial momentum",
152  "Physics","G4BestUnit","G4ThreeVector");
153 
154  G4String iMag("IMag");
155  (*store)[iMag] =
156  G4AttDef(iMag, "Initial momentum magnitude",
157  "Physics","G4BestUnit","G4double");
158 
159  G4String ntp("NTP");
160  (*store)[ntp] = G4AttDef(ntp,"No. of points","Physics","","G4int");
161 
162  }
163  return store;
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167 std::vector<G4AttValue>* RE04Trajectory::CreateAttValues() const
168 {
169  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
170 
171  values->push_back
172  (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
173 
174  values->push_back
175  (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
176 
177  values->push_back(G4AttValue("PN",fParticleName,""));
178 
179  values->push_back
180  (G4AttValue("Ch",G4UIcommand::ConvertToString(fPDGCharge),""));
181 
182  values->push_back
183  (G4AttValue("PDG",G4UIcommand::ConvertToString(fPDGEncoding),""));
184 
185  values->push_back
186  (G4AttValue("IKE",G4BestUnit(fInitialKineticEnergy,"Energy"),""));
187 
188  values->push_back
189  (G4AttValue("IMom",G4BestUnit(fInitialMomentum,"Energy"),""));
190 
191  values->push_back
192  (G4AttValue("IMag",G4BestUnit(fInitialMomentum.mag(),"Energy"),""));
193 
194  values->push_back
196 
197 #ifdef G4ATTDEBUG
198  G4cout << G4AttCheck(values,GetAttDefs());
199 #endif
200 
201  return values;
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206 {
207  fPositionRecord->push_back( new RE04TrajectoryPoint(
208  aStep->GetPostStepPoint()->GetPosition(),
209  aStep->GetPreStepPoint()->GetMaterial() ));
210 }
211 
212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
214 {
215  return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
216 }
217 
218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 {
221  if(!secondTrajectory) return;
222 
223  RE04Trajectory* seco = (RE04Trajectory*)secondTrajectory;
224  G4int ent = seco->GetPointEntries();
225  for(G4int i=1;i<ent;i++) // initial point of the second trajectory
226  // should not be merged
227  {
228  fPositionRecord->push_back((*(seco->fPositionRecord))[i]);
229  // fPositionRecord->push_back(seco->fPositionRecord->removeAt(1));
230  }
231  delete (*seco->fPositionRecord)[0];
232  seco->fPositionRecord->clear();
233 }
234 
235 
G4ParticleDefinition * GetDefinition() const
virtual ~RE04Trajectory()
G4int GetParentID() const
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:357
Definition of the RE04Trajectory class.
Definition of the RE04TrajectoryPoint class.
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4ThreadLocal G4Allocator< RE04Trajectory > * faTrajAllocator
virtual void DrawTrajectory() const
G4StepPoint * GetPreStepPoint() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetParticleDefinition()
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
virtual void DrawTrajectory() const
Definition: G4Step.hh:76
G4int GetTrackID() const
G4Material * GetMaterial() const
virtual int GetPointEntries() const
G4ThreeVector GetMomentum() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4ParticleTable * GetParticleTable()
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
G4StepPoint * GetPostStepPoint() const
virtual void AppendStep(const G4Step *aStep)
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4double GetPDGCharge() const
double mag() const
virtual std::vector< G4AttValue > * CreateAttValues() const