Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GenericTrap Class Reference

#include <G4GenericTrap.hh>

+ Inheritance diagram for G4GenericTrap:

Public Member Functions

 G4GenericTrap (const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
 
 ~G4GenericTrap ()
 
G4double GetZHalfLength () const
 
G4int GetNofVertices () const
 
G4TwoVector GetVertex (G4int index) const
 
const std::vector< G4TwoVector > & GetVertices () const
 
G4double GetTwistAngle (G4int index) const
 
G4bool IsTwisted () const
 
G4int GetVisSubdivisions () const
 
void SetVisSubdivisions (G4int subdiv)
 
EInside Inside (const G4ThreeVector &p) const
 
G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const
 
G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const
 
G4double DistanceToIn (const G4ThreeVector &p) const
 
G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
 
G4double DistanceToOut (const G4ThreeVector &p) const
 
G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
 
G4GeometryType GetEntityType () const
 
G4VSolidClone () const
 
std::ostream & StreamInfo (std::ostream &os) const
 
G4ThreeVector GetPointOnSurface () const
 
G4double GetCubicVolume ()
 
G4double GetSurfaceArea ()
 
G4PolyhedronGetPolyhedron () const
 
void DescribeYourselfTo (G4VGraphicsScene &scene) const
 
G4VisExtent GetExtent () const
 
G4PolyhedronCreatePolyhedron () const
 
G4NURBSCreateNURBS () const
 
 G4GenericTrap (__void__ &)
 
 G4GenericTrap (const G4GenericTrap &rhs)
 
G4GenericTrapoperator= (const G4GenericTrap &rhs)
 
- Public Member Functions inherited from G4VSolid
 G4VSolid (const G4String &name)
 
virtual ~G4VSolid ()
 
G4bool operator== (const G4VSolid &s) const
 
G4String GetName () const
 
void SetName (const G4String &name)
 
G4double GetTolerance () const
 
virtual G4bool CalculateExtent (const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
 
virtual EInside Inside (const G4ThreeVector &p) const =0
 
virtual G4ThreeVector SurfaceNormal (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p, const G4ThreeVector &v) const =0
 
virtual G4double DistanceToIn (const G4ThreeVector &p) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
 
virtual G4double DistanceToOut (const G4ThreeVector &p) const =0
 
virtual void ComputeDimensions (G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
 
virtual G4double GetCubicVolume ()
 
virtual G4double GetSurfaceArea ()
 
virtual G4GeometryType GetEntityType () const =0
 
virtual G4ThreeVector GetPointOnSurface () const
 
virtual G4VSolidClone () const
 
virtual std::ostream & StreamInfo (std::ostream &os) const =0
 
void DumpInfo () const
 
virtual void DescribeYourselfTo (G4VGraphicsScene &scene) const =0
 
virtual G4VisExtent GetExtent () const
 
virtual G4PolyhedronCreatePolyhedron () const
 
virtual G4NURBSCreateNURBS () const
 
virtual G4PolyhedronGetPolyhedron () const
 
virtual const G4VSolidGetConstituentSolid (G4int no) const
 
virtual G4VSolidGetConstituentSolid (G4int no)
 
virtual const G4DisplacedSolidGetDisplacedSolidPtr () const
 
virtual G4DisplacedSolidGetDisplacedSolidPtr ()
 
 G4VSolid (__void__ &)
 
 G4VSolid (const G4VSolid &rhs)
 
G4VSolidoperator= (const G4VSolid &rhs)
 

Protected Attributes

G4PolyhedronfpPolyhedron
 
- Protected Attributes inherited from G4VSolid
G4double kCarTolerance
 

Additional Inherited Members

- Protected Member Functions inherited from G4VSolid
void CalculateClippedPolygonExtent (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipCrossSection (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipBetweenSections (G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
 
void ClipPolygon (G4ThreeVectorList &pPolygon, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis) const
 
G4double EstimateCubicVolume (G4int nStat, G4double epsilon) const
 
G4double EstimateSurfaceArea (G4int nStat, G4double ell) const
 

Detailed Description

Definition at line 79 of file G4GenericTrap.hh.

Constructor & Destructor Documentation

◆ G4GenericTrap() [1/3]

G4GenericTrap::G4GenericTrap ( const G4String name,
G4double  halfZ,
const std::vector< G4TwoVector > &  vertices 
)

Definition at line 68 of file G4GenericTrap.cc.

70 : G4VSolid(name),
71 fpPolyhedron(0),
72 fDz(halfZ),
73 fVertices(),
74 fIsTwisted(false),
75 fTessellatedSolid(0),
76 fMinBBoxVector(G4ThreeVector(0,0,0)),
77 fMaxBBoxVector(G4ThreeVector(0,0,0)),
78 fVisSubdivisions(0),
79 fSurfaceArea(0.),
80 fCubicVolume(0.)
81
82{
83 // General constructor
84 const G4double min_length=5*1.e-6;
85 G4double length = 0.;
86 G4int k=0;
87 G4String errorDescription = "InvalidSetup in \" ";
88 errorDescription += name;
89 errorDescription += "\"";
90
91 // Check vertices size
92
93 if ( G4int(vertices.size()) != fgkNofVertices )
94 {
95 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
96 FatalErrorInArgument, "Number of vertices != 8");
97 }
98
99 // Check dZ
100 //
101 if (halfZ < kCarTolerance)
102 {
103 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids0002",
104 FatalErrorInArgument, "dZ is too small or negative");
105 }
106
107 // Check Ordering and Copy vertices
108 //
109 if(CheckOrder(vertices))
110 {
111 for (G4int i=0; i<fgkNofVertices; ++i) {fVertices.push_back(vertices[i]);}
112 }
113 else
114 {
115 for (G4int i=0; i <4; ++i) {fVertices.push_back(vertices[3-i]);}
116 for (G4int i=0; i <4; ++i) {fVertices.push_back(vertices[7-i]);}
117 }
118
119 // Check length of segments and Adjust
120 //
121 for (G4int j=0; j < 2; j++)
122 {
123 for (G4int i=1; i<4; ++i)
124 {
125 k = j*4+i;
126 length = (fVertices[k]-fVertices[k-1]).mag();
127 if ( ( length < min_length) && ( length > kCarTolerance ) )
128 {
129 std::ostringstream message;
130 message << "Length segment is too small." << G4endl
131 << "Distance between " << fVertices[k-1] << " and "
132 << fVertices[k] << " is only " << length << " mm !";
133 G4Exception("G4GenericTrap::G4GenericTrap()", "GeomSolids1001",
134 JustWarning, message, "Vertices will be collapsed.");
135 fVertices[k]=fVertices[k-1];
136 }
137 }
138 }
139
140 // Compute Twist
141 //
142 for( G4int i=0; i<4; i++) { fTwist[i]=0.; }
143 fIsTwisted = ComputeIsTwisted();
144
145 // Compute Bounding Box
146 //
147 ComputeBBox();
148
149 // If not twisted - create tessellated solid
150 // (an alternative implementation for testing)
151 //
152#ifdef G4TESS_TEST
153 if ( !fIsTwisted ) { fTessellatedSolid = CreateTessellatedSolid(); }
154#endif
155}
@ JustWarning
@ FatalErrorInArgument
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4Polyhedron * fpPolyhedron
G4double kCarTolerance
Definition: G4VSolid.hh:307
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ ~G4GenericTrap()

G4GenericTrap::~G4GenericTrap ( )

Definition at line 180 of file G4GenericTrap.cc.

181{
182 // Destructor
183 delete fTessellatedSolid;
184}

◆ G4GenericTrap() [2/3]

G4GenericTrap::G4GenericTrap ( __void__ &  a)

Definition at line 159 of file G4GenericTrap.cc.

160 : G4VSolid(a),
161 fpPolyhedron(0),
162 fDz(0.),
163 fVertices(),
164 fIsTwisted(false),
165 fTessellatedSolid(0),
166 fMinBBoxVector(G4ThreeVector(0,0,0)),
167 fMaxBBoxVector(G4ThreeVector(0,0,0)),
168 fVisSubdivisions(0),
169 fSurfaceArea(0.),
170 fCubicVolume(0.)
171{
172 // Fake default constructor - sets only member data and allocates memory
173 // for usage restricted to object persistency.
174
175 for (size_t i=0; i<4; ++i) { fTwist[i]=0.; }
176}

◆ G4GenericTrap() [3/3]

G4GenericTrap::G4GenericTrap ( const G4GenericTrap rhs)

Definition at line 188 of file G4GenericTrap.cc.

189 : G4VSolid(rhs),
190 fpPolyhedron(0), fDz(rhs.fDz), fVertices(rhs.fVertices),
191 fIsTwisted(rhs.fIsTwisted), fTessellatedSolid(0),
192 fMinBBoxVector(rhs.fMinBBoxVector), fMaxBBoxVector(rhs.fMaxBBoxVector),
193 fVisSubdivisions(rhs.fVisSubdivisions),
194 fSurfaceArea(rhs.fSurfaceArea), fCubicVolume(rhs.fCubicVolume)
195{
196 for (size_t i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
197#ifdef G4TESS_TEST
198 if (rhs.fTessellatedSolid && !fIsTwisted )
199 { fTessellatedSolid = CreateTessellatedSolid(); }
200#endif
201}

Member Function Documentation

◆ CalculateExtent()

G4bool G4GenericTrap::CalculateExtent ( const EAxis  pAxis,
const G4VoxelLimits pVoxelLimit,
const G4AffineTransform pTransform,
G4double pmin,
G4double pmax 
) const
virtual

Implements G4VSolid.

Definition at line 1169 of file G4GenericTrap.cc.

1173{
1174#ifdef G4TESS_TEST
1175 if ( fTessellatedSolid )
1176 {
1177 return fTessellatedSolid->CalculateExtent(pAxis, pVoxelLimit,
1178 pTransform, pMin, pMax);
1179 }
1180#endif
1181
1182 // Computes bounding vectors for a shape
1183 //
1184 G4double Dx,Dy;
1185 G4ThreeVector minVec = GetMinimumBBox();
1186 G4ThreeVector maxVec = GetMaximumBBox();
1187 Dx = 0.5*(maxVec.x()- minVec.x());
1188 Dy = 0.5*(maxVec.y()- minVec.y());
1189
1190 if (!pTransform.IsRotated())
1191 {
1192 // Special case handling for unrotated shapes
1193 // Compute x/y/z mins and maxs respecting limits, with early returns
1194 // if outside limits. Then switch() on pAxis
1195 //
1196 G4double xoffset,xMin,xMax;
1197 G4double yoffset,yMin,yMax;
1198 G4double zoffset,zMin,zMax;
1199
1200 xoffset=pTransform.NetTranslation().x();
1201 xMin=xoffset-Dx;
1202 xMax=xoffset+Dx;
1203 if (pVoxelLimit.IsXLimited())
1204 {
1205 if ( (xMin>pVoxelLimit.GetMaxXExtent()+kCarTolerance)
1206 || (xMax<pVoxelLimit.GetMinXExtent()-kCarTolerance) )
1207 {
1208 return false;
1209 }
1210 else
1211 {
1212 if (xMin<pVoxelLimit.GetMinXExtent())
1213 {
1214 xMin=pVoxelLimit.GetMinXExtent();
1215 }
1216 if (xMax>pVoxelLimit.GetMaxXExtent())
1217 {
1218 xMax=pVoxelLimit.GetMaxXExtent();
1219 }
1220 }
1221 }
1222
1223 yoffset=pTransform.NetTranslation().y();
1224 yMin=yoffset-Dy;
1225 yMax=yoffset+Dy;
1226 if (pVoxelLimit.IsYLimited())
1227 {
1228 if ( (yMin>pVoxelLimit.GetMaxYExtent()+kCarTolerance)
1229 || (yMax<pVoxelLimit.GetMinYExtent()-kCarTolerance) )
1230 {
1231 return false;
1232 }
1233 else
1234 {
1235 if (yMin<pVoxelLimit.GetMinYExtent())
1236 {
1237 yMin=pVoxelLimit.GetMinYExtent();
1238 }
1239 if (yMax>pVoxelLimit.GetMaxYExtent())
1240 {
1241 yMax=pVoxelLimit.GetMaxYExtent();
1242 }
1243 }
1244 }
1245
1246 zoffset=pTransform.NetTranslation().z();
1247 zMin=zoffset-fDz;
1248 zMax=zoffset+fDz;
1249 if (pVoxelLimit.IsZLimited())
1250 {
1251 if ( (zMin>pVoxelLimit.GetMaxZExtent()+kCarTolerance)
1252 || (zMax<pVoxelLimit.GetMinZExtent()-kCarTolerance) )
1253 {
1254 return false;
1255 }
1256 else
1257 {
1258 if (zMin<pVoxelLimit.GetMinZExtent())
1259 {
1260 zMin=pVoxelLimit.GetMinZExtent();
1261 }
1262 if (zMax>pVoxelLimit.GetMaxZExtent())
1263 {
1264 zMax=pVoxelLimit.GetMaxZExtent();
1265 }
1266 }
1267 }
1268
1269 switch (pAxis)
1270 {
1271 case kXAxis:
1272 pMin = xMin;
1273 pMax = xMax;
1274 break;
1275 case kYAxis:
1276 pMin = yMin;
1277 pMax = yMax;
1278 break;
1279 case kZAxis:
1280 pMin = zMin;
1281 pMax = zMax;
1282 break;
1283 default:
1284 break;
1285 }
1286 pMin-=kCarTolerance;
1287 pMax+=kCarTolerance;
1288
1289 return true;
1290 }
1291 else
1292 {
1293 // General rotated case - create and clip mesh to boundaries
1294
1295 G4bool existsAfterClip=false;
1296 G4ThreeVectorList *vertices;
1297
1298 pMin=+kInfinity;
1299 pMax=-kInfinity;
1300
1301 // Calculate rotated vertex coordinates
1302 //
1303 vertices=CreateRotatedVertices(pTransform);
1304 ClipCrossSection(vertices,0,pVoxelLimit,pAxis,pMin,pMax);
1305 ClipCrossSection(vertices,4,pVoxelLimit,pAxis,pMin,pMax);
1306 ClipBetweenSections(vertices,0,pVoxelLimit,pAxis,pMin,pMax);
1307
1308 if ( (pMin!=kInfinity) || (pMax!=-kInfinity) )
1309 {
1310 existsAfterClip=true;
1311
1312 // Add 2*tolerance to avoid precision troubles
1313 //
1314 pMin-=kCarTolerance;
1315 pMax+=kCarTolerance;
1316 }
1317 else
1318 {
1319 // Check for case where completely enveloping clipping volume.
1320 // If point inside then we are confident that the solid completely
1321 // envelopes the clipping volume. Hence set min/max extents according
1322 // to clipping volume extents along the specified axis.
1323 //
1324 G4ThreeVector clipCentre(
1325 (pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
1326 (pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
1327 (pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
1328
1329 if (Inside(pTransform.Inverse().TransformPoint(clipCentre))!=kOutside)
1330 {
1331 existsAfterClip=true;
1332 pMin=pVoxelLimit.GetMinExtent(pAxis);
1333 pMax=pVoxelLimit.GetMaxExtent(pAxis);
1334 }
1335 }
1336 delete vertices;
1337 return existsAfterClip;
1338 }
1339}
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
double z() const
double x() const
double y() const
G4bool IsRotated() const
G4AffineTransform Inverse() const
G4ThreeVector NetTranslation() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
EInside Inside(const G4ThreeVector &p) const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:376
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:345
G4double GetMinExtent(const EAxis pAxis) const
G4bool IsYLimited() const
G4double GetMinZExtent() const
G4bool IsXLimited() const
G4double GetMaxExtent(const EAxis pAxis) const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4double GetMaxXExtent() const
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54
@ kOutside
Definition: geomdefs.hh:58

◆ Clone()

G4VSolid * G4GenericTrap::Clone ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1396 of file G4GenericTrap.cc.

1397{
1398 return new G4GenericTrap(*this);
1399}

◆ CreateNURBS()

G4NURBS * G4GenericTrap::CreateNURBS ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 2219 of file G4GenericTrap.cc.

2220{
2221#ifdef G4TESS_TEST
2222 if ( fTessellatedSolid )
2223 {
2224 return fTessellatedSolid->CreateNURBS();
2225 }
2226#endif
2227
2228 // Computes bounding vectors for the shape
2229 //
2230 G4double Dx,Dy;
2231 G4ThreeVector minVec = GetMinimumBBox();
2232 G4ThreeVector maxVec = GetMaximumBBox();
2233 Dx = 0.5*(maxVec.x()- minVec.y());
2234 Dy = 0.5*(maxVec.y()- minVec.y());
2235
2236 return new G4NURBSbox (Dx, Dy, fDz);
2237}
virtual G4NURBS * CreateNURBS() const

◆ CreatePolyhedron()

G4Polyhedron * G4GenericTrap::CreatePolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 2123 of file G4GenericTrap.cc.

2124{
2125
2126#ifdef G4TESS_TEST
2127 if ( fTessellatedSolid )
2128 {
2129 return fTessellatedSolid->CreatePolyhedron();
2130 }
2131#endif
2132
2133 // Approximation of Twisted Side
2134 // Construct extra Points, if Twisted Side
2135 //
2136 G4PolyhedronArbitrary* polyhedron;
2137 size_t nVertices, nFacets;
2138
2139 G4int subdivisions=0;
2140 G4int i;
2141 if(fIsTwisted)
2142 {
2143 if ( GetVisSubdivisions()!= 0 )
2144 {
2145 subdivisions=GetVisSubdivisions();
2146 }
2147 else
2148 {
2149 // Estimation of Number of Subdivisions for smooth visualisation
2150 //
2151 G4double maxTwist=0.;
2152 for(i=0; i<4; i++)
2153 {
2154 if(GetTwistAngle(i)>maxTwist) { maxTwist=GetTwistAngle(i); }
2155 }
2156
2157 // Computes bounding vectors for the shape
2158 //
2159 G4double Dx,Dy;
2160 G4ThreeVector minVec = GetMinimumBBox();
2161 G4ThreeVector maxVec = GetMaximumBBox();
2162 Dx = 0.5*(maxVec.x()- minVec.y());
2163 Dy = 0.5*(maxVec.y()- minVec.y());
2164 if (Dy > Dx) { Dx=Dy; }
2165
2166 subdivisions=8*G4int(maxTwist/(Dx*Dx*Dx)*fDz);
2167 if (subdivisions<4) { subdivisions=4; }
2168 if (subdivisions>30) { subdivisions=30; }
2169 }
2170 }
2171 G4int sub4=4*subdivisions;
2172 nVertices = 8+subdivisions*4;
2173 nFacets = 6+subdivisions*4;
2174 G4double cf=1./(subdivisions+1);
2175 polyhedron = new G4PolyhedronArbitrary (nVertices, nFacets);
2176
2177 // Add Vertex
2178 //
2179 for (i=0;i<4;i++)
2180 {
2181 polyhedron->AddVertex(G4ThreeVector(fVertices[i].x(),
2182 fVertices[i].y(),-fDz));
2183 }
2184 for( i=0;i<subdivisions;i++)
2185 {
2186 for(G4int j=0;j<4;j++)
2187 {
2188 G4TwoVector u=fVertices[j]+cf*(i+1)*( fVertices[j+4]-fVertices[j]);
2189 polyhedron->AddVertex(G4ThreeVector(u.x(),u.y(),-fDz+cf*2*fDz*(i+1)));
2190 }
2191 }
2192 for (i=4;i<8;i++)
2193 {
2194 polyhedron->AddVertex(G4ThreeVector(fVertices[i].x(),
2195 fVertices[i].y(),fDz));
2196 }
2197
2198 // Add Facets
2199 //
2200 polyhedron->AddFacet(1,4,3,2); //Z-plane
2201 for (i=0;i<subdivisions+1;i++)
2202 {
2203 G4int is=i*4;
2204 polyhedron->AddFacet(5+is,8+is,4+is,1+is);
2205 polyhedron->AddFacet(8+is,7+is,3+is,4+is);
2206 polyhedron->AddFacet(7+is,6+is,2+is,3+is);
2207 polyhedron->AddFacet(6+is,5+is,1+is,2+is);
2208 }
2209 polyhedron->AddFacet(5+sub4,6+sub4,7+sub4,8+sub4); //Z-plane
2210
2211 polyhedron->SetReferences();
2212 polyhedron->InvertFacets();
2213
2214 return (G4Polyhedron*) polyhedron;
2215}
double x() const
double y() const
G4double GetTwistAngle(G4int index) const
G4int GetVisSubdivisions() const
void AddVertex(const G4ThreeVector v)
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
virtual G4Polyhedron * CreatePolyhedron() const

Referenced by GetPolyhedron().

◆ DescribeYourselfTo()

void G4GenericTrap::DescribeYourselfTo ( G4VGraphicsScene scene) const
virtual

Implements G4VSolid.

Definition at line 2086 of file G4GenericTrap.cc.

2087{
2088
2089#ifdef G4TESS_TEST
2090 if ( fTessellatedSolid )
2091 {
2092 return fTessellatedSolid->DescribeYourselfTo(scene);
2093 }
2094#endif
2095
2096 scene.AddSolid(*this);
2097}
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
virtual void AddSolid(const G4Box &)=0

◆ DistanceToIn() [1/2]

G4double G4GenericTrap::DistanceToIn ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 789 of file G4GenericTrap.cc.

790{
791 // Computes the closest distance from given point to this shape
792
793#ifdef G4TESS_TEST
794 if ( fTessellatedSolid )
795 {
796 return fTessellatedSolid->DistanceToIn(p);
797 }
798#endif
799
800 G4double safz = std::fabs(p.z())-fDz;
801 if(safz<0) { safz=0; }
802
803 G4int iseg;
804 G4double safe = safz;
805 G4double safxy = safz;
806
807 for (iseg=0; iseg<4; iseg++)
808 {
809 safxy = SafetyToFace(p,iseg);
810 if (safxy>safe) { safe=safxy; }
811 }
812
813 return safe;
814}
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const

◆ DistanceToIn() [2/2]

G4double G4GenericTrap::DistanceToIn ( const G4ThreeVector p,
const G4ThreeVector v 
) const
virtual

Implements G4VSolid.

Definition at line 719 of file G4GenericTrap.cc.

721{
722#ifdef G4TESS_TEST
723 if ( fTessellatedSolid )
724 {
725 return fTessellatedSolid->DistanceToIn(p, v);
726 }
727#endif
728
729 static const G4double halfCarTolerance=kCarTolerance*0.5;
730
731 G4double dist[5];
733
734 // Check lateral faces
735 //
736 G4int i;
737 for (i=0; i<4; i++)
738 {
739 dist[i]=DistToPlane(p, v, i);
740 }
741
742 // Check Z planes
743 //
744 dist[4]=kInfinity;
745 if (std::fabs(p.z())>fDz-halfCarTolerance)
746 {
747 if (v.z())
748 {
749 G4ThreeVector pt;
750 if (p.z()>0)
751 {
752 dist[4] = (fDz-p.z())/v.z();
753 }
754 else
755 {
756 dist[4] = (-fDz-p.z())/v.z();
757 }
758 if (dist[4]<-halfCarTolerance)
759 {
760 dist[4]=kInfinity;
761 }
762 else
763 {
764 if(dist[4]<halfCarTolerance)
765 {
766 if(p.z()>0) { n=G4ThreeVector(0,0,1); }
767 else { n=G4ThreeVector(0,0,-1); }
768 if (n.dot(v)<0) { dist[4]=0.; }
769 else { dist[4]=kInfinity; }
770 }
771 pt=p+dist[4]*v;
772 if (Inside(pt)==kOutside) { dist[4]=kInfinity; }
773 }
774 }
775 }
776 G4double distmin = dist[0];
777 for (i=1;i<5;i++)
778 {
779 if (dist[i] < distmin) { distmin = dist[i]; }
780 }
781
782 if (distmin<halfCarTolerance) { distmin=0.; }
783
784 return distmin;
785}

◆ DistanceToOut() [1/2]

G4double G4GenericTrap::DistanceToOut ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 1142 of file G4GenericTrap.cc.

1143{
1144
1145#ifdef G4TESS_TEST
1146 if ( fTessellatedSolid )
1147 {
1148 return fTessellatedSolid->DistanceToOut(p);
1149 }
1150#endif
1151
1152 G4double safz = fDz-std::fabs(p.z());
1153 if (safz<0) { safz = 0; }
1154
1155 G4double safe = safz;
1156 G4double safxy = safz;
1157
1158 for (G4int iseg=0; iseg<4; iseg++)
1159 {
1160 safxy = std::fabs(SafetyToFace(p,iseg));
1161 if (safxy < safe) { safe = safxy; }
1162 }
1163
1164 return safe;
1165}
virtual G4double DistanceToOut(const G4ThreeVector &p) const

◆ DistanceToOut() [2/2]

G4double G4GenericTrap::DistanceToOut ( const G4ThreeVector p,
const G4ThreeVector v,
const G4bool  calcNorm = false,
G4bool validNorm = 0,
G4ThreeVector n = 0 
) const
virtual

Implements G4VSolid.

Definition at line 895 of file G4GenericTrap.cc.

900{
901#ifdef G4TESS_TEST
902 if ( fTessellatedSolid )
903 {
904 return fTessellatedSolid->DistanceToOut(p, v, calcNorm, validNorm, n);
905 }
906#endif
907
908 static const G4double halfCarTolerance=kCarTolerance*0.5;
909
910 G4double distmin;
911 G4bool lateral_cross = false;
912 ESide side = kUndefined;
913
914 if (calcNorm) { *validNorm=true; } // All normals are valid
915
916 if (v.z() < 0)
917 {
918 distmin=(-fDz-p.z())/v.z();
919 if (calcNorm) { side=kMZ; *n=G4ThreeVector(0,0,-1); }
920 }
921 else
922 {
923 if (v.z() > 0)
924 {
925 distmin = (fDz-p.z())/v.z();
926 if (calcNorm) { side=kPZ; *n=G4ThreeVector(0,0,1); }
927 }
928 else { distmin = kInfinity; }
929 }
930
931 G4double dz2 =0.5/fDz;
932 G4double xa,xb,xc,xd;
933 G4double ya,yb,yc,yd;
934
935 for (G4int ipl=0; ipl<4; ipl++)
936 {
937 G4int j = (ipl+1)%4;
938 xa=fVertices[ipl].x();
939 ya=fVertices[ipl].y();
940 xb=fVertices[ipl+4].x();
941 yb=fVertices[ipl+4].y();
942 xc=fVertices[j].x();
943 yc=fVertices[j].y();
944 xd=fVertices[4+j].x();
945 yd=fVertices[4+j].y();
946
947 if ( ((std::fabs(xb-xd)+std::fabs(yb-yd))<halfCarTolerance)
948 || ((std::fabs(xa-xc)+std::fabs(ya-yc))<halfCarTolerance) )
949 {
950 G4double q=DistToTriangle(p,v,ipl) ;
951 if ( (q>=0) && (q<distmin) )
952 {
953 distmin=q;
954 lateral_cross=true;
955 side=ESide(ipl+1);
956 }
957 continue;
958 }
959 G4double tx1 =dz2*(xb-xa);
960 G4double ty1 =dz2*(yb-ya);
961 G4double tx2 =dz2*(xd-xc);
962 G4double ty2 =dz2*(yd-yc);
963 G4double dzp =fDz+p.z();
964 G4double xs1 =xa+tx1*dzp;
965 G4double ys1 =ya+ty1*dzp;
966 G4double xs2 =xc+tx2*dzp;
967 G4double ys2 =yc+ty2*dzp;
968 G4double dxs =xs2-xs1;
969 G4double dys =ys2-ys1;
970 G4double dtx =tx2-tx1;
971 G4double dty =ty2-ty1;
972 G4double a = (dtx*v.y()-dty*v.x()+(tx1*ty2-tx2*ty1)*v.z())*v.z();
973 G4double b = dxs*v.y()-dys*v.x()+(dtx*p.y()-dty*p.x()+ty2*xs1-ty1*xs2
974 + tx1*ys2-tx2*ys1)*v.z();
975 G4double c=dxs*p.y()-dys*p.x()+xs1*ys2-xs2*ys1;
976 G4double q=kInfinity;
977
978 if (std::fabs(a) < kCarTolerance)
979 {
980 if (std::fabs(b) < kCarTolerance) { continue; }
981 q=-c/b;
982
983 // Check for Point on the Surface
984 //
985 if ((q > -halfCarTolerance) && (q < distmin))
986 {
987 if (q < halfCarTolerance)
988 {
989 if (NormalToPlane(p,ipl).dot(v)<0.) { continue; }
990 }
991 distmin =q;
992 lateral_cross=true;
993 side=ESide(ipl+1);
994 }
995 continue;
996 }
997 G4double d=b*b-4*a*c;
998 if (d >= 0.)
999 {
1000 if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1001 else { q=0.5*(-b+std::sqrt(d))/a; }
1002
1003 // Check for Point on the Surface
1004 //
1005 if (q > -halfCarTolerance )
1006 {
1007 if (q < distmin)
1008 {
1009 if(q < halfCarTolerance)
1010 {
1011 if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1012 {
1013 if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1014 else { q=0.5*(-b-std::sqrt(d))/a; }
1015 if (( q > halfCarTolerance) && (q < distmin))
1016 {
1017 distmin=q;
1018 lateral_cross = true;
1019 side=ESide(ipl+1);
1020 }
1021 continue;
1022 }
1023 }
1024 distmin = q;
1025 lateral_cross = true;
1026 side=ESide(ipl+1);
1027 }
1028 }
1029 else
1030 {
1031 if (a > 0) { q=0.5*(-b+std::sqrt(d))/a; }
1032 else { q=0.5*(-b-std::sqrt(d))/a; }
1033
1034 // Check for Point on the Surface
1035 //
1036 if ((q > -halfCarTolerance) && (q < distmin))
1037 {
1038 if (q < halfCarTolerance)
1039 {
1040 if (NormalToPlane(p,ipl).dot(v)<0.) // Check second root
1041 {
1042 if (a > 0) { q=0.5*(-b-std::sqrt(d))/a; }
1043 else { q=0.5*(-b+std::sqrt(d))/a; }
1044 if ( ( q > halfCarTolerance) && (q < distmin) )
1045 {
1046 distmin=q;
1047 lateral_cross = true;
1048 side=ESide(ipl+1);
1049 }
1050 continue;
1051 }
1052 }
1053 distmin =q;
1054 lateral_cross = true;
1055 side=ESide(ipl+1);
1056 }
1057 }
1058 }
1059 }
1060 if (!lateral_cross) // Make sure that track crosses the top or bottom
1061 {
1062 if (distmin >= kInfinity) { distmin=kCarTolerance; }
1063 G4ThreeVector pt=p+distmin*v;
1064
1065 // Check if propagated point is in the polygon
1066 //
1067 G4int i=0;
1068 if (v.z()>0.) { i=4; }
1069 std::vector<G4TwoVector> xy;
1070 for ( G4int j=0; j<4; j++) { xy.push_back(fVertices[i+j]); }
1071
1072 // Check Inside
1073 //
1074 if (InsidePolygone(pt,xy)==kOutside)
1075 {
1076 if(calcNorm)
1077 {
1078 if (v.z()>0) {side= kPZ; *n = G4ThreeVector(0,0,1);}
1079 else { side=kMZ; *n = G4ThreeVector(0,0,-1);}
1080 }
1081 return 0.;
1082 }
1083 else
1084 {
1085 if(v.z()>0) {side=kPZ;}
1086 else {side=kMZ;}
1087 }
1088 }
1089
1090 if (calcNorm)
1091 {
1092 G4ThreeVector pt=p+v*distmin;
1093 switch (side)
1094 {
1095 case kXY0:
1096 *n=NormalToPlane(pt,0);
1097 break;
1098 case kXY1:
1099 *n=NormalToPlane(pt,1);
1100 break;
1101 case kXY2:
1102 *n=NormalToPlane(pt,2);
1103 break;
1104 case kXY3:
1105 *n=NormalToPlane(pt,3);
1106 break;
1107 case kPZ:
1108 *n=G4ThreeVector(0,0,1);
1109 break;
1110 case kMZ:
1111 *n=G4ThreeVector(0,0,-1);
1112 break;
1113 default:
1114 DumpInfo();
1115 std::ostringstream message;
1116 G4int oldprc = message.precision(16);
1117 message << "Undefined side for valid surface normal to solid." << G4endl
1118 << "Position:" << G4endl
1119 << " p.x() = " << p.x()/mm << " mm" << G4endl
1120 << " p.y() = " << p.y()/mm << " mm" << G4endl
1121 << " p.z() = " << p.z()/mm << " mm" << G4endl
1122 << "Direction:" << G4endl
1123 << " v.x() = " << v.x() << G4endl
1124 << " v.y() = " << v.y() << G4endl
1125 << " v.z() = " << v.z() << G4endl
1126 << "Proposed distance :" << G4endl
1127 << " distmin = " << distmin/mm << " mm";
1128 message.precision(oldprc);
1129 G4Exception("G4GenericTrap::DistanceToOut(p,v,..)",
1130 "GeomSolids1002", JustWarning, message);
1131 break;
1132 }
1133 }
1134
1135 if (distmin<halfCarTolerance) { distmin=0.; }
1136
1137 return distmin;
1138}
ESide
Definition: G4Cons.cc:68
double dot(const Hep3Vector &) const
void DumpInfo() const

◆ GetCubicVolume()

G4double G4GenericTrap::GetCubicVolume ( )
virtual

Reimplemented from G4VSolid.

Definition at line 1514 of file G4GenericTrap.cc.

1515{
1516 if(fCubicVolume != 0.) {;}
1517 else { fCubicVolume = G4VSolid::GetCubicVolume(); }
1518 return fCubicVolume;
1519}
virtual G4double GetCubicVolume()
Definition: G4VSolid.cc:188

◆ GetEntityType()

G4GeometryType G4GenericTrap::GetEntityType ( ) const
virtual

Implements G4VSolid.

Definition at line 1389 of file G4GenericTrap.cc.

1390{
1391 return G4String("G4GenericTrap");
1392}

Referenced by StreamInfo().

◆ GetExtent()

G4VisExtent G4GenericTrap::GetExtent ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 2101 of file G4GenericTrap.cc.

2102{
2103 // Computes bounding vectors for the shape
2104
2105#ifdef G4TESS_TEST
2106 if ( fTessellatedSolid )
2107 {
2108 return fTessellatedSolid->GetExtent();
2109 }
2110#endif
2111
2112 G4double Dx,Dy;
2113 G4ThreeVector minVec = GetMinimumBBox();
2114 G4ThreeVector maxVec = GetMaximumBBox();
2115 Dx = 0.5*(maxVec.x()- minVec.x());
2116 Dy = 0.5*(maxVec.y()- minVec.y());
2117
2118 return G4VisExtent (-Dx, Dx, -Dy, Dy, -fDz, fDz);
2119}
virtual G4VisExtent GetExtent() const

◆ GetNofVertices()

G4int G4GenericTrap::GetNofVertices ( ) const
inline

◆ GetPointOnSurface()

G4ThreeVector G4GenericTrap::GetPointOnSurface ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 1426 of file G4GenericTrap.cc.

1427{
1428
1429#ifdef G4TESS_TEST
1430 if ( fTessellatedSolid )
1431 {
1432 return fTessellatedSolid->GetPointOnSurface();
1433 }
1434#endif
1435
1436 G4ThreeVector point;
1437 G4TwoVector u,v,w;
1438 G4double rand,area,chose,cf,lambda0,lambda1,alfa,beta,zp;
1439 G4int ipl,j;
1440
1441 std::vector<G4ThreeVector> vertices;
1442 for (G4int i=0; i<4;i++)
1443 {
1444 vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz));
1445 }
1446 for (G4int i=4; i<8;i++)
1447 {
1448 vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),fDz));
1449 }
1450
1451 // Surface Area of Planes(only estimation for twisted)
1452 //
1453 G4double Surface0=GetFaceSurfaceArea(vertices[0],vertices[1],
1454 vertices[2],vertices[3]);//-fDz plane
1455 G4double Surface1=GetFaceSurfaceArea(vertices[0],vertices[1],
1456 vertices[5],vertices[4]);// Lat plane
1457 G4double Surface2=GetFaceSurfaceArea(vertices[3],vertices[0],
1458 vertices[4],vertices[7]);// Lat plane
1459 G4double Surface3=GetFaceSurfaceArea(vertices[2],vertices[3],
1460 vertices[7],vertices[6]);// Lat plane
1461 G4double Surface4=GetFaceSurfaceArea(vertices[2],vertices[1],
1462 vertices[5],vertices[6]);// Lat plane
1463 G4double Surface5=GetFaceSurfaceArea(vertices[4],vertices[5],
1464 vertices[6],vertices[7]);// fDz plane
1465 rand = G4UniformRand();
1466 area = Surface0+Surface1+Surface2+Surface3+Surface4+Surface5;
1467 chose = rand*area;
1468
1469 if ( ( chose < Surface0)
1470 || ( chose > (Surface0+Surface1+Surface2+Surface3+Surface4)) )
1471 { // fDz or -fDz Plane
1472 ipl = G4int(G4UniformRand()*4);
1473 j = (ipl+1)%4;
1474 if(chose < Surface0)
1475 {
1476 zp = -fDz;
1477 u = fVertices[ipl]; v = fVertices[j];
1478 w = fVertices[(ipl+3)%4];
1479 }
1480 else
1481 {
1482 zp = fDz;
1483 u = fVertices[ipl+4]; v = fVertices[j+4];
1484 w = fVertices[(ipl+3)%4+4];
1485 }
1486 alfa = G4UniformRand();
1487 beta = G4UniformRand();
1488 lambda1=alfa*beta;
1489 lambda0=alfa-lambda1;
1490 v = v-u;
1491 w = w-u;
1492 v = u+lambda0*v+lambda1*w;
1493 }
1494 else // Lateral Plane Twisted or Not
1495 {
1496 if (chose < Surface0+Surface1) { ipl=0; }
1497 else if (chose < Surface0+Surface1+Surface2) { ipl=1; }
1498 else if (chose < Surface0+Surface1+Surface2+Surface3) { ipl=2; }
1499 else { ipl=3; }
1500 j = (ipl+1)%4;
1501 zp = -fDz+G4UniformRand()*2*fDz;
1502 cf = 0.5*(fDz-zp)/fDz;
1503 u = fVertices[ipl+4]+cf*( fVertices[ipl]-fVertices[ipl+4]);
1504 v = fVertices[j+4]+cf*(fVertices[j]-fVertices[j+4]);
1505 v = u+(v-u)*G4UniformRand();
1506 }
1507 point=G4ThreeVector(v.x(),v.y(),zp);
1508
1509 return point;
1510}
#define G4UniformRand()
Definition: Randomize.hh:53
virtual G4ThreeVector GetPointOnSurface() const

◆ GetPolyhedron()

G4Polyhedron * G4GenericTrap::GetPolyhedron ( ) const
virtual

Reimplemented from G4VSolid.

Definition at line 2064 of file G4GenericTrap.cc.

2065{
2066
2067#ifdef G4TESS_TEST
2068 if ( fTessellatedSolid )
2069 {
2070 return fTessellatedSolid->GetPolyhedron();
2071 }
2072#endif
2073
2074 if ( (!fpPolyhedron)
2077 {
2078 delete fpPolyhedron;
2080 }
2081 return fpPolyhedron;
2082}
G4Polyhedron * CreatePolyhedron() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
virtual G4Polyhedron * GetPolyhedron() const
static G4int GetNumberOfRotationSteps()

◆ GetSurfaceArea()

G4double G4GenericTrap::GetSurfaceArea ( )
virtual

Reimplemented from G4VSolid.

Definition at line 1523 of file G4GenericTrap.cc.

1524{
1525 if(fSurfaceArea != 0.) {;}
1526 else
1527 {
1528 std::vector<G4ThreeVector> vertices;
1529 for (G4int i=0; i<4;i++)
1530 {
1531 vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),-fDz));
1532 }
1533 for (G4int i=4; i<8;i++)
1534 {
1535 vertices.push_back(G4ThreeVector(fVertices[i].x(),fVertices[i].y(),fDz));
1536 }
1537
1538 // Surface Area of Planes(only estimation for twisted)
1539 //
1540 G4double fSurface0=GetFaceSurfaceArea(vertices[0],vertices[1],
1541 vertices[2],vertices[3]);//-fDz plane
1542 G4double fSurface1=GetFaceSurfaceArea(vertices[0],vertices[1],
1543 vertices[5],vertices[4]);// Lat plane
1544 G4double fSurface2=GetFaceSurfaceArea(vertices[3],vertices[0],
1545 vertices[4],vertices[7]);// Lat plane
1546 G4double fSurface3=GetFaceSurfaceArea(vertices[2],vertices[3],
1547 vertices[7],vertices[6]);// Lat plane
1548 G4double fSurface4=GetFaceSurfaceArea(vertices[2],vertices[1],
1549 vertices[5],vertices[6]);// Lat plane
1550 G4double fSurface5=GetFaceSurfaceArea(vertices[4],vertices[5],
1551 vertices[6],vertices[7]);// fDz plane
1552
1553 // Total Surface Area
1554 //
1555 if(!fIsTwisted)
1556 {
1557 fSurfaceArea = fSurface0+fSurface1+fSurface2
1558 + fSurface3+fSurface4+fSurface5;
1559 }
1560 else
1561 {
1562 fSurfaceArea = G4VSolid::GetSurfaceArea();
1563 }
1564 }
1565 return fSurfaceArea;
1566}
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:248

◆ GetTwistAngle()

G4double G4GenericTrap::GetTwistAngle ( G4int  index) const
inline

Referenced by CreatePolyhedron(), and SurfaceNormal().

◆ GetVertex()

G4TwoVector G4GenericTrap::GetVertex ( G4int  index) const
inline

◆ GetVertices()

const std::vector< G4TwoVector > & G4GenericTrap::GetVertices ( ) const
inline

◆ GetVisSubdivisions()

G4int G4GenericTrap::GetVisSubdivisions ( ) const
inline

Referenced by CreatePolyhedron().

◆ GetZHalfLength()

G4double G4GenericTrap::GetZHalfLength ( ) const
inline

◆ Inside()

EInside G4GenericTrap::Inside ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 304 of file G4GenericTrap.cc.

305{
306 // Test if point is inside this shape
307
308#ifdef G4TESS_TEST
309 if ( fTessellatedSolid )
310 {
311 return fTessellatedSolid->Inside(p);
312 }
313#endif
314
315 static const G4double halfCarTolerance=kCarTolerance*0.5;
316 EInside innew=kOutside;
317 std::vector<G4TwoVector> xy;
318
319 if (std::fabs(p.z()) <= fDz+halfCarTolerance) // First check Z range
320 {
321 // Compute intersection between Z plane containing point and the shape
322 //
323 G4double cf = 0.5*(fDz-p.z())/fDz;
324 for (G4int i=0; i<4; i++)
325 {
326 xy.push_back(fVertices[i+4]+cf*( fVertices[i]-fVertices[i+4]));
327 }
328
329 innew=InsidePolygone(p,xy);
330
331 if( (innew==kInside) || (innew==kSurface) )
332 {
333 if(std::fabs(p.z()) > fDz-halfCarTolerance) { innew=kSurface; }
334 }
335 }
336 return innew;
337}
virtual EInside Inside(const G4ThreeVector &p) const
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58

Referenced by CalculateExtent(), and DistanceToIn().

◆ IsTwisted()

G4bool G4GenericTrap::IsTwisted ( ) const
inline

◆ operator=()

G4GenericTrap & G4GenericTrap::operator= ( const G4GenericTrap rhs)

Definition at line 205 of file G4GenericTrap.cc.

206{
207 // Check assignment to self
208 //
209 if (this == &rhs) { return *this; }
210
211 // Copy base class data
212 //
214
215 // Copy data
216 //
217 fpPolyhedron = 0; fDz = rhs.fDz; fVertices = rhs.fVertices;
218 fIsTwisted = rhs.fIsTwisted; fTessellatedSolid = 0;
219 fMinBBoxVector = rhs.fMinBBoxVector; fMaxBBoxVector = rhs.fMaxBBoxVector;
220 fVisSubdivisions = rhs.fVisSubdivisions;
221 fSurfaceArea = rhs.fSurfaceArea; fCubicVolume = rhs.fCubicVolume;
222
223 for (size_t i=0; i<4; ++i) { fTwist[i] = rhs.fTwist[i]; }
224#ifdef G4TESS_TEST
225 if (rhs.fTessellatedSolid && !fIsTwisted )
226 { delete fTessellatedSolid; fTessellatedSolid = CreateTessellatedSolid(); }
227#endif
228
229 return *this;
230}
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110

◆ SetVisSubdivisions()

void G4GenericTrap::SetVisSubdivisions ( G4int  subdiv)
inline

◆ StreamInfo()

std::ostream & G4GenericTrap::StreamInfo ( std::ostream &  os) const
virtual

Implements G4VSolid.

Definition at line 1403 of file G4GenericTrap.cc.

1404{
1405 G4int oldprc = os.precision(16);
1406 os << "-----------------------------------------------------------\n"
1407 << " *** Dump for solid - " << GetName() << " *** \n"
1408 << " =================================================== \n"
1409 << " Solid geometry type: " << GetEntityType() << G4endl
1410 << " half length Z: " << fDz/mm << " mm \n"
1411 << " list of vertices:\n";
1412
1413 for ( G4int i=0; i<fgkNofVertices; ++i )
1414 {
1415 os << std::setw(5) << "#" << i
1416 << " vx = " << fVertices[i].x()/mm << " mm"
1417 << " vy = " << fVertices[i].y()/mm << " mm" << G4endl;
1418 }
1419 os.precision(oldprc);
1420
1421 return os;
1422}
G4GeometryType GetEntityType() const
G4String GetName() const

◆ SurfaceNormal()

G4ThreeVector G4GenericTrap::SurfaceNormal ( const G4ThreeVector p) const
virtual

Implements G4VSolid.

Definition at line 341 of file G4GenericTrap.cc.

342{
343 // Calculate side nearest to p, and return normal
344 // If two sides are equidistant, sum of the Normal is returned
345
346#ifdef G4TESS_TEST
347 if ( fTessellatedSolid )
348 {
349 return fTessellatedSolid->SurfaceNormal(p);
350 }
351#endif
352
353 static const G4double halfCarTolerance=kCarTolerance*0.5;
354
355 G4ThreeVector lnorm, sumnorm(0.,0.,0.), apprnorm(0.,0.,1.),
356 p0, p1, p2, r1, r2, r3, r4;
357 G4int noSurfaces = 0;
358 G4double distxy,distz;
359 G4bool zPlusSide=false;
360
361 distz = fDz-std::fabs(p.z());
362 if (distz < halfCarTolerance)
363 {
364 if(p.z()>0)
365 {
366 zPlusSide=true;
367 sumnorm=G4ThreeVector(0,0,1);
368 }
369 else
370 {
371 sumnorm=G4ThreeVector(0,0,-1);
372 }
373 noSurfaces ++;
374 }
375
376 // Check lateral planes
377 //
378 std:: vector<G4TwoVector> vertices;
379 G4double cf = 0.5*(fDz-p.z())/fDz;
380 for (G4int i=0; i<4; i++)
381 {
382 vertices.push_back(fVertices[i+4]+cf*(fVertices[i]-fVertices[i+4]));
383 }
384
385 // Compute distance for lateral planes
386 //
387 for (G4int q=0; q<4; q++)
388 {
389 p0=G4ThreeVector(vertices[q].x(),vertices[q].y(),p.z());
390 if(zPlusSide)
391 {
392 p1=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
393 }
394 else
395 {
396 p1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
397 }
398 p2=G4ThreeVector(vertices[(q+1)%4].x(),vertices[(q+1)%4].y(),p.z());
399
400 // Collapsed vertices
401 //
402 if ( (p2-p0).mag2() < kCarTolerance )
403 {
404 if ( std::fabs(p.z()+fDz) > kCarTolerance )
405 {
406 p2=G4ThreeVector(fVertices[(q+1)%4].x(),fVertices[(q+1)%4].y(),-fDz);
407 }
408 else
409 {
410 p2=G4ThreeVector(fVertices[(q+1)%4+4].x(),fVertices[(q+1)%4+4].y(),fDz);
411 }
412 }
413 lnorm = (p1-p0).cross(p2-p0);
414 lnorm = lnorm.unit();
415 if(zPlusSide) { lnorm=-lnorm; }
416
417 // Adjust Normal for Twisted Surface
418 //
419 if ( (fIsTwisted) && (GetTwistAngle(q)!=0) )
420 {
421 G4double normP=(p2-p0).mag();
422 if(normP)
423 {
424 G4double proj=(p-p0).dot(p2-p0)/normP;
425 if(proj<0) { proj=0; }
426 if(proj>normP) { proj=normP; }
427 G4int j=(q+1)%4;
428 r1=G4ThreeVector(fVertices[q+4].x(),fVertices[q+4].y(),fDz);
429 r2=G4ThreeVector(fVertices[j+4].x(),fVertices[j+4].y(),fDz);
430 r3=G4ThreeVector(fVertices[q].x(),fVertices[q].y(),-fDz);
431 r4=G4ThreeVector(fVertices[j].x(),fVertices[j].y(),-fDz);
432 r1=r1+proj*(r2-r1)/normP;
433 r3=r3+proj*(r4-r3)/normP;
434 r2=r1-r3;
435 r4=r2.cross(p2-p0); r4=r4.unit();
436 lnorm=r4;
437 }
438 } // End if fIsTwisted
439
440 distxy=std::fabs((p0-p).dot(lnorm));
441 if ( distxy<halfCarTolerance )
442 {
443 noSurfaces ++;
444
445 // Negative sign for Normal is taken for Outside Normal
446 //
447 sumnorm=sumnorm+lnorm;
448 }
449
450 // For ApproxSurfaceNormal
451 //
452 if (distxy<distz)
453 {
454 distz=distxy;
455 apprnorm=lnorm;
456 }
457 } // End for loop
458
459 // Calculate final Normal, add Normal in the Corners and Touching Sides
460 //
461 if ( noSurfaces == 0 )
462 {
463 G4Exception("G4GenericTrap::SurfaceNormal(p)", "GeomSolids1002",
464 JustWarning, "Point p is not on surface !?" );
465 sumnorm=apprnorm;
466 // Add Approximative Surface Normal Calculation?
467 }
468 else if ( noSurfaces == 1 ) { sumnorm = sumnorm; }
469 else { sumnorm = sumnorm.unit(); }
470
471 return sumnorm ;
472}
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const

Member Data Documentation

◆ fpPolyhedron

G4Polyhedron* G4GenericTrap::fpPolyhedron
mutableprotected

Definition at line 192 of file G4GenericTrap.hh.

Referenced by GetPolyhedron(), and operator=().


The documentation for this class was generated from the following files: