Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Trajectory.cc
이 파일의 문서화 페이지로 가기
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 // $Id: G4Trajectory.cc 110262 2018-05-17 14:25:55Z gcosmo $
28 //
29 // ---------------------------------------------------------------
30 //
31 // G4Trajectory.cc
32 //
33 // Contact:
34 // Questions and comments to this code should be sent to
35 // Katsuya Amako (e-mail: Katsuya.Amako@kek.jp)
36 // Makoto Asai (e-mail: asai@kekvax.kek.jp)
37 // Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
38 //
39 // ---------------------------------------------------------------
40 
41 #include "G4Trajectory.hh"
42 #include "G4TrajectoryPoint.hh"
43 #include "G4ParticleTable.hh"
44 #include "G4AttDefStore.hh"
45 #include "G4AttDef.hh"
46 #include "G4AttValue.hh"
47 #include "G4UIcommand.hh"
48 #include "G4UnitsTable.hh"
49 
50 //#define G4ATTDEBUG
51 #ifdef G4ATTDEBUG
52 #include "G4AttCheck.hh"
53 #endif
54 
56 {
58  return _instance;
59 }
60 
62 : positionRecord(0), fTrackID(0), fParentID(0),
63  PDGEncoding( 0 ), PDGCharge(0.0), ParticleName(""),
64  initialKineticEnergy( 0. ), initialMomentum( G4ThreeVector() )
65 {;}
66 
68 {
69  G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
70  ParticleName = fpParticleDefinition->GetParticleName();
71  PDGCharge = fpParticleDefinition->GetPDGCharge();
72  PDGEncoding = fpParticleDefinition->GetPDGEncoding();
73  fTrackID = aTrack->GetTrackID();
74  fParentID = aTrack->GetParentID();
76  initialMomentum = aTrack->GetMomentum();
78  // Following is for the first trajectory point
79  positionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
80 }
81 
83 {
84  ParticleName = right.ParticleName;
85  PDGCharge = right.PDGCharge;
86  PDGEncoding = right.PDGEncoding;
87  fTrackID = right.fTrackID;
88  fParentID = right.fParentID;
92 
93  for(size_t i=0;i<right.positionRecord->size();i++)
94  {
95  G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right.positionRecord))[i]);
96  positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
97  }
98 }
99 
101 {
102  if (positionRecord) {
103  size_t i;
104  for(i=0;i<positionRecord->size();i++){
105  delete (*positionRecord)[i];
106  }
107  positionRecord->clear();
108  delete positionRecord;
109  }
110 }
111 
112 void G4Trajectory::ShowTrajectory(std::ostream& os) const
113 {
114  // Invoke the default implementation in G4VTrajectory...
116  // ... or override with your own code here.
117 }
118 
120 {
121  // Invoke the default implementation in G4VTrajectory...
123  // ... or override with your own code here.
124 }
125 
126 const std::map<G4String,G4AttDef>* G4Trajectory::GetAttDefs() const
127 {
128  G4bool isNew;
129  std::map<G4String,G4AttDef>* store
130  = G4AttDefStore::GetInstance("G4Trajectory",isNew);
131  if (isNew) {
132 
133  G4String ID("ID");
134  (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int");
135 
136  G4String PID("PID");
137  (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int");
138 
139  G4String PN("PN");
140  (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
141 
142  G4String Ch("Ch");
143  (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
144 
145  G4String PDG("PDG");
146  (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
147 
148  G4String IKE("IKE");
149  (*store)[IKE] =
150  G4AttDef(IKE, "Initial kinetic energy",
151  "Physics","G4BestUnit","G4double");
152 
153  G4String IMom("IMom");
154  (*store)[IMom] = G4AttDef(IMom, "Initial momentum",
155  "Physics","G4BestUnit","G4ThreeVector");
156 
157  G4String IMag("IMag");
158  (*store)[IMag] =
159  G4AttDef(IMag, "Initial momentum magnitude",
160  "Physics","G4BestUnit","G4double");
161 
162  G4String NTP("NTP");
163  (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
164 
165  }
166  return store;
167 }
168 
169 std::vector<G4AttValue>* G4Trajectory::CreateAttValues() const
170 {
171  std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
172 
173  values->push_back
175 
176  values->push_back
178 
179  values->push_back(G4AttValue("PN",ParticleName,""));
180 
181  values->push_back
183 
184  values->push_back
186 
187  values->push_back
188  (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
189 
190  values->push_back
191  (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
192 
193  values->push_back
194  (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
195 
196  values->push_back
198 
199 #ifdef G4ATTDEBUG
200  G4cout << G4AttCheck(values,GetAttDefs());
201 #endif
202 
203  return values;
204 }
205 
207 {
208  positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
209  GetPosition() ));
210 }
211 
213 {
214  return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
215 }
216 
218 {
219  if(!secondTrajectory) return;
220 
221  G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
222  G4int ent = seco->GetPointEntries();
223  for(G4int i=1;i<ent;i++) // initial point of the second trajectory should not be merged
224  {
225  positionRecord->push_back((*(seco->positionRecord))[i]);
226  }
227  delete (*seco->positionRecord)[0];
228  seco->positionRecord->clear();
229 }
virtual ~G4Trajectory()
G4double GetKineticEnergy() const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
static G4ParticleTable * GetParticleTable()
const G4String & GetParticleName() const
G4double initialKineticEnergy
G4int GetTrackID() const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4double GetPDGCharge() const
virtual int GetPointEntries() const
G4double PDGCharge
TrajectoryPointContainer * positionRecord
G4ParticleDefinition * GetDefinition() const
#define G4ThreadLocalStatic
Definition: tls.hh:68
virtual void DrawTrajectory() const
G4ThreeVector GetMomentum() const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetPosition() const
G4ThreeVector initialMomentum
Definition: G4Step.hh:76
G4StepPoint * GetPostStepPoint() const
virtual std::vector< G4AttValue > * CreateAttValues() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:374
std::map< G4String, G4AttDef > * GetInstance(const G4String &storeKey, G4bool &isNew)
double mag() const
int G4int
Definition: G4Types.hh:78
virtual void DrawTrajectory() const
G4ParticleDefinition * GetParticleDefinition()
G4GLOB_DLL std::ostream G4cout
G4TRACKING_DLL G4Allocator< G4Trajectory > *& aTrajectoryAllocator()
Definition: G4Trajectory.cc:55
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
G4String ParticleName
virtual void AppendStep(const G4Step *aStep)
G4int GetParentID() const