Geant4 User's Guide
For Application Developers Geometry |
The STEP standard supports multiple solid representations. Constructive
Solid Geometry (CSG) representations and Boundary Represented Solids (BREPs)
are available. Different representations are suitable for different
purposes, applications, required complexity, and levels of detail.
CSG representations are easy to use and normally give superior performance,
but they cannot reproduce complex solids such as those used in CAD systems.
BREP representations can handle more extended topologies and reproduce the
most complex solids, thus allowing the exchange of models with CAD systems.
All constructed solids can stream out their contents via appropriate methods
and streaming operators.
For all solids it is possible to estimate the geometrical volume and the surface area by invoking the methods:
G4double GetCubicVolume() G4double GetSurfaceArea()which return an estimate of the solid volume and total area in internal units respectively. For elementary solids the functions compute the exact geometrical quantities, while for composite or complex solids an estimate is made using Monte Carlo techniques.
For all solids it is also possible to generate pseudo-random points lying on their surfaces, by invoking the method
G4ThreeVector GetPointOnSurface() constwhich returns the generated point in local coordinates relative to the solid.
G4Box* aBox = new G4Box("BoxA", 1.0*cm, 3.0*cm, 5.0*cm);
Similarly to create a cylindrical section or tube,
one would use the constructor:
G4Tubs(const G4String& pName,
In the picture: pRMin = 10, pRMax = 15, pDz = 20
pSPhi = 0*Degree, pDPhi = 90*Degree |
pRMin | Inner radius | pRMax | Outer radius |
pDz | half length in z | pSPhi | the starting phi angle in radians |
pDPhi | the angle of the segment in radians | &n |
Similarly to create a cone, or conical section,
one would use the constructor:
G4Cons(const G4String& pName,
In the picture: pRmin1 = 5, pRmax1 = 10,
pRmin2 = 20, pRmax2 = 25 pDz = 40, pSPhi = 0, pDPhi = 4/3*Pi |
pRmin1 inside radius at -pDz
| pRmax1 | outside radius at -pDz
| pRmin2 | inside radius at +pDz
| pRmax2 | outside radius at +pDz
| pDz | half length in z
| pSPhi | starting angle of the segment in radians
| pDPhi | the angle of the segment in radians
| |
| |
A parallelepiped is constructed using:
G4Para(const G4String& pName,
In the picture: dx = 30, dy = 40, dz = 60
alpha = 10*Degree, theta = 20*Degree, phi = 5*Degree |
dx,dy,dz | Half-length in x,y,z |
alpha | Angle formed by the y axis and by the plane joining the centre of the faces parallel to the z-x plane at -dy and +dy |
theta | Polar angle of the line joining the centres of the faces at -dz and +dz in z |
phi | Azimuthal angle of the line joining the centres of the faces at -dz and +dz in z |
To construct a trapezoid use:
G4Trd(const G4String& pName,
In the picture: dx1 = 30, dx2 = 10
dy1 = 40, dy2 = 15 dz = 60 |
dx1 | Half-length along x at the surface positioned at -dz |
dx2 | Half-length along x at the surface positioned at +dz |
dy1 | Half-length along y at the surface positioned at -dz |
dy2 | Half-length along y at the surface positioned at +dz |
dz | Half-length along z axis |
To build a generic trapezoid, the G4Trap class is provided.
Here are the two costructors for a Right Angular Wedge and for the general
trapezoid for it:
G4Trap(const G4String& pName,
In the picture: pDx1 = 30, pDx2 = 40, pDy1 = 40
pDx3 = 10, pDx4 = 14, pDy2 = 16 pDz = 60, pTheta = 20*Degree pDphi = 5*Degree, pAlph1 = pAlph2 = 10*Degree |
pZ | Length along z |
pY | Length along y |
pX | Length along x at the wider side |
pLTX | Length along x at the narrower side (plTX<=pX) |
or to obtain the general trapezoid (see the Software Reference Manual):
pDx1 | Half x length at y=-pDy |
pDx2 | Half x length at y=+pDy |
pDy | Half y length |
pDz | Half z length |
pTheta | Polar angle of the line joining the centres of the faces at -/+pDz |
pDy1 | Half y length at -pDz |
pDx1 | Half x length at -pDz, y=-pDy1 |
pDx2 | Half x length at -pDz, y=+pDy1 |
pDy2 | Half y length at +pDz |
pDx3 | Half x length at +pDz, y=-pDy2 |
pDx4 | Half x length at +pDz, y=+pDy2 |
pAlph1 | Angle with respect to the y axis from the centre of the side (lower endcap) |
pAlph2 | Angle with respect to the y axis from the centre of the side (upper endcap) |
Note on pAlph1/2: the two angles have to be the same due to the planarity condition.
To build a sphere, or a spherical shell section, use:
G4Sphere(const G4String& pName,
In the picture: pRmin = 100, pRmax = 120 pSPhi = 0*Degree, pDPhi = 180*Degree pSTheta = 0 Degree, pDTheta = 180*Degree |
pRmin | Inner radius |
pRmax | Outer radius |
pSPhi | Starting Phi angle of the segment in radians |
pDPhi | Delta Phi angle of the segment in radians |
pSTheta | Starting Theta angle of the segment in radians |
pDTheta | Delta Theta angle of the segment in radians |
To build a full solid sphere use:
G4Orb(const G4String& pName,
In the picture: pRmax = 100
The Orb can be obtained from a Sphere with:
|
To build a torus use:
G4Torus(const G4String& pName,
In the picture: pRmin = 40, pRmax = 60, pRtor = 200
pSPhi = 0, pDPhi = 90*Degree |
pRmin | Inside radius |
pRmax | Outside radius |
pRtor | Swept radius of torus |
pSPhi | Starting Phi angle in radians (fSPhi+fDPhi<=2PI, fSPhi>-2PI) |
pDPhi | Delta angle of the segment in radians |
In addition, the Geant4 Design Documentation shows in the Solids Class Diagram the complete list of CSG classes, and the STEP documentation contains a detailed EXPRESS description of each CSG solid.
Specific CSG SolidsPolycons (PCON) are implemented in Geant4 through the G4Polycon class:
G4Polycone(const G4String& pName, G4double phiStart, G4double phiTotal, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[]) G4Polycone(const G4String& pName, G4double phiStart, G4double phiTotal, G4int numRZ, const G4double r[], const G4double z[])
In the picture: phiStart = 0*Degree, phiTotal = 2*Pi
numZPlanes = 9 rInner = { 0, 0, 0, 0, 0, 0, 0, 0, 0} rOuter = { 0, 10, 10, 5 , 5, 10 , 10 , 2, 2} z = { 5, 7, 9, 11, 25, 27, 29, 31, 35 } |
phiStart | Initial Phi starting angle |
phiTotal | Total Phi angle |
numZPlanes | Number of z planes |
numRZ | Number of corners in r,z space |
zPlane | Position of z planes |
rInner | Tangent distance to inner surface |
rOuter | Tangent distance to outer surface |
r | r coordinate of corners |
z | z coordinate of corners |
Polyhedra (PGON) are implemented through G4Polyhedra:
G4Polyhedra(const G4String& pName, G4double phiStart, G4double phiTotal, G4int numSide, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[] ) G4Polyhedra(const G4String& pName, G4double phiStart, G4double phiTotal, G4int numSide, G4int numRZ, const G4double r[], const G4double z[])
In the picture: phiStart = 0, phiTotal= 2 Pi
numSide = 5, nunZPlanes = 7 rInner = { 0, 0, 0, 0, 0, 0, 0 } rOuter = { 0, 15, 15, 4, 4, 10, 10 } z = { 0, 5, 8, 13 , 30, 32, 35 } |
phiStart | Initial Phi starting angle |
phiTotal | Total Phi angle |
numSide | Number of sides |
numZPlanes | Number of z planes |
numRZ | Number of corners in r,z space |
zPlane | Position of z planes |
rInner | Tangent distance to inner surface |
rOuter | Tangent distance to outer surface |
r | r coordinate of corners |
z | z coordinate of corners |
A tube with an elliptical cross section (ELTU) can be defined
as follows:
G4EllipticalTube(const G4String& pName, The equation of the surface in x/y is 1.0 = (x/dx)**2 + (y/dy)**2
In the picture: Dx = 5, Dy = 10, Dz = 20
|
The general ellipsoid with possible cut in Z can be
defined as follows:
G4Ellipsoid(const G4String& pName,
In the picture: pxSemiAxis = 10, pySemiAxis = 20, pzSemiAxis = 50
pzBottomCut = -10, pzTopCut = 40 |
1.0 = (x/pxSemiAxis)**2 + (y/pySemiAxis)**2 + (z/pzSemiAxis)**2
where:
pxSemiAxis | Semiaxis in X |
pySemiAxis | Semiaxis in Y |
pzSemiAxis | Semiaxis in Z |
pzBottomCut | lower cut plane level, z |
pzTopCut | upper cut plane level, z |
A cone with an elliptical cross section can be defined as follows:
G4EllipticalCone(const G4String& pName, G4double pxSemiAxis, G4double pySemiAxis, G4double zMax, G4double pzTopCut)
In the picture: pxSemiAxis = 30, pySemiAxis = 60
zMax = 50, pzTopCut = 25 |
pxSemiAxis | Semiaxis in X |
pySemiAxis | Semiaxis in Y |
zMax | Height of elliptical cone |
pzTopCut | upper cut plane level |
An elliptical cone of height zMax, semiaxis pxSemiAxis, and semiaxis pySemiAxis is given by the parametric equations:
x = pxSemiAxis * ( zMax - u ) / u * Cos(v)
y = pySemiAxis * ( zMax - u ) / u * Sin(v)
z = u
Where v is between 0 and 2*Pi, and u between 0 and h respectively.
A tube with a hyperbolic profile (HYPE) can be defined as follows:
G4Hype(const G4String& pName,
In the picture: innerStereo = 0.7, outerStereo = 0.7
halfLenZ = 50 innerRadius = 20, outerRadius = 30 |
G4Hype is shaped with curved sides parallel to the z-axis,
has a specified half-length along the z axis about which it is
centred, and a given minimum and maximum radius.
A minimum radius of 0 defines a filled Hype (with hyperbolic
inner surface), i.e. inner radius = 0 AND inner stereo angle = 0.
The inner and outer hyperbolic surfaces can have different
stereo angles. A stereo angle of 0 gives a cylindrical surface:
innerRadius | Inner radius |
outerRadius | Outer radius |
innerStereo | Inner stereo angle in radians |
outerStereo | Outer stereo angle in radians |
halfLenZ | Half length in Z |
A tetrahedra solid can be defined as follows:
G4Tet(const G4String& pName,
In the picture: anchor = {0, 0, sqrt(3)}
p2 = { 0, 2*sqrt(2/3), -1/sqrt(3) } p3 = { -sqrt(2), -sqrt(2/3),-1/sqrt(3) } p4 = { sqrt(2), -sqrt(2/3) , -1/sqrt(3) } |
anchor | Anchor point |
p2 | Point 2 |
p3 | Point 3 |
p4 | Point 4 |
degeneracyFlag | Flag indicating degeneracy of points |
A box twisted along one axis can be defined as follows:
G4TwistedBox(const G4String& pName,
In the picture: twistedangle = 30*Degree
pDx = 30, pDy =40, pDz = 60 |
G4TwistedBox is a box twisted along the z-axis. The twist angle cannot be greater than 90 degrees:
twistedangle | Twist angle |
pDx | Half x length |
pDy | Half y length |
pDz | Half z length |
A trapezoid twisted along one axis can be defined as follows:
G4TwistedTrap(const G4String& pName, G4double twistedangle, G4double pDxx1, G4double pDxx2, G4double pDy, G4double pDz) G4TwistedTrap(const G4String& pName, G4double twistedangle, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlph)
In the picture: pDx1 = 30, pDx2 = 40, pDy1 = 40
pDx3 = 10, pDx4 = 14, pDy2 = 16 pDz = 60 pTheta = 20*Degree, pDphi = 5*Degree pAlph = 10*Degree, twistedangle = 30*Degree |
The first constructor of G4TwistedTrap produces a regular trapezoid
twisted along the z-axis, where the caps of the trapezoid are of the
same shape and size.
The second constructor produces a generic trapezoid with
polar, azimuthal and tilt angles.
The twist angle cannot be greater than 90 degrees:
twistedangle | Twisted angle |
pDx1 | Half x length at y=-pDy |
pDx2 | Half x length at y=+pDy |
pDy | Half y length |
pDz | Half z length |
pTheta | Polar angle of the line joining the centres of the faces at -/+pDz |
pDy1 | Half y length at -pDz |
pDx1 | Half x length at -pDz, y=-pDy1 |
pDx2 | Half x length at -pDz, y=+pDy1 |
pDy2 | Half y length at +pDz |
pDx3 | Half x length at +pDz, y=-pDy2 |
pDx4 | Half x length at +pDz, y=+pDy2 |
pAlph | Angle with respect to the y axis from the centre of the side |
A twisted trapezoid with the x and y dimensions
varying along z can be defined as follows:
G4TwistedTrd(const G4String& pName,
In the picture: dx1 = 30, dx2 = 10
dy1 = 40, dy2 = 15 dz = 60 twistedangle = 30*Degree |
pDx1 | Half x length at the surface positioned at -dz |
pDx2 | Half x length at the surface positioned at +dz |
pDy1 | Half y length at the surface positioned at -dz |
pDy2 | Half y length at the surface positioned at +dz |
pDz | Half z length |
twistedangle | Twisted angle |
A tube section twisted along its axis can be defined as follows:
G4TwistedTubs(const G4String& pName,
In the picture: endinnerrad = 10, endouterrad = 15
halfzlen = 20, dphi = 90*Degree twistedangle = 60*Degree |
G4TwistedTubs is a sort of twisted cylinder which, placed along
the z-axis and divided into phi-segments is shaped like an
hyperboloid, where each of its segmented pieces can be tilted with a stereo
angle.
It can have inner and outer surfaces with the same stereo angle:
twistedangle | Twisted angle |
endinnerrad | Inner radius at endcap |
endouterrad | Outer radius at endcap |
halfzlen | Half z length |
dphi | Phi angle of a segment |
Additional constructors are provided, allowing the shape to be specified either as:
Creating such a new Boolean solid, requires:
The solids used should be either CSG solids (for examples a box, a spherical shell, or a tube) or another Boolean solid: the product of a previous Boolean operation. An important purpose of Boolean solids is to allow the description of solids with peculiar shapes in a simple and intuitive way, still allowing an efficient geometrical navigation inside them.
Note: The solids used can actually be of any type. However, in order to fully support the export of a Geant4 solid model via STEP to CAD systems, we restrict the use of Boolean operations to this subset of solids. But this subset contains all the most interesting use cases.
Note: The tracking cost for navigating in a Boolean solid in the current implementation, is proportional to the number of constituent solids. So care must be taken to avoid extensive, unecessary use of Boolean solids in performance-critical areas of a geometry description, where each solid is created from Boolean combinations of many other solids.
Examples of the creation of the simplest Boolean solids are given below:
G4Box* box = new G4Box("Box",20*mm,30*mm,40*mm); G4Tubs* cyl = new G4Tubs("Cylinder",0,50*mm,50*mm,0,twopi); // r: 0 mm -> 50 mm // z: -50 mm -> 50 mm // phi: 0 -> 2 pi G4UnionSolid* union = new G4UnionSolid("Box+Cylinder", box, cyl); G4IntersectionSolid* intersection = new G4IntersectionSolid("Box*Cylinder", box, cyl); G4SubtractionSolid* subtraction = new G4SubtractionSolid("Box-Cylinder", box, cyl);where the union, intersection and subtraction of a box and cylinder are constructed.
The more useful case where one of the solids is displaced from the origin of coordinates also exists. In this case the second solid is positioned relative to the coordinate system (and thus relative to the first). This can be done in two ways:
In the first case, the translation is applied first to move the origin of coordinates. Then the rotation is used to rotate the coordinate system of the second solid to the coordinate system of the first.
G4RotationMatrix* yRot = new G4RotationMatrix; // Rotates X and Z axes only yRot->rotateY(M_PI/4.*rad); // Rotates 45 degrees G4ThreeVector zTrans(0, 0, 50); G4UnionSolid* unionMoved = new G4UnionSolid("Box+CylinderMoved", box, cyl, yRot, zTrans); // // The new coordinate system of the cylinder is translated so that // its centre is at +50 on the original Z axis, and it is rotated // with its X axis halfway between the original X and Z axes. // Now we build the same solid using the alternative method // G4RotationMatrix invRot = *(yRot->invert()); G4Transform3D transform(invRot, zTrans); G4UnionSolid* unionMoved = new G4UnionSolid("Box+CylinderMoved", box, cyl, transform);Note that the first constructor that takes a pointer to the rotation-matrix (G4RotationMatrix*), does NOT copy it. Therefore once used a rotation-matrix to construct a Boolean solid, it must NOT be modified.
When positioning a volume associated to a Boolean solid, the relative center of coordinates considered for the positioning is the one related to the first of the two constituent solids.
In addition, the boundary surfaces can be made of Bezier surfaces and B-Splines,
or of NURBS (Non-Uniform-Rational-B-Splines) surfaces.
The resulting solids are Advanced BREPS.
Note - Currently, the implementation for surfaces generated by Beziers, B-Splines
or NURBS is only at the level of prototype and not fully functional.
Extensions are foreseen in the near future, also to allow exchange of geometrical
models with CAD systems.
We have defined a few simple Elementary BREPS, that can be instantiated simply by a user in a manner similar to the construction of Constructed Solids (CSGs). We summarize their capabilities in the following section.
However most BREPS Solids are defined by creating each surface separately and tying them together. Though it is possible to do this using code, it is potentially error prone. So it is generally much more productive to utilize a tool to create these volumes, and the tools of choice are CAD systems. In future, it will be possible to import/export models created in CAD systems, utilizing the STEP standard.
Specific BREP Solids
We have defined one polygonal and one polyconical shape using BREPS. The polycone provides a shape defined by a series of conical sections with the same axis, contiguous along it.
The polyconical solid G4BREPSolidPCone is a shape defined by a set of inner and outer conical or cylindrical surface sections and two planes perpendicular to the Z axis. Each conical surface is defined by its radius at two different planes perpendicular to the Z-axis. Inner and outer conical surfaces are defined using common Z planes.
G4BREPSolidPCone( const G4String& pName, G4double start_angle, G4double opening_angle, G4int num_z_planes, // sections, G4double z_start, const G4double z_values[], const G4double RMIN[], const G4double RMAX[] )The conical sections do not need to fill 360 degrees, but can have a common start and opening angle.
start_angle | starting angle |
opening_angle | opening angle |
num_z_planes | number of planes perpendicular to the z-axis used. |
z_start | starting value of z |
z_values | z coordinates of each plane |
RMIN | radius of inner cone at each plane |
RMAX | radius of outer cone at each plane |
The polygonal solid G4BREPSolidPolyhedra is a shape defined by an inner and outer polygonal surface and two planes perpendicular to the Z axis. Each polygonal surface is created by linking a series of polygons created at different planes perpendicular to the Z-axis. All these polygons all have the same number of sides (sides) and are defined at the same Z planes for both inner and outer polygonal surfaces.
The polygons do not need to fill 360 degrees, but have a start and opening angle.
The constructor takes the following parameters:
G4BREPSolidPolyhedra( const G4String& pName, G4double start_angle, G4double opening_angle, G4int sides, G4int num_z_planes, G4double z_start, const G4double z_values[], const G4double RMIN[], const G4double RMAX[] )which in addition to its name have the following meaning:
start_angle | starting angle |
opening_angle | opening angle |
sides | number of sides of each polygon in the x-y plane |
num_z_planes | number of planes perpendicular to the z-axis used. |
z_start | starting value of z |
z_values | z coordinates of each plane |
RMIN | radius of inner polygon at each corner |
RMAX | radius of outer polygon at each corner |
the shape is defined by the number of sides sides of the polygon in the plane perpendicular to the z-axis.
Figure 4.1.1 Example of geometries imported from CAD system and converted to tessellated solids. |
They can also be used to generate a solid bounded with a generic surface made
of planar facets. It is important that the supplied facets shall form a fully
enclose space to represent the solid.
Two types of facet can be used for the construction of a
G4TessellatedSolid: a triangular facet (G4TriangularFacet)
and a quadrangular facet (G4QuadrangularFacet).
An example on how to generate a simple tessellated shape is given below.
Example:
// First declare a tessellated solid // G4TessellatedSolid solidTarget = new G4TessellatedSolid("Solid_name"); // Define the facets which form the solid // G4double targetSize = 10*cm ; G4TriangularFacet *facet1 = new G4TriangularFacet (G4ThreeVector(-targetSize,-targetSize, 0.0), G4ThreeVector(+targetSize,-targetSize, 0.0), G4ThreeVector( 0.0, 0.0,+targetSize), ABSOLUTE); G4TriangularFacet *facet2 = new G4TriangularFacet (G4ThreeVector(+targetSize,-targetSize, 0.0), G4ThreeVector(+targetSize,+targetSize, 0.0), G4ThreeVector( 0.0, 0.0,+targetSize), ABSOLUTE); G4TriangularFacet *facet3 = new G4TriangularFacet (G4ThreeVector(+targetSize,+targetSize, 0.0), G4ThreeVector(-targetSize,+targetSize, 0.0), G4ThreeVector( 0.0, 0.0,+targetSize), ABSOLUTE); G4TriangularFacet *facet4 = new G4TriangularFacet (G4ThreeVector(-targetSize,+targetSize, 0.0), G4ThreeVector(-targetSize,-targetSize, 0.0), G4ThreeVector( 0.0, 0.0,+targetSize), ABSOLUTE); G4QuadrangularFacet *facet5 = new G4QuadrangularFacet (G4ThreeVector(-targetSize,-targetSize, 0.0), G4ThreeVector(-targetSize,+targetSize, 0.0), G4ThreeVector(+targetSize,+targetSize, 0.0), G4ThreeVector(+targetSize,-targetSize, 0.0), ABSOLUTE); // Now add the facets to the solid // solidTarget->AddFacet((G4VFacet*) facet1); solidTarget->AddFacet((G4VFacet*) facet2); solidTarget->AddFacet((G4VFacet*) facet3); solidTarget->AddFacet((G4VFacet*) facet4); solidTarget->AddFacet((G4VFacet*) facet5); Finally declare the solid is complete // solidTarget->SetSolidClosed(true); |
Source listing 4.1.1 An example of a simple tessellated solid with G4TessellatedSolid. |
The G4TriangularFacet class is used for the contruction of G4TessellatedSolid. It is defined by three vertices, which shall be supplied in anti-clockwise order looking from the outside of the solid where it belongs. Its constructor looks like:
G4TriangularFacet ( const G4ThreeVector Pt0, const G4ThreeVector vt1, const G4ThreeVector vt2, G4FacetVertexType fType )i.e., it takes 4 parameters to define the three vertices:
G4FacetVertexType | ABSOLUTE in which case Pt0, vt1 and vt2 are the three vertices in anti-clockwise order looking from the outside. |
G4FacetVertexType | RELATIVE in which case the first vertex is Pt0, the second vertex is Pt0+vt1 and the third vertex is Pt0+vt2, all in anti-clockwise order when looking from the outside. |
The G4QuadrangularFacet class can be used for the contruction of G4TessellatedSolid as well. It is defined by four vertices, which shall be in the same plane and be supplied in anti-clockwise order looking from the outside of the solid where it belongs. Its constructor looks like:
G4QuadrangularFacet ( const G4ThreeVector Pt0, const G4ThreeVector vt1, const G4ThreeVector vt2, const G4ThreeVector vt3, G4FacetVertexType fType )i.e., it takes 5 parameters to define the four vertices:
G4FacetVertexType | ABSOLUTE in which case Pt0, vt1, vt2 and vt3 are the four vertices required in anti-clockwise order when looking from the outside. |
G4FacetVertexType | RELATIVE in which case the first vertex is Pt0, the second vertex is Pt0+vt, the third vertex is Pt0+vt2 and the fourth vertex is Pt0+vt3, in anti-clockwise order when looking from the outside. |