Geant4 Cross Reference |
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 /// \file optical/wls/src/WLSTrajectory.cc 28 /// \brief Implementation of the WLSTrajectory class 29 // 30 // 31 32 #include "WLSTrajectory.hh" 33 34 #include "WLSTrajectoryPoint.hh" 35 36 #include "G4AttDef.hh" 37 #include "G4AttDefStore.hh" 38 #include "G4AttValue.hh" 39 #include "G4Colour.hh" 40 #include "G4ParticleTable.hh" 41 #include "G4ParticleTypes.hh" 42 #include "G4Polyline.hh" 43 #include "G4Polymarker.hh" 44 #include "G4UIcommand.hh" 45 #include "G4UnitsTable.hh" 46 #include "G4VVisManager.hh" 47 #include "G4VisAttributes.hh" 48 49 G4ThreadLocal G4Allocator<WLSTrajectory>* WLSTrajectoryAllocator = nullptr; 50 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 52 53 WLSTrajectory::WLSTrajectory(const G4Track* aTrack) 54 { 55 fParticleDefinition = aTrack->GetDefinition(); 56 fParticleName = fParticleDefinition->GetParticleName(); 57 fPDGCharge = fParticleDefinition->GetPDGCharge(); 58 fPDGEncoding = fParticleDefinition->GetPDGEncoding(); 59 fTrackID = aTrack->GetTrackID(); 60 fParentID = aTrack->GetParentID(); 61 fInitialMomentum = aTrack->GetMomentum(); 62 fpPointsContainer = new WLSTrajectoryPointContainer(); 63 // Following is for the first trajectory point 64 fpPointsContainer->push_back(new WLSTrajectoryPoint(aTrack)); 65 } 66 67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 68 69 WLSTrajectory::WLSTrajectory(WLSTrajectory& right) : G4VTrajectory() 70 { 71 fParticleDefinition = right.fParticleDefinition; 72 fParticleName = right.fParticleName; 73 fPDGCharge = right.fPDGCharge; 74 fPDGEncoding = right.fPDGEncoding; 75 fTrackID = right.fTrackID; 76 fParentID = right.fParentID; 77 fInitialMomentum = right.fInitialMomentum; 78 fpPointsContainer = new WLSTrajectoryPointContainer(); 79 80 for (size_t i = 0; i < right.fpPointsContainer->size(); ++i) { 81 auto rightPoint = (WLSTrajectoryPoint*)((*(right.fpPointsContainer))[i]); 82 fpPointsContainer->push_back(new WLSTrajectoryPoint(*rightPoint)); 83 } 84 } 85 86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 87 88 WLSTrajectory::~WLSTrajectory() 89 { 90 for (size_t i = 0; i < fpPointsContainer->size(); ++i) { 91 delete (*fpPointsContainer)[i]; 92 } 93 fpPointsContainer->clear(); 94 delete fpPointsContainer; 95 } 96 97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 98 99 void WLSTrajectory::ShowTrajectory(std::ostream& os) const 100 { 101 // Invoke the default implementation in G4VTrajectory... 102 G4VTrajectory::ShowTrajectory(os); 103 // ... or override with your own code here. 104 } 105 106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 107 108 void WLSTrajectory::AppendStep(const G4Step* aStep) 109 { 110 fpPointsContainer->push_back(new WLSTrajectoryPoint(aStep)); 111 } 112 113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 114 115 G4ParticleDefinition* WLSTrajectory::GetParticleDefinition() 116 { 117 return (G4ParticleTable::GetParticleTable()->FindParticle(fParticleName)); 118 } 119 120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 121 122 void WLSTrajectory::MergeTrajectory(G4VTrajectory* secondTrajectory) 123 { 124 if (!secondTrajectory) return; 125 126 auto second = (WLSTrajectory*)secondTrajectory; 127 G4int ent = second->GetPointEntries(); 128 // initial point of the second trajectory should not be merged 129 for (G4int i = 1; i < ent; ++i) { 130 fpPointsContainer->push_back((*(second->fpPointsContainer))[i]); 131 } 132 delete (*second->fpPointsContainer)[0]; 133 second->fpPointsContainer->clear(); 134 } 135 136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 137 138 const std::map<G4String, G4AttDef>* WLSTrajectory::GetAttDefs() const 139 { 140 G4bool isNew; 141 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("Trajectory", isNew); 142 143 if (isNew) { 144 G4String ID("ID"); 145 (*store)[ID] = G4AttDef(ID, "Track ID", "Bookkeeping", "", "G4int"); 146 147 G4String PID("PID"); 148 (*store)[PID] = G4AttDef(PID, "Parent ID", "Bookkeeping", "", "G4int"); 149 150 G4String PN("PN"); 151 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String"); 152 153 G4String Ch("Ch"); 154 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double"); 155 156 G4String PDG("PDG"); 157 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int"); 158 159 G4String IMom("IMom"); 160 (*store)[IMom] = G4AttDef(IMom, "Momentum of track at start of trajectory", "Physics", 161 "G4BestUnit", "G4ThreeVector"); 162 163 G4String IMag("IMag"); 164 (*store)[IMag] = G4AttDef(IMag, "Magnitude of momentum of track at start of trajectory", 165 "Physics", "G4BestUnit", "G4double"); 166 167 G4String NTP("NTP"); 168 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Bookkeeping", "", "G4int"); 169 } 170 return store; 171 } 172 173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 174 175 std::vector<G4AttValue>* WLSTrajectory::CreateAttValues() const 176 { 177 auto values = new std::vector<G4AttValue>; 178 179 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), "")); 180 181 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), "")); 182 183 values->push_back(G4AttValue("PN", fParticleName, "")); 184 185 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(fPDGCharge), "")); 186 187 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(fPDGEncoding), "")); 188 189 values->push_back(G4AttValue("IMom", G4BestUnit(fInitialMomentum, "Energy"), "")); 190 191 values->push_back(G4AttValue("IMag", G4BestUnit(fInitialMomentum.mag(), "Energy"), "")); 192 193 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), "")); 194 195 return values; 196 } 197