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 // G4ClonedRichTrajectoryPoint class implementation 27 // 28 // Makoto Asai (JLab) - Oct.2024 29 // -------------------------------------------------------------------- 30 31 #include "G4ClonedRichTrajectoryPoint.hh" 32 #include "G4RichTrajectoryPoint.hh" 33 34 #include "G4AttDef.hh" 35 #include "G4AttDefStore.hh" 36 #include "G4AttValue.hh" 37 #include "G4Step.hh" 38 #include "G4Track.hh" 39 #include "G4UIcommand.hh" 40 #include "G4UnitsTable.hh" 41 #include "G4VProcess.hh" 42 43 // #define G4ATTDEBUG 44 #ifdef G4ATTDEBUG 45 # include "G4AttCheck.hh" 46 #endif 47 48 #include <sstream> 49 50 G4Allocator<G4ClonedRichTrajectoryPoint>*& aClonedRichTrajectoryPointAllocator() 51 { 52 static G4Allocator<G4ClonedRichTrajectoryPoint>* _instance = nullptr; 53 return _instance; 54 } 55 56 G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint(const G4Track* aTrack) 57 : fPosition(aTrack->GetPosition()), 58 fPreStepPointGlobalTime(aTrack->GetGlobalTime()), 59 fPostStepPointGlobalTime(aTrack->GetGlobalTime()), 60 fpPreStepPointVolume(aTrack->GetTouchableHandle()), 61 fpPostStepPointVolume(aTrack->GetNextTouchableHandle()), 62 fPreStepPointWeight(aTrack->GetWeight()), 63 fPostStepPointWeight(aTrack->GetWeight()) 64 {} 65 66 G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint(const G4Step* aStep) 67 : fPosition(aStep->GetPostStepPoint()->GetPosition()), 68 fpAuxiliaryPointVector(aStep->GetPointerToVectorOfAuxiliaryPoints()), 69 fTotEDep(aStep->GetTotalEnergyDeposit()) 70 { 71 G4StepPoint* preStepPoint = aStep->GetPreStepPoint(); 72 G4StepPoint* postStepPoint = aStep->GetPostStepPoint(); 73 if (aStep->GetTrack()->GetCurrentStepNumber() <= 0) // First step 74 { 75 fRemainingEnergy = aStep->GetTrack()->GetKineticEnergy(); 76 } 77 else { 78 fRemainingEnergy = preStepPoint->GetKineticEnergy() - fTotEDep; 79 } 80 fpProcess = postStepPoint->GetProcessDefinedStep(); 81 fPreStepPointStatus = preStepPoint->GetStepStatus(); 82 fPostStepPointStatus = postStepPoint->GetStepStatus(); 83 fPreStepPointGlobalTime = preStepPoint->GetGlobalTime(); 84 fPostStepPointGlobalTime = postStepPoint->GetGlobalTime(); 85 fpPreStepPointVolume = preStepPoint->GetTouchableHandle(); 86 fpPostStepPointVolume = postStepPoint->GetTouchableHandle(); 87 fPreStepPointWeight = preStepPoint->GetWeight(); 88 fPostStepPointWeight = postStepPoint->GetWeight(); 89 } 90 91 G4ClonedRichTrajectoryPoint::G4ClonedRichTrajectoryPoint(const G4RichTrajectoryPoint& right) 92 { 93 fPosition = right.fPosition; 94 fTotEDep = right.fTotEDep; 95 fRemainingEnergy = right.fRemainingEnergy; 96 fpProcess = right.fpProcess; 97 fPreStepPointStatus = right.fPreStepPointStatus; 98 fPostStepPointStatus = right.fPostStepPointStatus; 99 fPreStepPointGlobalTime = right.fPreStepPointGlobalTime; 100 fPostStepPointGlobalTime = right.fPostStepPointGlobalTime; 101 fpPreStepPointVolume = right.fpPreStepPointVolume; 102 fpPostStepPointVolume = right.fpPostStepPointVolume; 103 fPreStepPointWeight = right.fPreStepPointWeight; 104 fPostStepPointWeight = right.fPostStepPointWeight; 105 if(right.fpAuxiliaryPointVector!=nullptr && right.fpAuxiliaryPointVector->size()>0) 106 { 107 fpAuxiliaryPointVector = new std::vector<G4ThreeVector>; 108 for(auto& i : *(right.fpAuxiliaryPointVector)) 109 { fpAuxiliaryPointVector->push_back(i); } 110 } 111 } 112 113 G4ClonedRichTrajectoryPoint::~G4ClonedRichTrajectoryPoint() { delete fpAuxiliaryPointVector; } 114 115 const std::map<G4String, G4AttDef>* G4ClonedRichTrajectoryPoint::GetAttDefs() const 116 { 117 G4bool isNew; 118 std::map<G4String, G4AttDef>* store = G4AttDefStore::GetInstance("G4RichTrajectoryPoint", isNew); 119 if (isNew) { 120 G4String ID; 121 122 ID = "Pos"; 123 (*store)[ID] = G4AttDef(ID, "Position", "Physics", "G4BestUnit", "G4ThreeVector"); 124 ID = "Aux"; 125 (*store)[ID] = 126 G4AttDef(ID, "Auxiliary Point Position", "Physics", "G4BestUnit", "G4ThreeVector"); 127 ID = "TED"; 128 (*store)[ID] = G4AttDef(ID, "Total Energy Deposit", "Physics", "G4BestUnit", "G4double"); 129 ID = "RE"; 130 (*store)[ID] = G4AttDef(ID, "Remaining Energy", "Physics", "G4BestUnit", "G4double"); 131 ID = "PDS"; 132 (*store)[ID] = G4AttDef(ID, "Process Defined Step", "Physics", "", "G4String"); 133 ID = "PTDS"; 134 (*store)[ID] = G4AttDef(ID, "Process Type Defined Step", "Physics", "", "G4String"); 135 ID = "PreStatus"; 136 (*store)[ID] = G4AttDef(ID, "Pre-step-point status", "Physics", "", "G4String"); 137 ID = "PostStatus"; 138 (*store)[ID] = G4AttDef(ID, "Post-step-point status", "Physics", "", "G4String"); 139 ID = "PreT"; 140 (*store)[ID] = G4AttDef(ID, "Pre-step-point global time", "Physics", "G4BestUnit", "G4double"); 141 ID = "PostT"; 142 (*store)[ID] = G4AttDef(ID, "Post-step-point global time", "Physics", "G4BestUnit", "G4double"); 143 ID = "PreVPath"; 144 (*store)[ID] = G4AttDef(ID, "Pre-step Volume Path", "Physics", "", "G4String"); 145 ID = "PostVPath"; 146 (*store)[ID] = G4AttDef(ID, "Post-step Volume Path", "Physics", "", "G4String"); 147 ID = "PreW"; 148 (*store)[ID] = G4AttDef(ID, "Pre-step-point weight", "Physics", "", "G4double"); 149 ID = "PostW"; 150 (*store)[ID] = G4AttDef(ID, "Post-step-point weight", "Physics", "", "G4double"); 151 } 152 return store; 153 } 154 155 static G4String Status(G4StepStatus stps) 156 { 157 G4String status; 158 switch (stps) { 159 case fWorldBoundary: 160 status = "fWorldBoundary"; 161 break; 162 case fGeomBoundary: 163 status = "fGeomBoundary"; 164 break; 165 case fAtRestDoItProc: 166 status = "fAtRestDoItProc"; 167 break; 168 case fAlongStepDoItProc: 169 status = "fAlongStepDoItProc"; 170 break; 171 case fPostStepDoItProc: 172 status = "fPostStepDoItProc"; 173 break; 174 case fUserDefinedLimit: 175 status = "fUserDefinedLimit"; 176 break; 177 case fExclusivelyForcedProc: 178 status = "fExclusivelyForcedProc"; 179 break; 180 case fUndefined: 181 status = "fUndefined"; 182 break; 183 default: 184 status = "Not recognised"; 185 break; 186 } 187 return status; 188 } 189 190 static G4String Path(const G4TouchableHandle& th) 191 { 192 std::ostringstream oss; 193 G4int depth = th->GetHistoryDepth(); 194 for (G4int i = depth; i >= 0; --i) { 195 oss << th->GetVolume(i)->GetName() << ':' << th->GetCopyNumber(i); 196 if (i != 0) oss << '/'; 197 } 198 return oss.str(); 199 } 200 201 std::vector<G4AttValue>* G4ClonedRichTrajectoryPoint::CreateAttValues() const 202 { 203 // Create base class att values... 204 205 auto values = new std::vector<G4AttValue>; 206 values->push_back(G4AttValue("Pos", G4BestUnit(fPosition, "Length"), "")); 207 208 #ifdef G4ATTDEBUG 209 G4cout << G4AttCheck(values, GetAttDefs()); 210 #endif 211 212 if (fpAuxiliaryPointVector != nullptr) { 213 for (const auto& iAux : *fpAuxiliaryPointVector) { 214 values->push_back(G4AttValue("Aux", G4BestUnit(iAux, "Length"), "")); 215 } 216 } 217 218 values->push_back(G4AttValue("TED", G4BestUnit(fTotEDep, "Energy"), "")); 219 values->push_back(G4AttValue("RE", G4BestUnit(fRemainingEnergy, "Energy"), "")); 220 221 if (fpProcess != nullptr) { 222 values->push_back(G4AttValue("PDS", fpProcess->GetProcessName(), "")); 223 values->push_back( 224 G4AttValue("PTDS", G4VProcess::GetProcessTypeName(fpProcess->GetProcessType()), "")); 225 } 226 else { 227 values->push_back(G4AttValue("PDS", "None", "")); 228 values->push_back(G4AttValue("PTDS", "None", "")); 229 } 230 231 values->push_back(G4AttValue("PreStatus", Status(fPreStepPointStatus), "")); 232 233 values->push_back(G4AttValue("PostStatus", Status(fPostStepPointStatus), "")); 234 235 values->push_back(G4AttValue("PreT", G4BestUnit(fPreStepPointGlobalTime, "Time"), "")); 236 237 values->push_back(G4AttValue("PostT", G4BestUnit(fPostStepPointGlobalTime, "Time"), "")); 238 239 if (fpPreStepPointVolume && (fpPreStepPointVolume->GetVolume() != nullptr)) { 240 values->push_back(G4AttValue("PreVPath", Path(fpPreStepPointVolume), "")); 241 } 242 else { 243 values->push_back(G4AttValue("PreVPath", "None", "")); 244 } 245 246 if (fpPostStepPointVolume && (fpPostStepPointVolume->GetVolume() != nullptr)) { 247 values->push_back(G4AttValue("PostVPath", Path(fpPostStepPointVolume), "")); 248 } 249 else { 250 values->push_back(G4AttValue("PostVPath", "None", "")); 251 } 252 253 std::ostringstream oss1; 254 oss1 << fPreStepPointWeight; 255 values->push_back(G4AttValue("PreW", oss1.str(), "")); 256 257 std::ostringstream oss2; 258 oss2 << fPostStepPointWeight; 259 values->push_back(G4AttValue("PostW", oss2.str(), "")); 260 261 #ifdef G4ATTDEBUG 262 G4cout << G4AttCheck(values, GetAttDefs()); 263 #endif 264 265 return values; 266 } 267