19#include "G4AffineTransform.hh"
20#include "G4Polyhedron.hh"
21#include "G4VGraphicsScene.hh"
22#include "G4VPVParameterisation.hh"
23#include "G4Version.hh"
24#include "G4VoxelLimits.hh"
25#if G4VERSION_NUMBER < 1000
27#include "G4NURBSbox.hh"
30#include "G4VisExtent.hh"
38#include <G4SystemOfUnits.hh>
43static const Double_t
gCm = 1. / cm;
54 const G4VoxelLimits& ,
const G4AffineTransform& ,
55 G4double& , G4double& )
const
60 G4cout <<
"Warning: TG4RootSolid::CalculateExtent() not implemented"
75 Bool_t in =
fShape->Contains(pt);
78 G4double safety =
fShape->Safety(pt, in) * cm;
79 if (TMath::Abs(safety) < 0.5 * kCarTolerance)
return kSurface;
80 if (in)
return kInside;
89 Double_t pt[3], dir[3], norm[3];
96 fShape->ComputeNormal(pt, dir, norm);
97 pt[0] += norm[0] * 2. * kCarTolerance;
98 pt[1] += norm[1] * 2. * kCarTolerance;
99 pt[2] += norm[2] * 2. * kCarTolerance;
101 G4ThreeVector n(norm[0], norm[1], norm[2]);
102 Bool_t in = fShape->Contains(pt);
104 n.set(-norm[0], -norm[1], -norm[2]);
110 const G4ThreeVector& p,
const G4ThreeVector& v)
const
117 Double_t pt[3], dir[3];
124 G4double dist =
fShape->DistFromOutside(pt, dir, 3) * cm;
125 if (dist < TGeoShape::Big())
return dist;
138 G4double safety =
fShape->Safety(pt, kFALSE) * cm;
144 const G4ThreeVector& v,
const G4bool calcNorm, G4bool* validNorm,
145 G4ThreeVector* n)
const
162 Double_t pt[3], dir[3], norm[3];
169 G4double dist =
fShape->DistFromInside(pt, dir, 3) * cm;
170 if (calcNorm) *validNorm =
true;
171 if (dist < 0.5 * kCarTolerance)
174 pt[0] += dist * dir[0];
175 pt[1] += dist * dir[1];
176 pt[2] += dist * dir[2];
179 fShape->ComputeNormal(pt, dir, norm);
180 *n = G4ThreeVector(norm[0], norm[1], norm[2]);
194 G4double safety =
fShape->Safety(pt, kTRUE) * cm;
200 const G4int ,
const G4VPhysicalVolume* )
204 G4cout <<
"Warning: TG4RootSolid::ComputeDimensions() not implemented"
216 G4double capacity =
fShape->Capacity() * cm3;
225 return G4String(
fShape->ClassName());
232 G4cout <<
"Warning: TG4RootSolid::GetPointOnSurface() not implemented"
234 return G4ThreeVector(0., 0., 0.);
241 os <<
"-----------------------------------------------------------\n"
242 <<
" *** Dump for solid - " << GetName() <<
" ***\n"
243 <<
" ===================================================\n"
244 <<
" Solid type: ROOT solid / " <<
fShape->ClassName() <<
"\n"
245 <<
" Bounding box: \n"
246 <<
" half length X: " << ((TGeoBBox*)
fShape)->GetDX() * cm / mm
248 <<
" half length Y: " << ((TGeoBBox*)
fShape)->GetDY() * cm / mm
250 <<
" half length Z: " << ((TGeoBBox*)
fShape)->GetDZ() * cm / mm
252 <<
"-----------------------------------------------------------\n";
264 scene.AddSolid(*
this);
271 G4double dx = ((TGeoBBox*)
fShape)->GetDX() * cm;
272 G4double dy = ((TGeoBBox*)
fShape)->GetDY() * cm;
273 G4double dz = ((TGeoBBox*)
fShape)->GetDZ() * cm;
274 const Double_t* origin = ((TGeoBBox*)
fShape)->GetOrigin();
275 G4double ox = origin[0] * cm;
276 G4double oy = origin[1] * cm;
277 G4double oz = origin[2] * cm;
278 return G4VisExtent(-dx + ox, dx + ox, -dy + oy, dy + oy, -dz + oz, dz + oz);
285 G4double dx = ((TGeoBBox*)
fShape)->GetDX() * cm;
286 G4double dy = ((TGeoBBox*)
fShape)->GetDY() * cm;
287 G4double dz = ((TGeoBBox*)
fShape)->GetDZ() * cm;
288 return new G4PolyhedronBox(dx, dy, dz);
291#if G4VERSION_NUMBER < 1000
static const double gCm
constant for conversion cm <-> mm
static const Double_t gCm
constant for conversion cm <-> mm
Definition of the TG4RootSolid class.
virtual G4GeometryType GetEntityType() const
virtual const G4VSolid * GetConstituentSolid(G4int no) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual const G4DisplacedSolid * GetDisplacedSolidPtr() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
virtual std::ostream & StreamInfo(std::ostream &os) const
virtual G4Polyhedron * GetPolyhedron() const
TGeoShape * fShape
TGeo associated shape.
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
virtual G4double GetCubicVolume()
virtual void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
TG4RootSolid()
Default ctor.
virtual G4NURBS * CreateNURBS() const
virtual G4ThreeVector GetPointOnSurface() const
virtual G4VisExtent GetExtent() const
virtual G4Polyhedron * CreatePolyhedron() const
virtual EInside Inside(const G4ThreeVector &p) const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const