Geant4  10.00.p02
G4Tet.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the intellectual property of the *
19 // * Vanderbilt University Free Electron Laser Center *
20 // * Vanderbilt University, Nashville, TN, USA *
21 // * Development supported by: *
22 // * United States MFEL program under grant FA9550-04-1-0045 *
23 // * and NASA under contract number NNG04CT05P *
24 // * Written by Marcus H. Mendenhall and Robert A. Weller. *
25 // * *
26 // * Contributed to the Geant4 Core, January, 2005. *
27 // * *
28 // ********************************************************************
29 //
30 // $Id: G4Tet.cc 81641 2014-06-04 09:11:38Z gcosmo $
31 //
32 // class G4Tet
33 //
34 // Implementation for G4Tet class
35 //
36 // History:
37 //
38 // 20040903 - Marcus Mendenhall, created G4Tet
39 // 20041101 - Marcus Mendenhall, optimized constant dot products with
40 // fCdotNijk values
41 // 20041101 - MHM removed tracking error by clipping DistanceToOut to 0
42 // for surface cases
43 // 20041101 - MHM many speed optimizations in if statements
44 // 20041101 - MHM changed vdotn comparisons to 1e-12 instead of 0.0 to
45 // avoid nearly-parallel problems
46 // 20041102 - MHM Added extra distance into solid to DistanceToIn(p,v)
47 // hit testing
48 // 20041102 - MHM added ability to check for degeneracy without throwing
49 // G4Exception
50 // 20041103 - MHM removed many unused variables from class
51 // 20040803 - Dionysios Anninos, added GetPointOnSurface() method
52 // 20061112 - MHM added code for G4VSolid GetSurfaceArea()
53 // 20100920 - Gabriele Cosmo added copy-ctor and operator=()
54 //
55 // --------------------------------------------------------------------
56 
57 #include "G4Tet.hh"
58 
59 #if !defined(G4GEOM_USE_UTET)
60 
61 const char G4Tet::CVSVers[]="$Id: G4Tet.cc 81641 2014-06-04 09:11:38Z gcosmo $";
62 
63 #include "G4VoxelLimits.hh"
64 #include "G4AffineTransform.hh"
65 
66 #include "G4VPVParameterisation.hh"
67 
68 #include "Randomize.hh"
69 
70 #include "G4VGraphicsScene.hh"
71 #include "G4Polyhedron.hh"
72 #include "G4VisExtent.hh"
73 
74 #include "G4ThreeVector.hh"
75 
76 #include <cmath>
77 
78 using namespace CLHEP;
79 
81 //
82 // Constructor - create a tetrahedron
83 // This class is implemented separately from general polyhedra,
84 // because the simplex geometry can be computed very quickly,
85 // which may become important in situations imported from mesh generators,
86 // in which a very large number of G4Tets are created.
87 // A Tet has all of its geometrical information precomputed
88 
89 G4Tet::G4Tet(const G4String& pName,
90  G4ThreeVector anchor,
91  G4ThreeVector p2,
92  G4ThreeVector p3,
93  G4ThreeVector p4, G4bool *degeneracyFlag)
94  : G4VSolid(pName), fpPolyhedron(0), warningFlag(0)
95 {
96  // fV<x><y> is vector from vertex <y> to vertex <x>
97  //
98  G4ThreeVector fV21=p2-anchor;
99  G4ThreeVector fV31=p3-anchor;
100  G4ThreeVector fV41=p4-anchor;
101 
102  // make sure this is a correctly oriented set of points for the tetrahedron
103  //
104  G4double signed_vol=fV21.cross(fV31).dot(fV41);
105  if(signed_vol<0.0)
106  {
107  G4ThreeVector temp(p4);
108  p4=p3;
109  p3=temp;
110  temp=fV41;
111  fV41=fV31;
112  fV31=temp;
113  }
114  fCubicVolume = std::fabs(signed_vol) / 6.;
115 
116  G4ThreeVector fV24=p2-p4;
117  G4ThreeVector fV43=p4-p3;
118  G4ThreeVector fV32=p3-p2;
119 
120  fXMin=std::min(std::min(std::min(anchor.x(), p2.x()),p3.x()),p4.x());
121  fXMax=std::max(std::max(std::max(anchor.x(), p2.x()),p3.x()),p4.x());
122  fYMin=std::min(std::min(std::min(anchor.y(), p2.y()),p3.y()),p4.y());
123  fYMax=std::max(std::max(std::max(anchor.y(), p2.y()),p3.y()),p4.y());
124  fZMin=std::min(std::min(std::min(anchor.z(), p2.z()),p3.z()),p4.z());
125  fZMax=std::max(std::max(std::max(anchor.z(), p2.z()),p3.z()),p4.z());
126 
127  fDx=(fXMax-fXMin)*0.5; fDy=(fYMax-fYMin)*0.5; fDz=(fZMax-fZMin)*0.5;
128 
129  fMiddle=G4ThreeVector(fXMax+fXMin, fYMax+fYMin, fZMax+fZMin)*0.5;
130  fMaxSize=std::max(std::max(std::max((anchor-fMiddle).mag(),
131  (p2-fMiddle).mag()),
132  (p3-fMiddle).mag()),
133  (p4-fMiddle).mag());
134 
135  G4bool degenerate=std::fabs(signed_vol) < 1e-9*fMaxSize*fMaxSize*fMaxSize;
136 
137  if(degeneracyFlag) *degeneracyFlag=degenerate;
138  else if (degenerate)
139  {
140  G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException,
141  "Degenerate tetrahedron not allowed.");
142  }
143 
144  fTol=1e-9*(std::fabs(fXMin)+std::fabs(fXMax)+std::fabs(fYMin)
145  +std::fabs(fYMax)+std::fabs(fZMin)+std::fabs(fZMax));
146  //fTol=kCarTolerance;
147 
148  fAnchor=anchor;
149  fP2=p2;
150  fP3=p3;
151  fP4=p4;
152 
153  G4ThreeVector fCenter123=(anchor+p2+p3)*(1.0/3.0); // face center
154  G4ThreeVector fCenter134=(anchor+p4+p3)*(1.0/3.0);
155  G4ThreeVector fCenter142=(anchor+p4+p2)*(1.0/3.0);
156  G4ThreeVector fCenter234=(p2+p3+p4)*(1.0/3.0);
157 
158  // compute area of each triangular face by cross product
159  // and sum for total surface area
160 
161  G4ThreeVector normal123=fV31.cross(fV21);
162  G4ThreeVector normal134=fV41.cross(fV31);
163  G4ThreeVector normal142=fV21.cross(fV41);
164  G4ThreeVector normal234=fV32.cross(fV43);
165 
166  fSurfaceArea=(
167  normal123.mag()+
168  normal134.mag()+
169  normal142.mag()+
170  normal234.mag()
171  )/2.0;
172 
173  fNormal123=normal123.unit();
174  fNormal134=normal134.unit();
175  fNormal142=normal142.unit();
176  fNormal234=normal234.unit();
177 
178  fCdotN123=fCenter123.dot(fNormal123);
179  fCdotN134=fCenter134.dot(fNormal134);
180  fCdotN142=fCenter142.dot(fNormal142);
181  fCdotN234=fCenter234.dot(fNormal234);
182 }
183 
185 //
186 // Fake default constructor - sets only member data and allocates memory
187 // for usage restricted to object persistency.
188 //
189 G4Tet::G4Tet( __void__& a )
190  : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0),
191  fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0),
192  fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0),
193  fNormal234(0,0,0), warningFlag(0),
194  fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.),
195  fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.),
196  fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.)
197 {
198 }
199 
201 //
202 // Destructor
203 
205 {
206  delete fpPolyhedron;
207 }
208 
210 //
211 // Copy constructor
212 
213 G4Tet::G4Tet(const G4Tet& rhs)
214  : G4VSolid(rhs),
215  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
216  fpPolyhedron(0), fAnchor(rhs.fAnchor),
217  fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle),
218  fNormal123(rhs.fNormal123), fNormal142(rhs.fNormal142),
219  fNormal134(rhs.fNormal134), fNormal234(rhs.fNormal234),
220  warningFlag(rhs.warningFlag), fCdotN123(rhs.fCdotN123),
221  fCdotN142(rhs.fCdotN142), fCdotN134(rhs.fCdotN134),
222  fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax),
223  fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax),
224  fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol),
225  fMaxSize(rhs.fMaxSize)
226 {
228 }
229 
230 
232 //
233 // Assignment operator
234 
236 {
237  // Check assignment to self
238  //
239  if (this == &rhs) { return *this; }
240 
241  // Copy base class data
242  //
243  G4VSolid::operator=(rhs);
244 
245  // Copy data
246  //
248  fAnchor = rhs.fAnchor;
249  fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
253  fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
254  fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
255  fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
256  fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
257  fMaxSize = rhs.fMaxSize;
259 
260  return *this;
261 }
262 
264 //
265 // CheckDegeneracy
266 
268  G4ThreeVector p2,
269  G4ThreeVector p3,
270  G4ThreeVector p4 )
271 {
272  G4bool result;
273  G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
274  delete object;
275  return result;
276 }
277 
279 //
280 // Dispatch to parameterisation for replication mechanism dimension
281 // computation & modification.
282 
284  const G4int ,
285  const G4VPhysicalVolume* )
286 {
287 }
288 
290 //
291 // Calculate extent under transform and specified limit
292 
294  const G4VoxelLimits& pVoxelLimit,
295  const G4AffineTransform& pTransform,
296  G4double& pMin, G4double& pMax) const
297 {
298  G4double xMin,xMax;
299  G4double yMin,yMax;
300  G4double zMin,zMax;
301 
302  if (pTransform.IsRotated())
303  {
304  G4ThreeVector pp0=pTransform.TransformPoint(fAnchor);
305  G4ThreeVector pp1=pTransform.TransformPoint(fP2);
306  G4ThreeVector pp2=pTransform.TransformPoint(fP3);
307  G4ThreeVector pp3=pTransform.TransformPoint(fP4);
308 
309  xMin = std::min(std::min(std::min(pp0.x(), pp1.x()),pp2.x()),pp3.x());
310  xMax = std::max(std::max(std::max(pp0.x(), pp1.x()),pp2.x()),pp3.x());
311  yMin = std::min(std::min(std::min(pp0.y(), pp1.y()),pp2.y()),pp3.y());
312  yMax = std::max(std::max(std::max(pp0.y(), pp1.y()),pp2.y()),pp3.y());
313  zMin = std::min(std::min(std::min(pp0.z(), pp1.z()),pp2.z()),pp3.z());
314  zMax = std::max(std::max(std::max(pp0.z(), pp1.z()),pp2.z()),pp3.z());
315 
316  }
317  else
318  {
319  G4double xoffset = pTransform.NetTranslation().x() ;
320  xMin = xoffset + fXMin;
321  xMax = xoffset + fXMax;
322  G4double yoffset = pTransform.NetTranslation().y() ;
323  yMin = yoffset + fYMin;
324  yMax = yoffset + fYMax;
325  G4double zoffset = pTransform.NetTranslation().z() ;
326  zMin = zoffset + fZMin;
327  zMax = zoffset + fZMax;
328  }
329 
330  if (pVoxelLimit.IsXLimited())
331  {
332  if ( (xMin > pVoxelLimit.GetMaxXExtent()+fTol) ||
333  (xMax < pVoxelLimit.GetMinXExtent()-fTol) ) { return false; }
334  else
335  {
336  xMin = std::max(xMin, pVoxelLimit.GetMinXExtent());
337  xMax = std::min(xMax, pVoxelLimit.GetMaxXExtent());
338  }
339  }
340 
341  if (pVoxelLimit.IsYLimited())
342  {
343  if ( (yMin > pVoxelLimit.GetMaxYExtent()+fTol) ||
344  (yMax < pVoxelLimit.GetMinYExtent()-fTol) ) { return false; }
345  else
346  {
347  yMin = std::max(yMin, pVoxelLimit.GetMinYExtent());
348  yMax = std::min(yMax, pVoxelLimit.GetMaxYExtent());
349  }
350  }
351 
352  if (pVoxelLimit.IsZLimited())
353  {
354  if ( (zMin > pVoxelLimit.GetMaxZExtent()+fTol) ||
355  (zMax < pVoxelLimit.GetMinZExtent()-fTol) ) { return false; }
356  else
357  {
358  zMin = std::max(zMin, pVoxelLimit.GetMinZExtent());
359  zMax = std::min(zMax, pVoxelLimit.GetMaxZExtent());
360  }
361  }
362 
363  switch (pAxis)
364  {
365  case kXAxis:
366  pMin=xMin;
367  pMax=xMax;
368  break;
369  case kYAxis:
370  pMin=yMin;
371  pMax=yMax;
372  break;
373  case kZAxis:
374  pMin=zMin;
375  pMax=zMax;
376  break;
377  default:
378  break;
379  }
380 
381  return true;
382 }
383 
385 //
386 // Return whether point inside/outside/on surface, using tolerance
387 
389 {
390  G4double r123, r134, r142, r234;
391 
392  // this is written to allow if-statement truncation so the outside test
393  // (where most of the world is) can fail very quickly and efficiently
394 
395  if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
396  (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
397  (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
398  (r234=p.dot(fNormal234)-fCdotN234) > fTol )
399  {
400  return kOutside; // at least one is out!
401  }
402  else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
403  {
404  return kInside; // all are definitively inside
405  }
406  else
407  {
408  return kSurface; // too close to tell
409  }
410 }
411 
413 //
414 // Calculate side nearest to p, and return normal
415 // If two sides are equidistant, normal of first side (x/y/z)
416 // encountered returned.
417 // This assumes that we are looking from the inside!
418 
420 {
421  G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
422  G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
423  G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
424  G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
425 
426  const G4double delta = 0.5*kCarTolerance;
427  G4ThreeVector sumnorm(0., 0., 0.);
428  G4int noSurfaces=0;
429 
430  if (r123 <= delta)
431  {
432  noSurfaces ++;
433  sumnorm= fNormal123;
434  }
435 
436  if (r134 <= delta)
437  {
438  noSurfaces ++;
439  sumnorm += fNormal134;
440  }
441 
442  if (r142 <= delta)
443  {
444  noSurfaces ++;
445  sumnorm += fNormal142;
446  }
447  if (r234 <= delta)
448  {
449  noSurfaces ++;
450  sumnorm += fNormal234;
451  }
452 
453  if( noSurfaces > 0 )
454  {
455  if( noSurfaces == 1 )
456  {
457  return sumnorm;
458  }
459  else
460  {
461  return sumnorm.unit();
462  }
463  }
464  else // Approximative Surface Normal
465  {
466 
467  if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
468  else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; }
469  else if (r142 <= r234) { return fNormal142; }
470  return fNormal234;
471  }
472 }
474 //
475 // Calculate distance to box from an outside point
476 // - return kInfinity if no intersection.
477 // All this is very unrolled, for speed.
478 
480  const G4ThreeVector& v) const
481 {
482  G4ThreeVector vu(v.unit()), hp;
483  G4double vdotn, t, tmin=kInfinity;
484 
485  G4double extraDistance=10.0*fTol; // a little ways into the solid
486 
487  vdotn=-vu.dot(fNormal123);
488  if(vdotn > 1e-12)
489  { // this is a candidate face, since it is pointing at us
490  t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection
491  if( (t>=-fTol) && (t<tmin) )
492  { // if not true, we're going away from this face or it's not close
493  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
494  if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
495  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
496  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
497  {
498  tmin=t;
499  }
500  }
501  }
502 
503  vdotn=-vu.dot(fNormal134);
504  if(vdotn > 1e-12)
505  { // # this is a candidate face, since it is pointing at us
506  t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection
507  if( (t>=-fTol) && (t<tmin) )
508  { // if not true, we're going away from this face
509  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
510  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
511  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
512  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
513  {
514  tmin=t;
515  }
516  }
517  }
518 
519  vdotn=-vu.dot(fNormal142);
520  if(vdotn > 1e-12)
521  { // # this is a candidate face, since it is pointing at us
522  t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection
523  if( (t>=-fTol) && (t<tmin) )
524  { // if not true, we're going away from this face
525  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
526  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
527  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
528  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
529  {
530  tmin=t;
531  }
532  }
533  }
534 
535  vdotn=-vu.dot(fNormal234);
536  if(vdotn > 1e-12)
537  { // # this is a candidate face, since it is pointing at us
538  t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection
539  if( (t>=-fTol) && (t<tmin) )
540  { // if not true, we're going away from this face
541  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
542  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
543  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
544  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
545  {
546  tmin=t;
547  }
548  }
549  }
550 
551  return std::max(0.0,tmin);
552 }
553 
555 //
556 // Approximate distance to tet.
557 // returns distance to sphere centered on bounding box
558 // - If inside return 0
559 
561 {
562  G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
563  return std::max(0.0, dd);
564 }
565 
567 //
568 // Calcluate distance to surface of box from inside
569 // by calculating distances to box's x/y/z planes.
570 // Smallest distance is exact distance to exiting.
571 
573  const G4bool calcNorm,
574  G4bool *validNorm, G4ThreeVector *n) const
575 {
576  G4ThreeVector vu(v.unit());
577  G4double t1=kInfinity,t2=kInfinity,t3=kInfinity,t4=kInfinity, vdotn, tt;
578 
579  vdotn=vu.dot(fNormal123);
580  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
581  {
582  t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection
583  }
584 
585  vdotn=vu.dot(fNormal134);
586  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
587  {
588  t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection
589  }
590 
591  vdotn=vu.dot(fNormal142);
592  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
593  {
594  t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection
595  }
596 
597  vdotn=vu.dot(fNormal234);
598  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
599  {
600  t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection
601  }
602 
603  tt=std::min(std::min(std::min(t1,t2),t3),t4);
604 
605  if (warningFlag && (tt == kInfinity || tt < -fTol))
606  {
607  DumpInfo();
608  std::ostringstream message;
609  message << "No good intersection found or already outside!?" << G4endl
610  << "p = " << p / mm << "mm" << G4endl
611  << "v = " << v << G4endl
612  << "t1, t2, t3, t4 (mm) "
613  << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
614  G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
615  JustWarning, message);
616  if(validNorm)
617  {
618  *validNorm=false; // flag normal as meaningless
619  }
620  }
621  else if(calcNorm && n)
622  {
624  if(tt==t1) { normal=fNormal123; }
625  else if (tt==t2) { normal=fNormal134; }
626  else if (tt==t3) { normal=fNormal142; }
627  else if (tt==t4) { normal=fNormal234; }
628  *n=normal;
629  if(validNorm) { *validNorm=true; }
630  }
631 
632  return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
633  // if we are right on a face
634 }
635 
637 //
638 // Calculate exact shortest distance to any boundary from inside
639 // - If outside return 0
640 
642 {
643  G4double t1,t2,t3,t4;
644  t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside
645  t2=fCdotN134-p.dot(fNormal134); // distance to plane
646  t3=fCdotN142-p.dot(fNormal142); // distance to plane
647  t4=fCdotN234-p.dot(fNormal234); // distance to plane
648 
649  // if any one of these is negative, we are outside,
650  // so return zero in that case
651 
652  G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
653  return (tmin < fTol)? 0:tmin;
654 }
655 
657 //
658 // Create a List containing the transformed vertices
659 // Note: Caller has deletion responsibility
660 
663 {
664  G4ThreeVectorList* vertices = new G4ThreeVectorList();
665 
666  if (vertices)
667  {
668  vertices->reserve(4);
669  G4ThreeVector vertex0(fAnchor);
670  G4ThreeVector vertex1(fP2);
671  G4ThreeVector vertex2(fP3);
672  G4ThreeVector vertex3(fP4);
673 
674  vertices->push_back(pTransform.TransformPoint(vertex0));
675  vertices->push_back(pTransform.TransformPoint(vertex1));
676  vertices->push_back(pTransform.TransformPoint(vertex2));
677  vertices->push_back(pTransform.TransformPoint(vertex3));
678  }
679  else
680  {
681  DumpInfo();
682  G4Exception("G4Tet::CreateRotatedVertices()",
683  "GeomSolids0003", FatalException,
684  "Error in allocation of vertices. Out of memory !");
685  }
686  return vertices;
687 }
688 
690 //
691 // GetEntityType
692 
694 {
695  return G4String("G4Tet");
696 }
697 
699 //
700 // Make a clone of the object
701 
703 {
704  return new G4Tet(*this);
705 }
706 
708 //
709 // Stream object contents to an output stream
710 
711 std::ostream& G4Tet::StreamInfo(std::ostream& os) const
712 {
713  G4int oldprc = os.precision(16);
714  os << "-----------------------------------------------------------\n"
715  << " *** Dump for solid - " << GetName() << " ***\n"
716  << " ===================================================\n"
717  << " Solid type: G4Tet\n"
718  << " Parameters: \n"
719  << " anchor: " << fAnchor/mm << " mm \n"
720  << " p2: " << fP2/mm << " mm \n"
721  << " p3: " << fP3/mm << " mm \n"
722  << " p4: " << fP4/mm << " mm \n"
723  << " normal123: " << fNormal123 << " \n"
724  << " normal134: " << fNormal134 << " \n"
725  << " normal142: " << fNormal142 << " \n"
726  << " normal234: " << fNormal234 << " \n"
727  << "-----------------------------------------------------------\n";
728  os.precision(oldprc);
729 
730  return os;
731 }
732 
733 
735 //
736 // GetPointOnFace
737 //
738 // Auxiliary method for get point on surface
739 
741  G4ThreeVector p3, G4double& area) const
742 {
743  G4double lambda1,lambda2;
744  G4ThreeVector v, w;
745 
746  v = p3 - p1;
747  w = p1 - p2;
748 
749  lambda1 = RandFlat::shoot(0.,1.);
750  lambda2 = RandFlat::shoot(0.,lambda1);
751 
752  area = 0.5*(v.cross(w)).mag();
753 
754  return (p2 + lambda1*w + lambda2*v);
755 }
756 
758 //
759 // GetPointOnSurface
760 
762 {
763  G4double chose,aOne,aTwo,aThree,aFour;
764  G4ThreeVector p1, p2, p3, p4;
765 
766  p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
767  p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
768  p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
769  p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
770 
771  chose = RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
772  if( (chose>=0.) && (chose <aOne) ) {return p1;}
773  else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
774  else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
775  return p4;
776 }
777 
779 //
780 // GetVertices
781 
782 std::vector<G4ThreeVector> G4Tet::GetVertices() const
783 {
784  std::vector<G4ThreeVector> vertices(4);
785  vertices[0] = fAnchor;
786  vertices[1] = fP2;
787  vertices[2] = fP3;
788  vertices[3] = fP4;
789 
790  return vertices;
791 }
792 
794 //
795 // GetCubicVolume
796 
798 {
799  return fCubicVolume;
800 }
801 
803 //
804 // GetSurfaceArea
805 
807 {
808  return fSurfaceArea;
809 }
810 
811 // Methods for visualisation
812 
814 //
815 // DescribeYourselfTo
816 
818 {
819  scene.AddSolid (*this);
820 }
821 
823 //
824 // GetExtent
825 
827 {
828  return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
829 }
830 
832 //
833 // CreatePolyhedron
834 
836 {
837  G4Polyhedron *ph=new G4Polyhedron;
838  G4double xyz[4][3];
839  const G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
840  xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
841  xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
842  xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
843  xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
844 
845  ph->createPolyhedron(4,4,xyz,faces);
846 
847  return ph;
848 }
849 
851 //
852 // GetPolyhedron
853 
855 {
856  if (!fpPolyhedron ||
858  fpPolyhedron->GetNumberOfRotationSteps())
859  {
860  delete fpPolyhedron;
862  }
863  return fpPolyhedron;
864 }
865 
866 #endif
G4ThreeVector fP4
Definition: G4Tet.hh:164
G4double fTol
Definition: G4Tet.hh:171
G4double fYMin
Definition: G4Tet.hh:170
G4String GetName() const
G4double fXMax
Definition: G4Tet.hh:170
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tet.cc:572
G4ThreeVector fNormal142
Definition: G4Tet.hh:165
G4double fDz
Definition: G4Tet.hh:171
G4ThreeVector fP2
Definition: G4Tet.hh:164
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetPointOnSurface() const
Definition: G4Tet.cc:761
G4double fYMax
Definition: G4Tet.hh:170
G4bool IsYLimited() const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tet.cc:479
G4GeometryType GetEntityType() const
Definition: G4Tet.cc:693
G4bool IsRotated() const
G4ThreeVector fNormal123
Definition: G4Tet.hh:165
G4double fCdotN142
Definition: G4Tet.hh:169
G4ThreeVector fP3
Definition: G4Tet.hh:164
G4ThreeVector NetTranslation() const
G4bool IsXLimited() const
G4double a
Definition: TRTMaterials.hh:39
G4double fCdotN134
Definition: G4Tet.hh:169
virtual void AddSolid(const G4Box &)=0
int G4int
Definition: G4Types.hh:78
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tet.cc:293
void DumpInfo() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4double fMaxSize
Definition: G4Tet.hh:171
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:782
G4double GetMaxXExtent() const
G4ThreeVector fAnchor
Definition: G4Tet.hh:164
G4double GetMinZExtent() const
G4double fZMin
Definition: G4Tet.hh:170
G4double fXMin
Definition: G4Tet.hh:170
Definition: G4Tet.hh:65
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tet.cc:817
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector GetPointOnFace(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double &area) const
Definition: G4Tet.cc:740
static G4bool CheckDegeneracy(G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4)
Definition: G4Tet.cc:267
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
G4ThreeVector fMiddle
Definition: G4Tet.hh:164
const G4int n
G4double fSurfaceArea
Definition: G4Tet.hh:154
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:835
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double fCdotN234
Definition: G4Tet.hh:169
G4double GetMinXExtent() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetMaxZExtent() const
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tet.cc:711
G4double fDy
Definition: G4Tet.hh:171
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
G4ThreeVector fNormal134
Definition: G4Tet.hh:165
EAxis
Definition: geomdefs.hh:54
G4double GetCubicVolume()
Definition: G4Tet.cc:797
G4double fCdotN123
Definition: G4Tet.hh:169
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double fCubicVolume
Definition: G4Tet.hh:154
G4double fDx
Definition: G4Tet.hh:171
G4double GetSurfaceArea()
Definition: G4Tet.cc:806
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tet.cc:388
G4VisExtent GetExtent() const
Definition: G4Tet.cc:826
#define G4endl
Definition: G4ios.hh:61
G4double GetMaxYExtent() const
static const char CVSVers[]
Definition: G4Tet.hh:160
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double kCarTolerance
Definition: G4VSolid.hh:305
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4VSolid * Clone() const
Definition: G4Tet.cc:702
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tet.cc:662
G4ThreeVector fNormal234
Definition: G4Tet.hh:165
G4Tet & operator=(const G4Tet &rhs)
Definition: G4Tet.cc:235
double G4double
Definition: G4Types.hh:76
virtual ~G4Tet()
Definition: G4Tet.cc:204
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tet.cc:283
G4bool IsZLimited() const
static const double mm
Definition: G4SIunits.hh:102
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:156
G4Polyhedron * GetPolyhedron() const
Definition: G4Tet.cc:854
G4bool warningFlag
Definition: G4Tet.hh:167
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:89
G4double fZMax
Definition: G4Tet.hh:170
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tet.cc:419