Geant4 11.2.2
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
45#include <cmath>
46
48{
49 delete fpHeadPolyhedron;
50 delete fpShaftPolyhedron;
51}
52
54(G4double x1, G4double y1, G4double z1,
55 G4double x2, G4double y2, G4double z2,
56 G4double width, const G4Colour& colour,
57 const G4String& description,
58 G4int lineSegmentsPerCircle,
59 const G4Transform3D& transform)
60: fpShaftPolyhedron(nullptr)
61, fpHeadPolyhedron(nullptr)
62, fTransform(transform)
63{
64 fType = "G4ArrowModel";
66 fGlobalDescription = fType + ": " + description;
68 (std::min(x1,x2),
69 std::max(x1,x2),
70 std::min(y1,y2),
71 std::max(y1,y2),
72 std::min(z1,z2),
73 std::max(z1,z2));
74
75 // Force number of line segments per circle (aka number of rotation steps)
77 G4Polyhedron::SetNumberOfRotationSteps(lineSegmentsPerCircle);
78
79 // Make a cylinder slightly shorter than the arrow length so that it
80 // doesn't stick out of the head.
82
83 G4double totalLength = std::hypot(x2-x1, y2-y1, z2-z1);
84 if (totalLength < tolerance)
85 {totalLength = tolerance;}
86
87 G4double shaftRadius = width/6.;
88 if (shaftRadius < tolerance)
89 {shaftRadius = tolerance;}
90
91 // case 1 - arrow length >> width -> arrow head is width and 1.5x width tall
92 // case 2 - arrow length < width -> arrow head is made to be 0.5x length
93 G4double arrowLength = std::min(1.5*width, 0.5*totalLength);
94
95 G4double shaftLength = totalLength - arrowLength;
96 if (shaftLength < 2*tolerance)
97 {shaftLength = 2*tolerance;}
98
99 const G4Tubs shaft("shaft",0.,shaftRadius,0.5*shaftLength,0.,twopi);
100 fpShaftPolyhedron = shaft.CreatePolyhedron();
101 // translate the polyhedron down w.r.t. the centre of the whole arrow
102 if (fpShaftPolyhedron)
103 {fpShaftPolyhedron->Transform(G4Translate3D(0,0,-0.5*arrowLength));}
104
105 // Locate the head at +halfShaftLength.
106 const G4double zHi = 0.5*totalLength;
107 const G4double zLow = zHi - arrowLength;
108 const G4double rExt = 0.5*width;
109 const G4double xExt = std::sqrt(3.)*rExt/2.;
110 const G4Tet head("head",
111 G4ThreeVector(0.,0.,zHi),
112 G4ThreeVector(0.,rExt,zLow),
113 G4ThreeVector(xExt,-rExt/2.,zLow),
114 G4ThreeVector(-xExt,-rExt/2.,zLow));
115 fpHeadPolyhedron = head.CreatePolyhedron();
116
117 // Transform to position
118 const G4Vector3D arrowDirection = G4Vector3D(x2-x1,y2-y1,z2-z1).unit();
119 const G4double theta = arrowDirection.theta();
120 const G4double phi = arrowDirection.phi();
121 const G4Point3D arrowCentre(0.5*(x1+x2),0.5*(y1+y2),0.5*(z1+z2));
122 const G4Transform3D tr =
123 G4Translate3D(arrowCentre) * G4RotateZ3D(phi) * G4RotateY3D(theta);
124 if (fpShaftPolyhedron) fpShaftPolyhedron->Transform(tr);
125 if (fpHeadPolyhedron) fpHeadPolyhedron->Transform(tr);
126
128 va.SetColour(colour);
129 va.SetForceSolid(true);
130 if (fpShaftPolyhedron) fpShaftPolyhedron->SetVisAttributes(va);
131 if (fpHeadPolyhedron) fpHeadPolyhedron->SetVisAttributes(va);
132
133 // Restore number of line segments per circle
135}
136
138{
139 if (fpShaftPolyhedron && fpHeadPolyhedron) {
140 sceneHandler.BeginPrimitives(fTransform);
141 sceneHandler.AddPrimitive(*fpShaftPolyhedron);
142 sceneHandler.AddPrimitive(*fpHeadPolyhedron);
143 sceneHandler.EndPrimitives();
144 }
145}
CLHEP::Hep3Vector G4ThreeVector
HepGeom::RotateZ3D G4RotateZ3D
HepGeom::RotateY3D G4RotateY3D
HepGeom::Translate3D G4Translate3D
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
HepGeom::Vector3D< G4double > G4Vector3D
Definition G4Vector3D.hh:34
G4ArrowModel(G4double x1, G4double y1, G4double z1, G4double x2, G4double y2, G4double z2, G4double width, const G4Colour &colour, const G4String &description="", G4int lineSegmentsPerCircle=6, const G4Transform3D &transform=G4Transform3D())
virtual ~G4ArrowModel()
void DescribeYourselfTo(G4VGraphicsScene &) override
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
Definition G4Tet.hh:56
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tet.cc:689
G4Polyhedron * CreatePolyhedron() const override
Definition G4Tubs.cc:1737
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())=0
virtual void AddPrimitive(const G4Polyline &)=0
virtual void EndPrimitives()=0
G4VisExtent fExtent
Definition G4VModel.hh:97
G4String fGlobalDescription
Definition G4VModel.hh:96
G4String fType
Definition G4VModel.hh:94
G4String fGlobalTag
Definition G4VModel.hh:95
void SetColour(const G4Colour &)
void SetForceSolid(G4bool=true)
void SetVisAttributes(const G4VisAttributes *)
Definition G4Visible.cc:98
BasicVector3D< T > unit() const
static void SetNumberOfRotationSteps(G4int n)
static G4int GetNumberOfRotationSteps()
HepPolyhedron & Transform(const G4Transform3D &t)