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 // G4ClonedTrajectory class implementation 27 // 28 // Makoto Asai (JLab) - Oct.2024 29 // -------------------------------------------------------------------- 30 31 #include "G4ClonedTrajectory.hh" 32 #include "G4Trajectory.hh" 33 34 #include "G4AttDef.hh" 35 #include "G4AttDefStore.hh" 36 #include "G4AttValue.hh" 37 #include "G4ParticleTable.hh" 38 #include "G4ClonedTrajectoryPoint.hh" 39 #include "G4UIcommand.hh" 40 #include "G4UnitsTable.hh" 41 #include "G4AutoLock.hh" 42 43 // #define G4ATTDEBUG 44 #ifdef G4ATTDEBUG 45 # include "G4AttCheck.hh" 46 #endif 47 48 G4Allocator<G4ClonedTrajectory>*& aClonedTrajectoryAllocator() 49 { 50 static G4Allocator<G4ClonedTrajectory>* _instance = nullptr; 51 return _instance; 52 } 53 54 G4ClonedTrajectory::G4ClonedTrajectory(const G4Trajectory& right) 55 { 56 ParticleName = right.ParticleName; 57 PDGCharge = right.PDGCharge; 58 PDGEncoding = right.PDGEncoding; 59 fTrackID = right.fTrackID; 60 fParentID = right.fParentID; 61 initialKineticEnergy = right.initialKineticEnergy; 62 initialMomentum = right.initialMomentum; 63 positionRecord = new G4ClonedTrajectoryPointContainer(); 64 65 for (auto& i : *right.positionRecord) { 66 auto rightPoint = (G4ClonedTrajectoryPoint*)i; 67 positionRecord->push_back(new G4ClonedTrajectoryPoint(*rightPoint)); 68 } 69 } 70 71 G4ClonedTrajectory::~G4ClonedTrajectory() 72 { 73 if (positionRecord != nullptr) { 74 for (auto& i : *positionRecord) { 75 delete i; 76 } 77 positionRecord->clear(); 78 delete positionRecord; 79 } 80 } 81 82 void G4ClonedTrajectory::ShowTrajectory(std::ostream& os) const 83 { 84 // Invoke the default implementation in G4VTrajectory... 85 G4VTrajectory::ShowTrajectory(os); 86 87 // ... or override with your own code here. 88 } 89 90 void G4ClonedTrajectory::DrawTrajectory() const 91 { 92 // Invoke the default implementation in G4VTrajectory... 93 G4VTrajectory::DrawTrajectory(); 94 95 // ... or override with your own code here. 96 } 97 98 const std::map<G4String, G4AttDef>* G4ClonedTrajectory::GetAttDefs() const 99 { 100 G4bool isNew; 101 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4Trajectory", isNew); 102 if (isNew) { 103 G4String ID("ID"); 104 (*store)[ID] = G4AttDef(ID, "Track ID", "Physics", "", "G4int"); 105 106 G4String PID("PID"); 107 (*store)[PID] = G4AttDef(PID, "Parent ID", "Physics", "", "G4int"); 108 109 G4String PN("PN"); 110 (*store)[PN] = G4AttDef(PN, "Particle Name", "Physics", "", "G4String"); 111 112 G4String Ch("Ch"); 113 (*store)[Ch] = G4AttDef(Ch, "Charge", "Physics", "e+", "G4double"); 114 115 G4String PDG("PDG"); 116 (*store)[PDG] = G4AttDef(PDG, "PDG Encoding", "Physics", "", "G4int"); 117 118 G4String IKE("IKE"); 119 (*store)[IKE] = G4AttDef(IKE, "Initial kinetic energy", "Physics", "G4BestUnit", "G4double"); 120 121 G4String IMom("IMom"); 122 (*store)[IMom] = G4AttDef(IMom, "Initial momentum", "Physics", "G4BestUnit", "G4ThreeVector"); 123 124 G4String IMag("IMag"); 125 (*store)[IMag] = 126 G4AttDef(IMag, "Initial momentum magnitude", "Physics", "G4BestUnit", "G4double"); 127 128 G4String NTP("NTP"); 129 (*store)[NTP] = G4AttDef(NTP, "No. of points", "Physics", "", "G4int"); 130 } 131 return store; 132 } 133 134 std::vector<G4AttValue>* G4ClonedTrajectory::CreateAttValues() const 135 { 136 auto values = new std::vector<G4AttValue>; 137 138 values->push_back(G4AttValue("ID", G4UIcommand::ConvertToString(fTrackID), "")); 139 140 values->push_back(G4AttValue("PID", G4UIcommand::ConvertToString(fParentID), "")); 141 142 values->push_back(G4AttValue("PN", ParticleName, "")); 143 144 values->push_back(G4AttValue("Ch", G4UIcommand::ConvertToString(PDGCharge), "")); 145 146 values->push_back(G4AttValue("PDG", G4UIcommand::ConvertToString(PDGEncoding), "")); 147 148 values->push_back(G4AttValue("IKE", G4BestUnit(initialKineticEnergy, "Energy"), "")); 149 150 values->push_back(G4AttValue("IMom", G4BestUnit(initialMomentum, "Energy"), "")); 151 152 values->push_back(G4AttValue("IMag", G4BestUnit(initialMomentum.mag(), "Energy"), "")); 153 154 values->push_back(G4AttValue("NTP", G4UIcommand::ConvertToString(GetPointEntries()), "")); 155 156 #ifdef G4ATTDEBUG 157 G4cout << G4AttCheck(values, GetAttDefs()); 158 #endif 159 160 return values; 161 } 162 163 void G4ClonedTrajectory::AppendStep(const G4Step* aStep) 164 { 165 positionRecord->push_back(new G4ClonedTrajectoryPoint(aStep->GetPostStepPoint()->GetPosition())); 166 } 167 168 G4ParticleDefinition* G4ClonedTrajectory::GetParticleDefinition() 169 { 170 return (G4ParticleTable::GetParticleTable()->FindParticle(ParticleName)); 171 } 172 173 void G4ClonedTrajectory::MergeTrajectory(G4VTrajectory* secondTrajectory) 174 { 175 if (secondTrajectory == nullptr) return; 176 177 auto seco = (G4ClonedTrajectory*)secondTrajectory; 178 G4int ent = seco->GetPointEntries(); 179 for (G4int i = 1; i < ent; ++i) // initial pt of 2nd trajectory shouldn't be merged 180 { 181 positionRecord->push_back((*(seco->positionRecord))[i]); 182 } 183 delete (*seco->positionRecord)[0]; 184 seco->positionRecord->clear(); 185 } 186