Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 "G4VisManager.hh"
40#include "G4VisAttributes.hh"
41#include "G4Colour.hh"
43#include "G4ios.hh"
45
47{
49 return _instance;
50}
51
52G4RayTrajectory :: G4RayTrajectory()
53{
54 positionRecord = new std::vector<G4RayTrajectoryPoint*>;
55}
56
57G4RayTrajectory :: G4RayTrajectory(G4RayTrajectory & right)
59{
60 positionRecord = new std::vector<G4RayTrajectoryPoint*>;
61 for(size_t i=0;i<right.positionRecord->size();i++)
62 {
64 ((*(right.positionRecord))[i]);
65 positionRecord->push_back(new G4RayTrajectoryPoint(*rightPoint));
66 }
67}
68
69G4RayTrajectory :: ~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
79{
80 G4RayTrajectoryPoint* trajectoryPoint = new G4RayTrajectoryPoint();
81
82 const G4Step* aStep = theStep;
83 G4Navigator* theNavigator
85
86 // Take care of parallel world(s)
88 {
91 std::vector<G4Navigator*>::iterator iNav =
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
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
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
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
157void G4RayTrajectory::ShowTrajectory(std::ostream&) const
158{ }
159
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
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 void ShowTrajectory(std::ostream &) const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual void AppendStep(const G4Step *)
virtual G4int GetPointEntries() const
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:34
virtual G4int GetHistoryDepth() const
Definition: G4VTouchable.cc:74
G4VSceneHandler * GetCurrentSceneHandler() const
static G4VisManager * GetInstance()
#define G4ThreadLocalStatic
Definition: tls.hh:76