#include <G4Tet.hh>
Inheritance diagram for G4Tet:
Definition at line 56 of file G4Tet.hh.
G4Tet::G4Tet | ( | const G4String & | pName, | |
G4ThreeVector | anchor, | |||
G4ThreeVector | p2, | |||
G4ThreeVector | p3, | |||
G4ThreeVector | p4, | |||
G4bool * | degeneracyFlag = 0 | |||
) |
Definition at line 89 of file G4Tet.cc.
References FatalException, and G4Exception().
Referenced by CheckDegeneracy(), and Clone().
00094 : G4VSolid(pName), fpPolyhedron(0), warningFlag(0) 00095 { 00096 // fV<x><y> is vector from vertex <y> to vertex <x> 00097 // 00098 G4ThreeVector fV21=p2-anchor; 00099 G4ThreeVector fV31=p3-anchor; 00100 G4ThreeVector fV41=p4-anchor; 00101 00102 // make sure this is a correctly oriented set of points for the tetrahedron 00103 // 00104 G4double signed_vol=fV21.cross(fV31).dot(fV41); 00105 if(signed_vol<0.0) 00106 { 00107 G4ThreeVector temp(p4); 00108 p4=p3; 00109 p3=temp; 00110 temp=fV41; 00111 fV41=fV31; 00112 fV31=temp; 00113 } 00114 fCubicVolume = std::fabs(signed_vol) / 6.; 00115 00116 G4ThreeVector fV24=p2-p4; 00117 G4ThreeVector fV43=p4-p3; 00118 G4ThreeVector fV32=p3-p2; 00119 00120 fXMin=std::min(std::min(std::min(anchor.x(), p2.x()),p3.x()),p4.x()); 00121 fXMax=std::max(std::max(std::max(anchor.x(), p2.x()),p3.x()),p4.x()); 00122 fYMin=std::min(std::min(std::min(anchor.y(), p2.y()),p3.y()),p4.y()); 00123 fYMax=std::max(std::max(std::max(anchor.y(), p2.y()),p3.y()),p4.y()); 00124 fZMin=std::min(std::min(std::min(anchor.z(), p2.z()),p3.z()),p4.z()); 00125 fZMax=std::max(std::max(std::max(anchor.z(), p2.z()),p3.z()),p4.z()); 00126 00127 fDx=(fXMax-fXMin)*0.5; fDy=(fYMax-fYMin)*0.5; fDz=(fZMax-fZMin)*0.5; 00128 00129 fMiddle=G4ThreeVector(fXMax+fXMin, fYMax+fYMin, fZMax+fZMin)*0.5; 00130 fMaxSize=std::max(std::max(std::max((anchor-fMiddle).mag(), 00131 (p2-fMiddle).mag()), 00132 (p3-fMiddle).mag()), 00133 (p4-fMiddle).mag()); 00134 00135 G4bool degenerate=std::fabs(signed_vol) < 1e-9*fMaxSize*fMaxSize*fMaxSize; 00136 00137 if(degeneracyFlag) *degeneracyFlag=degenerate; 00138 else if (degenerate) 00139 { 00140 G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException, 00141 "Degenerate tetrahedron not allowed."); 00142 } 00143 00144 fTol=1e-9*(std::fabs(fXMin)+std::fabs(fXMax)+std::fabs(fYMin) 00145 +std::fabs(fYMax)+std::fabs(fZMin)+std::fabs(fZMax)); 00146 //fTol=kCarTolerance; 00147 00148 fAnchor=anchor; 00149 fP2=p2; 00150 fP3=p3; 00151 fP4=p4; 00152 00153 G4ThreeVector fCenter123=(anchor+p2+p3)*(1.0/3.0); // face center 00154 G4ThreeVector fCenter134=(anchor+p4+p3)*(1.0/3.0); 00155 G4ThreeVector fCenter142=(anchor+p4+p2)*(1.0/3.0); 00156 G4ThreeVector fCenter234=(p2+p3+p4)*(1.0/3.0); 00157 00158 // compute area of each triangular face by cross product 00159 // and sum for total surface area 00160 00161 G4ThreeVector normal123=fV31.cross(fV21); 00162 G4ThreeVector normal134=fV41.cross(fV31); 00163 G4ThreeVector normal142=fV21.cross(fV41); 00164 G4ThreeVector normal234=fV32.cross(fV43); 00165 00166 fSurfaceArea=( 00167 normal123.mag()+ 00168 normal134.mag()+ 00169 normal142.mag()+ 00170 normal234.mag() 00171 )/2.0; 00172 00173 fNormal123=normal123.unit(); 00174 fNormal134=normal134.unit(); 00175 fNormal142=normal142.unit(); 00176 fNormal234=normal234.unit(); 00177 00178 fCdotN123=fCenter123.dot(fNormal123); 00179 fCdotN134=fCenter134.dot(fNormal134); 00180 fCdotN142=fCenter142.dot(fNormal142); 00181 fCdotN234=fCenter234.dot(fNormal234); 00182 }
G4Tet::~G4Tet | ( | ) | [virtual] |
G4Tet::G4Tet | ( | __void__ & | ) |
Definition at line 189 of file G4Tet.cc.
00190 : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0), 00191 fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0), 00192 fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0), 00193 fNormal234(0,0,0), warningFlag(0), 00194 fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.), 00195 fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.), 00196 fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.) 00197 { 00198 }
G4Tet::G4Tet | ( | const G4Tet & | rhs | ) |
Definition at line 213 of file G4Tet.cc.
00214 : G4VSolid(rhs), 00215 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea), 00216 fpPolyhedron(0), fAnchor(rhs.fAnchor), 00217 fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle), 00218 fNormal123(rhs.fNormal123), fNormal142(rhs.fNormal142), 00219 fNormal134(rhs.fNormal134), fNormal234(rhs.fNormal234), 00220 warningFlag(rhs.warningFlag), fCdotN123(rhs.fCdotN123), 00221 fCdotN142(rhs.fCdotN142), fCdotN134(rhs.fCdotN134), 00222 fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax), 00223 fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax), 00224 fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol), 00225 fMaxSize(rhs.fMaxSize) 00226 { 00227 }
G4bool G4Tet::CalculateExtent | ( | const EAxis | pAxis, | |
const G4VoxelLimits & | pVoxelLimit, | |||
const G4AffineTransform & | pTransform, | |||
G4double & | pmin, | |||
G4double & | pmax | |||
) | const [virtual] |
Implements G4VSolid.
Definition at line 291 of file G4Tet.cc.
References G4VoxelLimits::GetMaxXExtent(), G4VoxelLimits::GetMaxYExtent(), G4VoxelLimits::GetMaxZExtent(), G4VoxelLimits::GetMinXExtent(), G4VoxelLimits::GetMinYExtent(), G4VoxelLimits::GetMinZExtent(), G4AffineTransform::IsRotated(), G4VoxelLimits::IsXLimited(), G4VoxelLimits::IsYLimited(), G4VoxelLimits::IsZLimited(), kXAxis, kYAxis, kZAxis, G4AffineTransform::NetTranslation(), and G4AffineTransform::TransformPoint().
00295 { 00296 G4double xMin,xMax; 00297 G4double yMin,yMax; 00298 G4double zMin,zMax; 00299 00300 if (pTransform.IsRotated()) 00301 { 00302 G4ThreeVector pp0=pTransform.TransformPoint(fAnchor); 00303 G4ThreeVector pp1=pTransform.TransformPoint(fP2); 00304 G4ThreeVector pp2=pTransform.TransformPoint(fP3); 00305 G4ThreeVector pp3=pTransform.TransformPoint(fP4); 00306 00307 xMin = std::min(std::min(std::min(pp0.x(), pp1.x()),pp2.x()),pp3.x()); 00308 xMax = std::max(std::max(std::max(pp0.x(), pp1.x()),pp2.x()),pp3.x()); 00309 yMin = std::min(std::min(std::min(pp0.y(), pp1.y()),pp2.y()),pp3.y()); 00310 yMax = std::max(std::max(std::max(pp0.y(), pp1.y()),pp2.y()),pp3.y()); 00311 zMin = std::min(std::min(std::min(pp0.z(), pp1.z()),pp2.z()),pp3.z()); 00312 zMax = std::max(std::max(std::max(pp0.z(), pp1.z()),pp2.z()),pp3.z()); 00313 00314 } 00315 else 00316 { 00317 G4double xoffset = pTransform.NetTranslation().x() ; 00318 xMin = xoffset + fXMin; 00319 xMax = xoffset + fXMax; 00320 G4double yoffset = pTransform.NetTranslation().y() ; 00321 yMin = yoffset + fYMin; 00322 yMax = yoffset + fYMax; 00323 G4double zoffset = pTransform.NetTranslation().z() ; 00324 zMin = zoffset + fZMin; 00325 zMax = zoffset + fZMax; 00326 } 00327 00328 if (pVoxelLimit.IsXLimited()) 00329 { 00330 if ( (xMin > pVoxelLimit.GetMaxXExtent()+fTol) || 00331 (xMax < pVoxelLimit.GetMinXExtent()-fTol) ) { return false; } 00332 else 00333 { 00334 xMin = std::max(xMin, pVoxelLimit.GetMinXExtent()); 00335 xMax = std::min(xMax, pVoxelLimit.GetMaxXExtent()); 00336 } 00337 } 00338 00339 if (pVoxelLimit.IsYLimited()) 00340 { 00341 if ( (yMin > pVoxelLimit.GetMaxYExtent()+fTol) || 00342 (yMax < pVoxelLimit.GetMinYExtent()-fTol) ) { return false; } 00343 else 00344 { 00345 yMin = std::max(yMin, pVoxelLimit.GetMinYExtent()); 00346 yMax = std::min(yMax, pVoxelLimit.GetMaxYExtent()); 00347 } 00348 } 00349 00350 if (pVoxelLimit.IsZLimited()) 00351 { 00352 if ( (zMin > pVoxelLimit.GetMaxZExtent()+fTol) || 00353 (zMax < pVoxelLimit.GetMinZExtent()-fTol) ) { return false; } 00354 else 00355 { 00356 zMin = std::max(zMin, pVoxelLimit.GetMinZExtent()); 00357 zMax = std::min(zMax, pVoxelLimit.GetMaxZExtent()); 00358 } 00359 } 00360 00361 switch (pAxis) 00362 { 00363 case kXAxis: 00364 pMin=xMin; 00365 pMax=xMax; 00366 break; 00367 case kYAxis: 00368 pMin=yMin; 00369 pMax=yMax; 00370 break; 00371 case kZAxis: 00372 pMin=zMin; 00373 pMax=zMax; 00374 break; 00375 default: 00376 break; 00377 } 00378 00379 return true; 00380 }
G4bool G4Tet::CheckDegeneracy | ( | G4ThreeVector | anchor, | |
G4ThreeVector | p2, | |||
G4ThreeVector | p3, | |||
G4ThreeVector | p4 | |||
) | [static] |
G4VSolid * G4Tet::Clone | ( | ) | const [virtual] |
void G4Tet::ComputeDimensions | ( | G4VPVParameterisation * | p, | |
const G4int | n, | |||
const G4VPhysicalVolume * | pRep | |||
) | [virtual] |
G4NURBS * G4Tet::CreateNURBS | ( | ) | const [virtual] |
Reimplemented from G4VSolid.
Definition at line 852 of file G4Tet.cc.
00853 { 00854 return new G4NURBSbox (fDx, fDy, fDz); 00855 }
G4Polyhedron * G4Tet::CreatePolyhedron | ( | ) | const [virtual] |
Reimplemented from G4VSolid.
Definition at line 833 of file G4Tet.cc.
Referenced by GetPolyhedron().
00834 { 00835 G4Polyhedron *ph=new G4Polyhedron; 00836 G4double xyz[4][3]; 00837 const G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}}; 00838 xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z(); 00839 xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z(); 00840 xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z(); 00841 xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z(); 00842 00843 ph->createPolyhedron(4,4,xyz,faces); 00844 00845 return ph; 00846 }
G4ThreeVectorList * G4Tet::CreateRotatedVertices | ( | const G4AffineTransform & | pTransform | ) | const [protected] |
Definition at line 660 of file G4Tet.cc.
References G4VSolid::DumpInfo(), FatalException, G4Exception(), and G4AffineTransform::TransformPoint().
00661 { 00662 G4ThreeVectorList* vertices = new G4ThreeVectorList(); 00663 00664 if (vertices) 00665 { 00666 vertices->reserve(4); 00667 G4ThreeVector vertex0(fAnchor); 00668 G4ThreeVector vertex1(fP2); 00669 G4ThreeVector vertex2(fP3); 00670 G4ThreeVector vertex3(fP4); 00671 00672 vertices->push_back(pTransform.TransformPoint(vertex0)); 00673 vertices->push_back(pTransform.TransformPoint(vertex1)); 00674 vertices->push_back(pTransform.TransformPoint(vertex2)); 00675 vertices->push_back(pTransform.TransformPoint(vertex3)); 00676 } 00677 else 00678 { 00679 DumpInfo(); 00680 G4Exception("G4Tet::CreateRotatedVertices()", 00681 "GeomSolids0003", FatalException, 00682 "Error in allocation of vertices. Out of memory !"); 00683 } 00684 return vertices; 00685 }
const char* G4Tet::CVSFileVers | ( | ) | [inline] |
const char* G4Tet::CVSHeaderVers | ( | ) | [inline] |
void G4Tet::DescribeYourselfTo | ( | G4VGraphicsScene & | scene | ) | const [virtual] |
Implements G4VSolid.
Definition at line 815 of file G4Tet.cc.
References G4VGraphicsScene::AddSolid().
00816 { 00817 scene.AddSolid (*this); 00818 }
G4double G4Tet::DistanceToIn | ( | const G4ThreeVector & | p | ) | const [virtual] |
G4double G4Tet::DistanceToIn | ( | const G4ThreeVector & | p, | |
const G4ThreeVector & | v | |||
) | const [virtual] |
Implements G4VSolid.
Definition at line 477 of file G4Tet.cc.
00479 { 00480 G4ThreeVector vu(v.unit()), hp; 00481 G4double vdotn, t, tmin=kInfinity; 00482 00483 G4double extraDistance=10.0*fTol; // a little ways into the solid 00484 00485 vdotn=-vu.dot(fNormal123); 00486 if(vdotn > 1e-12) 00487 { // this is a candidate face, since it is pointing at us 00488 t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection 00489 if( (t>=-fTol) && (t<tmin) ) 00490 { // if not true, we're going away from this face or it's not close 00491 hp=p+vu*(t+extraDistance); // a little beyond point of intersection 00492 if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) && 00493 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) && 00494 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) ) 00495 { 00496 tmin=t; 00497 } 00498 } 00499 } 00500 00501 vdotn=-vu.dot(fNormal134); 00502 if(vdotn > 1e-12) 00503 { // # this is a candidate face, since it is pointing at us 00504 t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection 00505 if( (t>=-fTol) && (t<tmin) ) 00506 { // if not true, we're going away from this face 00507 hp=p+vu*(t+extraDistance); // a little beyond point of intersection 00508 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) && 00509 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) && 00510 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) ) 00511 { 00512 tmin=t; 00513 } 00514 } 00515 } 00516 00517 vdotn=-vu.dot(fNormal142); 00518 if(vdotn > 1e-12) 00519 { // # this is a candidate face, since it is pointing at us 00520 t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection 00521 if( (t>=-fTol) && (t<tmin) ) 00522 { // if not true, we're going away from this face 00523 hp=p+vu*(t+extraDistance); // a little beyond point of intersection 00524 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) && 00525 ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) && 00526 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) ) 00527 { 00528 tmin=t; 00529 } 00530 } 00531 } 00532 00533 vdotn=-vu.dot(fNormal234); 00534 if(vdotn > 1e-12) 00535 { // # this is a candidate face, since it is pointing at us 00536 t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection 00537 if( (t>=-fTol) && (t<tmin) ) 00538 { // if not true, we're going away from this face 00539 hp=p+vu*(t+extraDistance); // a little beyond point of intersection 00540 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) && 00541 ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) && 00542 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) ) 00543 { 00544 tmin=t; 00545 } 00546 } 00547 } 00548 00549 return std::max(0.0,tmin); 00550 }
G4double G4Tet::DistanceToOut | ( | const G4ThreeVector & | p | ) | const [virtual] |
Implements G4VSolid.
Definition at line 639 of file G4Tet.cc.
00640 { 00641 G4double t1,t2,t3,t4; 00642 t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside 00643 t2=fCdotN134-p.dot(fNormal134); // distance to plane 00644 t3=fCdotN142-p.dot(fNormal142); // distance to plane 00645 t4=fCdotN234-p.dot(fNormal234); // distance to plane 00646 00647 // if any one of these is negative, we are outside, 00648 // so return zero in that case 00649 00650 G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4); 00651 return (tmin < fTol)? 0:tmin; 00652 }
G4double G4Tet::DistanceToOut | ( | const G4ThreeVector & | p, | |
const G4ThreeVector & | v, | |||
const G4bool | calcNorm = false , |
|||
G4bool * | validNorm = 0 , |
|||
G4ThreeVector * | n = 0 | |||
) | const [virtual] |
Implements G4VSolid.
Definition at line 570 of file G4Tet.cc.
References G4VSolid::DumpInfo(), G4endl, G4Exception(), and JustWarning.
00573 { 00574 G4ThreeVector vu(v.unit()); 00575 G4double t1=kInfinity,t2=kInfinity,t3=kInfinity,t4=kInfinity, vdotn, tt; 00576 00577 vdotn=vu.dot(fNormal123); 00578 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate 00579 { 00580 t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection 00581 } 00582 00583 vdotn=vu.dot(fNormal134); 00584 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate 00585 { 00586 t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection 00587 } 00588 00589 vdotn=vu.dot(fNormal142); 00590 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate 00591 { 00592 t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection 00593 } 00594 00595 vdotn=vu.dot(fNormal234); 00596 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate 00597 { 00598 t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection 00599 } 00600 00601 tt=std::min(std::min(std::min(t1,t2),t3),t4); 00602 00603 if (warningFlag && (tt == kInfinity || tt < -fTol)) 00604 { 00605 DumpInfo(); 00606 std::ostringstream message; 00607 message << "No good intersection found or already outside!?" << G4endl 00608 << "p = " << p / mm << "mm" << G4endl 00609 << "v = " << v << G4endl 00610 << "t1, t2, t3, t4 (mm) " 00611 << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm; 00612 G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002", 00613 JustWarning, message); 00614 if(validNorm) 00615 { 00616 *validNorm=false; // flag normal as meaningless 00617 } 00618 } 00619 else if(calcNorm && n) 00620 { 00621 G4ThreeVector normal; 00622 if(tt==t1) { normal=fNormal123; } 00623 else if (tt==t2) { normal=fNormal134; } 00624 else if (tt==t3) { normal=fNormal142; } 00625 else if (tt==t4) { normal=fNormal234; } 00626 *n=normal; 00627 if(validNorm) { *validNorm=true; } 00628 } 00629 00630 return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit 00631 // if we are right on a face 00632 }
G4double G4Tet::GetCubicVolume | ( | ) | [virtual] |
G4GeometryType G4Tet::GetEntityType | ( | ) | const [virtual] |
G4VisExtent G4Tet::GetExtent | ( | ) | const [virtual] |
Reimplemented from G4VSolid.
Definition at line 824 of file G4Tet.cc.
00825 { 00826 return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax); 00827 }
G4ThreeVector G4Tet::GetPointOnSurface | ( | ) | const [virtual] |
Reimplemented from G4VSolid.
Definition at line 759 of file G4Tet.cc.
00760 { 00761 G4double chose,aOne,aTwo,aThree,aFour; 00762 G4ThreeVector p1, p2, p3, p4; 00763 00764 p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne); 00765 p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo); 00766 p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree); 00767 p4 = GetPointOnFace(fP4,fP3,fP2,aFour); 00768 00769 chose = RandFlat::shoot(0.,aOne+aTwo+aThree+aFour); 00770 if( (chose>=0.) && (chose <aOne) ) {return p1;} 00771 else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;} 00772 else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;} 00773 return p4; 00774 }
G4Polyhedron * G4Tet::GetPolyhedron | ( | ) | const [virtual] |
Reimplemented from G4VSolid.
Definition at line 861 of file G4Tet.cc.
References CreatePolyhedron(), and G4Polyhedron::GetNumberOfRotationStepsAtTimeOfCreation().
00862 { 00863 if (!fpPolyhedron || 00864 fpPolyhedron->GetNumberOfRotationStepsAtTimeOfCreation() != 00865 fpPolyhedron->GetNumberOfRotationSteps()) 00866 { 00867 delete fpPolyhedron; 00868 fpPolyhedron = CreatePolyhedron(); 00869 } 00870 return fpPolyhedron; 00871 }
G4double G4Tet::GetSurfaceArea | ( | ) | [virtual] |
std::vector< G4ThreeVector > G4Tet::GetVertices | ( | ) | const |
Definition at line 780 of file G4Tet.cc.
Referenced by G4GDMLWriteSolids::TetWrite().
00781 { 00782 std::vector<G4ThreeVector> vertices(4); 00783 vertices[0] = fAnchor; 00784 vertices[1] = fP2; 00785 vertices[2] = fP3; 00786 vertices[3] = fP4; 00787 00788 return vertices; 00789 }
EInside G4Tet::Inside | ( | const G4ThreeVector & | p | ) | const [virtual] |
Implements G4VSolid.
Definition at line 386 of file G4Tet.cc.
References kInside, kOutside, and kSurface.
00387 { 00388 G4double r123, r134, r142, r234; 00389 00390 // this is written to allow if-statement truncation so the outside test 00391 // (where most of the world is) can fail very quickly and efficiently 00392 00393 if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol || 00394 (r134=p.dot(fNormal134)-fCdotN134) > fTol || 00395 (r142=p.dot(fNormal142)-fCdotN142) > fTol || 00396 (r234=p.dot(fNormal234)-fCdotN234) > fTol ) 00397 { 00398 return kOutside; // at least one is out! 00399 } 00400 else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) ) 00401 { 00402 return kInside; // all are definitively inside 00403 } 00404 else 00405 { 00406 return kSurface; // too close to tell 00407 } 00408 }
Definition at line 234 of file G4Tet.cc.
References fAnchor, fCdotN123, fCdotN134, fCdotN142, fCdotN234, fCubicVolume, fDx, fDy, fDz, fMaxSize, fMiddle, fNormal123, fNormal134, fNormal142, fNormal234, fP2, fP3, fP4, fSurfaceArea, fTol, fXMax, fXMin, fYMax, fYMin, fZMax, fZMin, G4VSolid::operator=(), and warningFlag.
00235 { 00236 // Check assignment to self 00237 // 00238 if (this == &rhs) { return *this; } 00239 00240 // Copy base class data 00241 // 00242 G4VSolid::operator=(rhs); 00243 00244 // Copy data 00245 // 00246 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea; 00247 fpPolyhedron = 0; fAnchor = rhs.fAnchor; 00248 fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle; 00249 fNormal123 = rhs.fNormal123; fNormal142 = rhs.fNormal142; 00250 fNormal134 = rhs.fNormal134; fNormal234 = rhs.fNormal234; 00251 warningFlag = rhs.warningFlag; fCdotN123 = rhs.fCdotN123; 00252 fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134; 00253 fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax; 00254 fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax; 00255 fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol; 00256 fMaxSize = rhs.fMaxSize; 00257 00258 return *this; 00259 }
void G4Tet::PrintWarnings | ( | G4bool | flag | ) | [inline] |
std::ostream & G4Tet::StreamInfo | ( | std::ostream & | os | ) | const [virtual] |
Implements G4VSolid.
Definition at line 709 of file G4Tet.cc.
References G4VSolid::GetName().
00710 { 00711 G4int oldprc = os.precision(16); 00712 os << "-----------------------------------------------------------\n" 00713 << " *** Dump for solid - " << GetName() << " ***\n" 00714 << " ===================================================\n" 00715 << " Solid type: G4Tet\n" 00716 << " Parameters: \n" 00717 << " anchor: " << fAnchor/mm << " mm \n" 00718 << " p2: " << fP2/mm << " mm \n" 00719 << " p3: " << fP3/mm << " mm \n" 00720 << " p4: " << fP4/mm << " mm \n" 00721 << " normal123: " << fNormal123 << " \n" 00722 << " normal134: " << fNormal134 << " \n" 00723 << " normal142: " << fNormal142 << " \n" 00724 << " normal234: " << fNormal234 << " \n" 00725 << "-----------------------------------------------------------\n"; 00726 os.precision(oldprc); 00727 00728 return os; 00729 }
G4ThreeVector G4Tet::SurfaceNormal | ( | const G4ThreeVector & | p | ) | const [virtual] |
Implements G4VSolid.
Definition at line 417 of file G4Tet.cc.
References G4VSolid::kCarTolerance.
00418 { 00419 G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123); 00420 G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134); 00421 G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142); 00422 G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234); 00423 00424 static const G4double delta = 0.5*kCarTolerance; 00425 G4ThreeVector sumnorm(0., 0., 0.); 00426 G4int noSurfaces=0; 00427 00428 if (r123 <= delta) 00429 { 00430 noSurfaces ++; 00431 sumnorm= fNormal123; 00432 } 00433 00434 if (r134 <= delta) 00435 { 00436 noSurfaces ++; 00437 sumnorm += fNormal134; 00438 } 00439 00440 if (r142 <= delta) 00441 { 00442 noSurfaces ++; 00443 sumnorm += fNormal142; 00444 } 00445 if (r234 <= delta) 00446 { 00447 noSurfaces ++; 00448 sumnorm += fNormal234; 00449 } 00450 00451 if( noSurfaces > 0 ) 00452 { 00453 if( noSurfaces == 1 ) 00454 { 00455 return sumnorm; 00456 } 00457 else 00458 { 00459 return sumnorm.unit(); 00460 } 00461 } 00462 else // Approximative Surface Normal 00463 { 00464 00465 if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; } 00466 else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; } 00467 else if (r142 <= r234) { return fNormal142; } 00468 return fNormal234; 00469 } 00470 }