Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VtkUtility.cc File Reference
#include "G4VVtkPipeline.hh"
#include "G4VtkUtility.hh"
#include "G4Transform3D.hh"
#include "G4Colour.hh"
#include <vtkMatrix4x4.h>
#include <vtkPlane.h>

Go to the source code of this file.

Functions

vtkSmartPointer< vtkMatrix4x4 > G4Transform3DToVtkMatrix4x4 (const G4Transform3D &g4Trans)
 
vtkSmartPointer< vtkPlane > G4Plane3DToVtkPlane (const G4Plane3D &g4plane)
 
G4Plane3D VtkPlaneToG4Plane3D (vtkPlane *vtkPlane)
 
void MaxBounds (G4double *maxBound, G4double *boundToCheck)
 
std::size_t MakeHash (const G4ThreeVector &v)
 
std::size_t MakeHash (const G4Colour &c)
 

Function Documentation

◆ G4Plane3DToVtkPlane()

vtkSmartPointer< vtkPlane > G4Plane3DToVtkPlane ( const G4Plane3D & g4plane)

Definition at line 50 of file G4VtkUtility.cc.

51{
52 auto plane = vtkSmartPointer<vtkPlane>::New();
53
54 auto g4normal = g4plane.normal();
55 auto g4point = g4plane.point();
56
57 plane->SetNormal(g4normal.x(), g4normal.y(), g4normal.z());
58 plane->SetOrigin(g4point.x(), g4point.y(), g4point.z());
59
60 return plane;
61}
Point3D< T > point(const Point3D< T > &p) const
Definition Plane3D.h:115
Normal3D< T > normal() const
Definition Plane3D.h:97

Referenced by G4VtkViewer::AddClipperPlaneWidget(), and G4VtkViewer::AddCutterPlaneWidget().

◆ G4Transform3DToVtkMatrix4x4()

vtkSmartPointer< vtkMatrix4x4 > G4Transform3DToVtkMatrix4x4 ( const G4Transform3D & g4Trans)

Definition at line 38 of file G4VtkUtility.cc.

39{
40 auto transform = vtkSmartPointer<vtkMatrix4x4>();
41
42 double transformArray[16] = {g4Trans.xx(), g4Trans.xy(), g4Trans.xy(), g4Trans.dx(),
43 g4Trans.yx(), g4Trans.yy(), g4Trans.yy(), g4Trans.dy(),
44 g4Trans.zx(), g4Trans.zy(), g4Trans.zy(), g4Trans.dz(),
45 0., 0., 0., 1.};
46 transform->DeepCopy(transformArray);
47 return transform;
48}
double dy() const
double dz() const
double dx() const
double xy() const
double zx() const
double yx() const
double zy() const
double xx() const
double yy() const

◆ MakeHash() [1/2]

std::size_t MakeHash ( const G4Colour & c)

Definition at line 99 of file G4VtkUtility.cc.

100{
101
102 std::size_t rhash = std::hash<G4double>{}(c.GetRed());
103 std::size_t ghash = std::hash<G4double>{}(c.GetGreen());
104 std::size_t bhash = std::hash<G4double>{}(c.GetBlue());
105 std::size_t ahash = std::hash<G4double>{}(c.GetAlpha());
106
107 std::hash_combine(rhash, ghash);
108 std::hash_combine(rhash, bhash);
109 std::hash_combine(rhash, ahash);
110
111 return rhash;
112}
G4double GetBlue() const
Definition G4Colour.hh:154
G4double GetAlpha() const
Definition G4Colour.hh:155
G4double GetRed() const
Definition G4Colour.hh:152
G4double GetGreen() const
Definition G4Colour.hh:153
void hash_combine(std::size_t)

◆ MakeHash() [2/2]

std::size_t MakeHash ( const G4ThreeVector & v)

Definition at line 86 of file G4VtkUtility.cc.

87{
88
89 std::size_t xhash = std::hash<G4double>{}(v.x());
90 std::size_t yhash = std::hash<G4double>{}(v.y());
91 std::size_t zhash = std::hash<G4double>{}(v.z());
92
93 std::hash_combine(xhash, yhash);
94 std::hash_combine(yhash, zhash);
95
96 return xhash;
97}
double z() const
double x() const
double y() const

◆ MaxBounds()

void MaxBounds ( G4double * maxBound,
G4double * boundToCheck )

Definition at line 73 of file G4VtkUtility.cc.

74{
75 if (boundToCheck[0] < maxBound[0]) maxBound[0] = boundToCheck[0];
76 if (boundToCheck[1] > maxBound[1]) maxBound[1] = boundToCheck[1];
77
78 if (boundToCheck[2] < maxBound[2]) maxBound[2] = boundToCheck[2];
79 if (boundToCheck[3] > maxBound[3]) maxBound[3] = boundToCheck[3];
80
81 if (boundToCheck[4] < maxBound[4]) maxBound[4] = boundToCheck[4];
82 if (boundToCheck[5] > maxBound[5]) maxBound[5] = boundToCheck[5];
83}

Referenced by G4VtkStore::GetBounds().

◆ VtkPlaneToG4Plane3D()

G4Plane3D VtkPlaneToG4Plane3D ( vtkPlane * vtkPlane)

Definition at line 63 of file G4VtkUtility.cc.

64{
65 auto point = vtkPlane->GetOrigin();
66 auto normal = vtkPlane->GetNormal();
67
68 G4Plane3D plane =
69 G4Plane3D(G4Normal3D(normal[0], normal[1], normal[2]), G4Point3D(point[0], point[1], point[2]));
70 return plane;
71}
HepGeom::Normal3D< G4double > G4Normal3D
Definition G4Normal3D.hh:34
HepGeom::Plane3D< G4double > G4Plane3D
Definition G4Plane3D.hh:34
HepGeom::Point3D< G4double > G4Point3D
Definition G4Point3D.hh:34

Referenced by vtkIPWCallback::Execute().