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