Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ArrowModel.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// John Allison 15th July 2012
30// Model that knows how to draw an arrow.
31
32#include "G4ArrowModel.hh"
33
35#include "G4VGraphicsScene.hh"
36#include "G4VisAttributes.hh"
37#include "G4Tubs.hh"
38#include "G4Tet.hh"
39#include "G4Polyhedron.hh"
40#include "G4Vector3D.hh"
41#include "G4Point3D.hh"
42#include "G4Transform3D.hh"
44
46{
47 delete fpHeadPolyhedron;
48 delete fpShaftPolyhedron;
49}
50
52(G4double x1, G4double y1, G4double z1,
53 G4double x2, G4double y2, G4double z2,
54 G4double width, const G4Colour& colour,
55 const G4String& description,
56 G4int lineSegmentsPerCircle)
57: fpShaftPolyhedron(nullptr)
58, fpHeadPolyhedron(nullptr)
59{
60 fType = "G4ArrowModel";
62 fGlobalDescription = fType + ": " + description;
64 (std::min(x1,x2),
65 std::max(x1,x2),
66 std::min(y1,y2),
67 std::max(y1,y2),
68 std::min(z1,z2),
69 std::max(z1,z2));
70
71 // Force number of line segments per circle (aka number of rotation steps)
73 G4Polyhedron::SetNumberOfRotationSteps(lineSegmentsPerCircle);
74
75 // Make a cylinder slightly shorter than the arrow length so that it
76 // doesn't stick out of the head.
78 G4double shaftLength = std::sqrt
79 (std::pow(x2-x1,2)+std::pow(y2-y1,2)+std::pow(z2-z1,2));
80 if (shaftLength < tolerance) shaftLength = tolerance;
81 G4double shaftRadius = width/2.;
82 if (shaftRadius > shaftLength/100.) shaftRadius = shaftLength/100.;
83 if (shaftRadius < tolerance) shaftRadius = tolerance;
84 const G4double halfShaftLength = shaftLength/2.;
85 const G4double halfReduction = 4.*shaftRadius;
86 G4double halfLength = halfShaftLength - halfReduction;
87 if (halfLength < tolerance) halfLength = tolerance;
88 const G4Tubs shaft("shaft",0.,shaftRadius,halfLength,0.,twopi);
89 fpShaftPolyhedron = shaft.CreatePolyhedron();
90 // Move it a little so that the tail is at z = -halfShaftLength.
91 if (fpShaftPolyhedron)
92 fpShaftPolyhedron->Transform(G4Translate3D(0,0,-halfReduction));
93
94 // Locate the head at +halfShaftLength.
95 const G4double zHi = halfShaftLength;
96 const G4double zLow = halfShaftLength - 12.*shaftRadius;
97 const G4double rExt = 8. * shaftRadius;
98 const G4double xExt = std::sqrt(3.)*rExt/2.;
99 const G4Tet head("head",
100 G4ThreeVector(0.,0.,zHi),
101 G4ThreeVector(0.,rExt,zLow),
102 G4ThreeVector(xExt,-rExt/2.,zLow),
103 G4ThreeVector(-xExt,-rExt/2.,zLow));
104 fpHeadPolyhedron = head.CreatePolyhedron();
105
106 // Transform to position
107 const G4Vector3D arrowDirection = G4Vector3D(x2-x1,y2-y1,z2-z1).unit();
108 const G4double theta = arrowDirection.theta();
109 const G4double phi = arrowDirection.phi();
110 const G4Point3D arrowCentre(0.5*(x1+x2),0.5*(y1+y2),0.5*(z1+z2));
111 const G4Transform3D tr =
112 G4Translate3D(arrowCentre) * G4RotateZ3D(phi) * G4RotateY3D(theta);
113 if (fpShaftPolyhedron) fpShaftPolyhedron->Transform(tr);
114 if (fpHeadPolyhedron) fpHeadPolyhedron->Transform(tr);
115
117 va.SetColour(colour);
118 va.SetForceSolid(true);
119 if (fpShaftPolyhedron) fpShaftPolyhedron->SetVisAttributes(va);
120 if (fpHeadPolyhedron) fpHeadPolyhedron->SetVisAttributes(va);
121
122 // Restore number of line segments per circle
124}
125
127{
128 sceneHandler.BeginPrimitives();
129 if (fpShaftPolyhedron) sceneHandler.AddPrimitive(*fpShaftPolyhedron);
130 if (fpHeadPolyhedron) sceneHandler.AddPrimitive(*fpHeadPolyhedron);
131 sceneHandler.EndPrimitives();
132}
CLHEP::Hep3Vector G4ThreeVector
HepGeom::RotateZ3D G4RotateZ3D
HepGeom::Translate3D G4Translate3D
HepGeom::RotateY3D G4RotateY3D
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
virtual ~G4ArrowModel()
Definition: G4ArrowModel.cc:45
G4ArrowModel(G4double x1, G4double y1, G4double z1, G4double x2, G4double y2, G4double z2, G4double width, const G4Colour &colour, const G4String &description="", G4int lineSegmentsPerCircle=6)
Definition: G4ArrowModel.cc:52
virtual void DescribeYourselfTo(G4VGraphicsScene &)
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
Definition: G4Tet.hh:56
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:646
Definition: G4Tubs.hh:75
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tubs.cc:1759
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual void EndPrimitives()=0
G4VisExtent fExtent
Definition: G4VModel.hh:113
G4String fGlobalDescription
Definition: G4VModel.hh:112
G4String fType
Definition: G4VModel.hh:110
G4String fGlobalTag
Definition: G4VModel.hh:111
void SetColour(const G4Colour &)
void SetForceSolid(G4bool=true)
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:79
BasicVector3D< T > unit() const
static void SetNumberOfRotationSteps(G4int n)
static G4int GetNumberOfRotationSteps()
HepPolyhedron & Transform(const G4Transform3D &t)