File indexing completed on 2026-04-02 07:37:25
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #include "WLSTrajectory.hh"
0030
0031 #include "WLSTrajectoryPoint.hh"
0032
0033 #include "G4AttDef.hh"
0034 #include "G4AttDefStore.hh"
0035 #include "G4AttValue.hh"
0036 #include "G4Colour.hh"
0037 #include "G4ParticleTable.hh"
0038 #include "G4ParticleTypes.hh"
0039 #include "G4Polyline.hh"
0040 #include "G4Polymarker.hh"
0041 #include "G4UIcommand.hh"
0042 #include "G4UnitsTable.hh"
0043 #include "G4VVisManager.hh"
0044 #include "G4VisAttributes.hh"
0045
0046 G4ThreadLocal G4Allocator<WLSTrajectory>* WLSTrajectoryAllocator = nullptr;
0047
0048
0049
0050 WLSTrajectory::WLSTrajectory(const G4Track* aTrack)
0051 {
0052 fParticleDefinition = aTrack->GetDefinition();
0053 fParticleName = fParticleDefinition->GetParticleName();
0054 fPDGCharge = fParticleDefinition->GetPDGCharge();
0055 fPDGEncoding = fParticleDefinition->GetPDGEncoding();
0056 fTrackID = aTrack->GetTrackID();
0057 fParentID = aTrack->GetParentID();
0058 fInitialMomentum = aTrack->GetMomentum();
0059 fpPointsContainer = new WLSTrajectoryPointContainer();
0060
0061 fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack));
0062 }
0063
0064
0065
0066 WLSTrajectory::WLSTrajectory(WLSTrajectory& right) : G4VTrajectory()
0067 {
0068 fParticleDefinition = right.fParticleDefinition;
0069 fParticleName = right.fParticleName;
0070 fPDGCharge = right.fPDGCharge;
0071 fPDGEncoding = right.fPDGEncoding;
0072 fTrackID = right.fTrackID;
0073 fParentID = right.fParentID;
0074 fInitialMomentum = right.fInitialMomentum;
0075 fpPointsContainer = new WLSTrajectoryPointContainer();
0076
0077 for (size_t i = 0; i < right.fpPointsContainer->size(); ++i) {
0078 auto rightPoint = (WLSTrajectoryPoint*)((*(right.fpPointsContainer))[i]);
0079 fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint));
0080 }
0081 }
0082
0083
0084
0085 WLSTrajectory::~WLSTrajectory()
0086 {
0087 for (size_t i = 0; i < fpPointsContainer->size(); ++i) {
0088 delete (*fpPointsContainer)[i];
0089 }
0090 fpPointsContainer->clear();
0091 delete fpPointsContainer;
0092 }
0093
0094
0095
0096 void WLSTrajectory::ShowTrajectory(std::ostream& os) const
0097 {
0098
0099 G4VTrajectory::ShowTrajectory(os);
0100
0101 }
0102
0103
0104
0105 void WLSTrajectory::AppendStep(const G4Step* aStep)
0106 {
0107 fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep));
0108 }
0109
0110
0111
0112 G4ParticleDefinition* WLSTrajectory::GetParticleDefinition()
0113 {
0114 return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName));
0115 }
0116
0117
0118
0119 void WLSTrajectory::MergeTrajectory(G4VTrajectory* secondTrajectory)
0120 {
0121 if (!secondTrajectory) return;
0122
0123 auto second = (WLSTrajectory*)secondTrajectory;
0124 G4int ent = second->GetPointEntries();
0125
0126 for (G4int i = 1; i < ent; ++i) {
0127 fpPointsContainer->push_back((*(second->fpPointsContainer))[i]);
0128 }
0129 delete (*second->fpPointsContainer)[0];
0130 second->fpPointsContainer->clear();
0131 }
0132
0133
0134
0135 const std::map<G4String, G4AttDef>* WLSTrajectory::GetAttDefs() const
0136 {
0137 G4bool isNew;
0138 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("Trajectory", isNew);
0139
0140 if (isNew) {
0141 G4String ID("ID");
0142 (*store)[ID] = G4AttDef(ID, "Track ID", "Bookkeeping", "", "G4int");
0143
0144 G4String PID("PID");
0145 (*store)[PID] = G4AttDef(PID, "Parent ID", "Bookkeeping", "", "G4int");
0146
0147 G4String PN("PN");
0148 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String");
0149
0150 G4String Ch("Ch");
0151 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double");
0152
0153 G4String PDG("PDG");
0154 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int");
0155
0156 G4String IMom("IMom");
0157 (*store)[IMom] = G4AttDef(IMom, "Momentum of track at start of trajectory", "Physics",
0158 "G4BestUnit", "G4ThreeVector");
0159
0160 G4String IMag("IMag");
0161 (*store)[IMag] = G4AttDef(IMag, "Magnitude of momentum of track at start of trajectory",
0162 "Physics", "G4BestUnit", "G4double");
0163
0164 G4String NTP("NTP");
0165 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Bookkeeping", "", "G4int");
0166 }
0167 return store;
0168 }
0169
0170
0171
0172 std::vector<G4AttValue>* WLSTrajectory::CreateAttValues() const
0173 {
0174 auto values = new std::vector<G4AttValue>;
0175
0176 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), ""));
0177
0178 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), ""));
0179
0180 values->push_back(G4AttValue("PN", fParticleName, ""));
0181
0182 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(fPDGCharge), ""));
0183
0184 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(fPDGEncoding), ""));
0185
0186 values->push_back(G4AttValue("IMom", G4BestUnit(fInitialMomentum, "Energy"), ""));
0187
0188 values->push_back(G4AttValue("IMag", G4BestUnit(fInitialMomentum.mag(), "Energy"), ""));
0189
0190 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), ""));
0191
0192 return values;
0193 }