00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041 #include "G4Trajectory.hh"
00042 #include "G4TrajectoryPoint.hh"
00043 #include "G4ParticleTable.hh"
00044 #include "G4AttDefStore.hh"
00045 #include "G4AttDef.hh"
00046 #include "G4AttValue.hh"
00047 #include "G4UIcommand.hh"
00048 #include "G4UnitsTable.hh"
00049
00050
00051 #ifdef G4ATTDEBUG
00052 #include "G4AttCheck.hh"
00053 #endif
00054
00055 G4Allocator<G4Trajectory> aTrajectoryAllocator;
00056
00057 G4Trajectory::G4Trajectory()
00058 : positionRecord(0), fTrackID(0), fParentID(0),
00059 PDGEncoding( 0 ), PDGCharge(0.0), ParticleName(""),
00060 initialKineticEnergy( 0. ), initialMomentum( G4ThreeVector() )
00061 {;}
00062
00063 G4Trajectory::G4Trajectory(const G4Track* aTrack)
00064 {
00065 G4ParticleDefinition * fpParticleDefinition = aTrack->GetDefinition();
00066 ParticleName = fpParticleDefinition->GetParticleName();
00067 PDGCharge = fpParticleDefinition->GetPDGCharge();
00068 PDGEncoding = fpParticleDefinition->GetPDGEncoding();
00069 fTrackID = aTrack->GetTrackID();
00070 fParentID = aTrack->GetParentID();
00071 initialKineticEnergy = aTrack->GetKineticEnergy();
00072 initialMomentum = aTrack->GetMomentum();
00073 positionRecord = new TrajectoryPointContainer();
00074
00075 positionRecord->push_back(new G4TrajectoryPoint(aTrack->GetPosition()));
00076 }
00077
00078 G4Trajectory::G4Trajectory(G4Trajectory & right):G4VTrajectory()
00079 {
00080 ParticleName = right.ParticleName;
00081 PDGCharge = right.PDGCharge;
00082 PDGEncoding = right.PDGEncoding;
00083 fTrackID = right.fTrackID;
00084 fParentID = right.fParentID;
00085 initialKineticEnergy = right.initialKineticEnergy;
00086 initialMomentum = right.initialMomentum;
00087 positionRecord = new TrajectoryPointContainer();
00088
00089 for(size_t i=0;i<right.positionRecord->size();i++)
00090 {
00091 G4TrajectoryPoint* rightPoint = (G4TrajectoryPoint*)((*(right.positionRecord))[i]);
00092 positionRecord->push_back(new G4TrajectoryPoint(*rightPoint));
00093 }
00094 }
00095
00096 G4Trajectory::~G4Trajectory()
00097 {
00098 if (positionRecord) {
00099
00100 size_t i;
00101 for(i=0;i<positionRecord->size();i++){
00102 delete (*positionRecord)[i];
00103 }
00104 positionRecord->clear();
00105 delete positionRecord;
00106 }
00107 }
00108
00109 void G4Trajectory::ShowTrajectory(std::ostream& os) const
00110 {
00111
00112 G4VTrajectory::ShowTrajectory(os);
00113
00114 }
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125 void G4Trajectory::DrawTrajectory(G4int i_mode) const
00126 {
00127
00128 G4VTrajectory::DrawTrajectory(i_mode);
00129
00130 }
00131
00132 const std::map<G4String,G4AttDef>* G4Trajectory::GetAttDefs() const
00133 {
00134 G4bool isNew;
00135 std::map<G4String,G4AttDef>* store
00136 = G4AttDefStore::GetInstance("G4Trajectory",isNew);
00137 if (isNew) {
00138
00139 G4String ID("ID");
00140 (*store)[ID] = G4AttDef(ID,"Track ID","Physics","","G4int");
00141
00142 G4String PID("PID");
00143 (*store)[PID] = G4AttDef(PID,"Parent ID","Physics","","G4int");
00144
00145 G4String PN("PN");
00146 (*store)[PN] = G4AttDef(PN,"Particle Name","Physics","","G4String");
00147
00148 G4String Ch("Ch");
00149 (*store)[Ch] = G4AttDef(Ch,"Charge","Physics","e+","G4double");
00150
00151 G4String PDG("PDG");
00152 (*store)[PDG] = G4AttDef(PDG,"PDG Encoding","Physics","","G4int");
00153
00154 G4String IKE("IKE");
00155 (*store)[IKE] =
00156 G4AttDef(IKE, "Initial kinetic energy",
00157 "Physics","G4BestUnit","G4double");
00158
00159 G4String IMom("IMom");
00160 (*store)[IMom] = G4AttDef(IMom, "Initial momentum",
00161 "Physics","G4BestUnit","G4ThreeVector");
00162
00163 G4String IMag("IMag");
00164 (*store)[IMag] =
00165 G4AttDef(IMag, "Initial momentum magnitude",
00166 "Physics","G4BestUnit","G4double");
00167
00168 G4String NTP("NTP");
00169 (*store)[NTP] = G4AttDef(NTP,"No. of points","Physics","","G4int");
00170
00171 }
00172 return store;
00173 }
00174
00175 std::vector<G4AttValue>* G4Trajectory::CreateAttValues() const
00176 {
00177 std::vector<G4AttValue>* values = new std::vector<G4AttValue>;
00178
00179 values->push_back
00180 (G4AttValue("ID",G4UIcommand::ConvertToString(fTrackID),""));
00181
00182 values->push_back
00183 (G4AttValue("PID",G4UIcommand::ConvertToString(fParentID),""));
00184
00185 values->push_back(G4AttValue("PN",ParticleName,""));
00186
00187 values->push_back
00188 (G4AttValue("Ch",G4UIcommand::ConvertToString(PDGCharge),""));
00189
00190 values->push_back
00191 (G4AttValue("PDG",G4UIcommand::ConvertToString(PDGEncoding),""));
00192
00193 values->push_back
00194 (G4AttValue("IKE",G4BestUnit(initialKineticEnergy,"Energy"),""));
00195
00196 values->push_back
00197 (G4AttValue("IMom",G4BestUnit(initialMomentum,"Energy"),""));
00198
00199 values->push_back
00200 (G4AttValue("IMag",G4BestUnit(initialMomentum.mag(),"Energy"),""));
00201
00202 values->push_back
00203 (G4AttValue("NTP",G4UIcommand::ConvertToString(GetPointEntries()),""));
00204
00205 #ifdef G4ATTDEBUG
00206 G4cout << G4AttCheck(values,GetAttDefs());
00207 #endif
00208
00209 return values;
00210 }
00211
00212 void G4Trajectory::AppendStep(const G4Step* aStep)
00213 {
00214 positionRecord->push_back( new G4TrajectoryPoint(aStep->GetPostStepPoint()->
00215 GetPosition() ));
00216 }
00217
00218 G4ParticleDefinition* G4Trajectory::GetParticleDefinition()
00219 {
00220 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName));
00221 }
00222
00223 void G4Trajectory::MergeTrajectory(G4VTrajectory* secondTrajectory)
00224 {
00225 if(!secondTrajectory) return;
00226
00227 G4Trajectory* seco = (G4Trajectory*)secondTrajectory;
00228 G4int ent = seco->GetPointEntries();
00229 for(G4int i=1;i<ent;i++)
00230 {
00231 positionRecord->push_back((*(seco->positionRecord))[i]);
00232
00233 }
00234 delete (*seco->positionRecord)[0];
00235 seco->positionRecord->clear();
00236 }
00237
00238