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 // 28 // 29 // 30 31 /////////////////// 32 //G4RayTrajectory.cc 33 /////////////////// 34 35 #include "G4RayTrajectory.hh" 36 #include "G4RayTrajectoryPoint.hh" 37 #include "G4RayTracerSceneHandler.hh" 38 #include "G4Step.hh" 39 #include "G4VisManager.hh" 40 #include "G4VisAttributes.hh" 41 #include "G4Colour.hh" 42 #include "G4TransportationManager.hh" 43 #include "G4ios.hh" 44 #include "G4ParallelWorldProcess.hh" 45 46 G4Allocator<G4RayTrajectory>*& rayTrajectoryAllocator() 47 { 48 G4ThreadLocalStatic G4Allocator<G4RayTrajectory>* _instance = nullptr; 49 return _instance; 50 } 51 52 G4RayTrajectory :: G4RayTrajectory() 53 { 54 positionRecord = new std::vector<G4RayTrajectoryPoint*>; 55 } 56 57 G4RayTrajectory :: G4RayTrajectory(G4RayTrajectory & right) 58 : G4VTrajectory() 59 { 60 positionRecord = new std::vector<G4RayTrajectoryPoint*>; 61 for(size_t i=0;i<right.positionRecord->size();i++) 62 { 63 G4RayTrajectoryPoint* rightPoint = (G4RayTrajectoryPoint*) 64 ((*(right.positionRecord))[i]); 65 positionRecord->push_back(new G4RayTrajectoryPoint(*rightPoint)); 66 } 67 } 68 69 G4RayTrajectory :: ~G4RayTrajectory() 70 { 71 //positionRecord->clearAndDestroy(); 72 for(size_t i=0;i<positionRecord->size();i++) 73 { delete (*positionRecord)[i]; } 74 positionRecord->clear(); 75 delete positionRecord; 76 } 77 78 void G4RayTrajectory::AppendStep(const G4Step* theStep) 79 { 80 G4RayTrajectoryPoint* trajectoryPoint = new G4RayTrajectoryPoint(); 81 82 const G4Step* aStep = theStep; 83 G4Navigator* theNavigator 84 = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking(); 85 86 // Take care of parallel world(s) 87 if(G4ParallelWorldProcess::GetHyperStep()) 88 { 89 aStep = G4ParallelWorldProcess::GetHyperStep(); 90 G4int navID = G4ParallelWorldProcess::GetHypNavigatorID(); 91 std::vector<G4Navigator*>::iterator iNav = 92 G4TransportationManager::GetTransportationManager()-> 93 GetActiveNavigatorsIterator(); 94 theNavigator = iNav[navID]; 95 } 96 97 trajectoryPoint->SetStepLength(aStep->GetStepLength()); 98 99 // Surface normal 100 G4bool valid; 101 G4ThreeVector theLocalNormal = theNavigator->GetLocalExitNormal(&valid); 102 if(valid) { theLocalNormal = -theLocalNormal; } 103 G4ThreeVector theGrobalNormal 104 = theNavigator->GetLocalToGlobalTransform().TransformAxis(theLocalNormal); 105 trajectoryPoint->SetSurfaceNormal(theGrobalNormal); 106 107 G4VisManager* visManager = G4VisManager::GetInstance(); 108 G4RayTracerSceneHandler* sceneHandler = 109 static_cast<G4RayTracerSceneHandler*>(visManager->GetCurrentSceneHandler()); 110 const auto& sceneVisAttsMap = sceneHandler->GetSceneVisAttsMap(); 111 112 // Make a path from the preStepPoint touchable 113 G4StepPoint* preStepPoint = aStep -> GetPreStepPoint(); 114 const G4VTouchable* preTouchable = preStepPoint->GetTouchable(); 115 G4int preDepth = preTouchable->GetHistoryDepth(); 116 G4ModelingParameters::PVPointerCopyNoPath localPrePVPointerCopyNoPath; 117 for (G4int i = preDepth; i >= 0; --i) { 118 localPrePVPointerCopyNoPath.push_back 119 (G4ModelingParameters::PVPointerCopyNo 120 (preTouchable->GetVolume(i),preTouchable->GetCopyNumber(i))); 121 } 122 123 // Pick up the vis atts, if any, from the scene handler 124 auto preIterator = sceneVisAttsMap.find(localPrePVPointerCopyNoPath); 125 const G4VisAttributes* preVisAtts; 126 if (preIterator != sceneVisAttsMap.end()) { 127 preVisAtts = &preIterator->second; 128 } else { 129 preVisAtts = 0; 130 } 131 trajectoryPoint->SetPreStepAtt(preVisAtts); 132 133 // Make a path from the postStepPoint touchable 134 G4StepPoint* postStepPoint = aStep -> GetPostStepPoint(); 135 const G4VTouchable* postTouchable = postStepPoint->GetTouchable(); 136 G4int postDepth = postTouchable->GetHistoryDepth(); 137 G4ModelingParameters::PVPointerCopyNoPath localPostPVPointerCopyNoPath; 138 for (G4int i = postDepth; i >= 0; --i) { 139 localPostPVPointerCopyNoPath.push_back 140 (G4ModelingParameters::PVPointerCopyNo 141 (postTouchable->GetVolume(i),postTouchable->GetCopyNumber(i))); 142 } 143 144 // Pick up the vis atts, if any, from the scene handler 145 auto postIterator = sceneVisAttsMap.find(localPostPVPointerCopyNoPath); 146 const G4VisAttributes* postVisAtts; 147 if (postIterator != sceneVisAttsMap.end()) { 148 postVisAtts = &postIterator->second; 149 } else { 150 postVisAtts = 0; 151 } 152 trajectoryPoint->SetPostStepAtt(postVisAtts); 153 154 positionRecord->push_back(trajectoryPoint); 155 } 156 157 void G4RayTrajectory::ShowTrajectory(std::ostream&) const 158 { } 159 160 void G4RayTrajectory::MergeTrajectory(G4VTrajectory* secondTrajectory) 161 { 162 if(!secondTrajectory) return; 163 164 G4RayTrajectory* seco = (G4RayTrajectory*)secondTrajectory; 165 G4int ent = seco->GetPointEntries(); 166 for(G4int i=0;i<ent;i++) 167 { positionRecord->push_back((G4RayTrajectoryPoint*)seco->GetPoint(i)); } 168 seco->positionRecord->clear(); 169 } 170 171