Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4Tet.cc
이 파일의 문서화 페이지로 가기
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 106603 2017-10-16 09:17:44Z 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 // 20160924 - Evgueni Tcherniaev, use G4BoundingEnvelope for CalculateExtent()
55 //
56 // --------------------------------------------------------------------
57 
58 #include "G4Tet.hh"
59 
60 //#if !defined(G4GEOM_USE_UTET)
61 
62 const char G4Tet::CVSVers[]="$Id: G4Tet.cc 106603 2017-10-16 09:17:44Z gcosmo $";
63 
64 #include "G4VoxelLimits.hh"
65 #include "G4AffineTransform.hh"
66 #include "G4BoundingEnvelope.hh"
67 
68 #include "G4VPVParameterisation.hh"
69 
70 #include "Randomize.hh"
71 
72 #include "G4VGraphicsScene.hh"
73 #include "G4Polyhedron.hh"
74 #include "G4VisExtent.hh"
75 
76 #include "G4ThreeVector.hh"
77 
78 #include <cmath>
79 
80 #include "G4AutoLock.hh"
81 
82 namespace
83 {
84  G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
85 }
86 
87 using namespace CLHEP;
88 
90 //
91 // Constructor - create a tetrahedron
92 // This class is implemented separately from general polyhedra,
93 // because the simplex geometry can be computed very quickly,
94 // which may become important in situations imported from mesh generators,
95 // in which a very large number of G4Tets are created.
96 // A Tet has all of its geometrical information precomputed
97 
98 G4Tet::G4Tet(const G4String& pName,
99  G4ThreeVector anchor,
100  G4ThreeVector p2,
101  G4ThreeVector p3,
102  G4ThreeVector p4, G4bool *degeneracyFlag)
103  : G4VSolid(pName), fRebuildPolyhedron(false), fpPolyhedron(0), warningFlag(0)
104 {
105  // fV<x><y> is vector from vertex <y> to vertex <x>
106  //
107  G4ThreeVector fV21=p2-anchor;
108  G4ThreeVector fV31=p3-anchor;
109  G4ThreeVector fV41=p4-anchor;
110 
111  // make sure this is a correctly oriented set of points for the tetrahedron
112  //
113  G4double signed_vol=fV21.cross(fV31).dot(fV41);
114  if(signed_vol<0.0)
115  {
116  G4ThreeVector temp(p4);
117  p4=p3;
118  p3=temp;
119  temp=fV41;
120  fV41=fV31;
121  fV31=temp;
122  }
123  fCubicVolume = std::fabs(signed_vol) / 6.;
124 
125  G4ThreeVector fV24=p2-p4;
126  G4ThreeVector fV43=p4-p3;
127  G4ThreeVector fV32=p3-p2;
128 
129  fXMin=std::min(std::min(std::min(anchor.x(), p2.x()),p3.x()),p4.x());
130  fXMax=std::max(std::max(std::max(anchor.x(), p2.x()),p3.x()),p4.x());
131  fYMin=std::min(std::min(std::min(anchor.y(), p2.y()),p3.y()),p4.y());
132  fYMax=std::max(std::max(std::max(anchor.y(), p2.y()),p3.y()),p4.y());
133  fZMin=std::min(std::min(std::min(anchor.z(), p2.z()),p3.z()),p4.z());
134  fZMax=std::max(std::max(std::max(anchor.z(), p2.z()),p3.z()),p4.z());
135 
136  fDx=(fXMax-fXMin)*0.5; fDy=(fYMax-fYMin)*0.5; fDz=(fZMax-fZMin)*0.5;
137 
138  fMiddle=G4ThreeVector(fXMax+fXMin, fYMax+fYMin, fZMax+fZMin)*0.5;
139  fMaxSize=std::max(std::max(std::max((anchor-fMiddle).mag(),
140  (p2-fMiddle).mag()),
141  (p3-fMiddle).mag()),
142  (p4-fMiddle).mag());
143 
144  G4bool degenerate=std::fabs(signed_vol) < 1e-9*fMaxSize*fMaxSize*fMaxSize;
145 
146  if(degeneracyFlag) *degeneracyFlag=degenerate;
147  else if (degenerate)
148  {
149  G4Exception("G4Tet::G4Tet()", "GeomSolids0002", FatalException,
150  "Degenerate tetrahedron not allowed.");
151  }
152 
153  fTol=1e-9*(std::fabs(fXMin)+std::fabs(fXMax)+std::fabs(fYMin)
154  +std::fabs(fYMax)+std::fabs(fZMin)+std::fabs(fZMax));
155  //fTol=kCarTolerance;
156 
157  fAnchor=anchor;
158  fP2=p2;
159  fP3=p3;
160  fP4=p4;
161 
162  G4ThreeVector fCenter123=(anchor+p2+p3)*(1.0/3.0); // face center
163  G4ThreeVector fCenter134=(anchor+p4+p3)*(1.0/3.0);
164  G4ThreeVector fCenter142=(anchor+p4+p2)*(1.0/3.0);
165  G4ThreeVector fCenter234=(p2+p3+p4)*(1.0/3.0);
166 
167  // compute area of each triangular face by cross product
168  // and sum for total surface area
169 
170  G4ThreeVector normal123=fV31.cross(fV21);
171  G4ThreeVector normal134=fV41.cross(fV31);
172  G4ThreeVector normal142=fV21.cross(fV41);
173  G4ThreeVector normal234=fV32.cross(fV43);
174 
175  fSurfaceArea=(
176  normal123.mag()+
177  normal134.mag()+
178  normal142.mag()+
179  normal234.mag()
180  )/2.0;
181 
182  fNormal123=normal123.unit();
183  fNormal134=normal134.unit();
184  fNormal142=normal142.unit();
185  fNormal234=normal234.unit();
186 
187  fCdotN123=fCenter123.dot(fNormal123);
188  fCdotN134=fCenter134.dot(fNormal134);
189  fCdotN142=fCenter142.dot(fNormal142);
190  fCdotN234=fCenter234.dot(fNormal234);
191 }
192 
194 //
195 // Fake default constructor - sets only member data and allocates memory
196 // for usage restricted to object persistency.
197 //
198 G4Tet::G4Tet( __void__& a )
199  : G4VSolid(a), fCubicVolume(0.), fSurfaceArea(0.),
200  fRebuildPolyhedron(false), fpPolyhedron(0),
201  fAnchor(0,0,0), fP2(0,0,0), fP3(0,0,0), fP4(0,0,0), fMiddle(0,0,0),
202  fNormal123(0,0,0), fNormal142(0,0,0), fNormal134(0,0,0),
203  fNormal234(0,0,0), warningFlag(0),
204  fCdotN123(0.), fCdotN142(0.), fCdotN134(0.), fCdotN234(0.),
205  fXMin(0.), fXMax(0.), fYMin(0.), fYMax(0.), fZMin(0.), fZMax(0.),
206  fDx(0.), fDy(0.), fDz(0.), fTol(0.), fMaxSize(0.)
207 {
208 }
209 
211 //
212 // Destructor
213 
215 {
216  delete fpPolyhedron; fpPolyhedron = 0;
217 }
218 
220 //
221 // Copy constructor
222 
223 G4Tet::G4Tet(const G4Tet& rhs)
224  : G4VSolid(rhs),
225  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
226  fRebuildPolyhedron(false), fpPolyhedron(0), fAnchor(rhs.fAnchor),
227  fP2(rhs.fP2), fP3(rhs.fP3), fP4(rhs.fP4), fMiddle(rhs.fMiddle),
228  fNormal123(rhs.fNormal123), fNormal142(rhs.fNormal142),
229  fNormal134(rhs.fNormal134), fNormal234(rhs.fNormal234),
230  warningFlag(rhs.warningFlag), fCdotN123(rhs.fCdotN123),
231  fCdotN142(rhs.fCdotN142), fCdotN134(rhs.fCdotN134),
232  fCdotN234(rhs.fCdotN234), fXMin(rhs.fXMin), fXMax(rhs.fXMax),
233  fYMin(rhs.fYMin), fYMax(rhs.fYMax), fZMin(rhs.fZMin), fZMax(rhs.fZMax),
234  fDx(rhs.fDx), fDy(rhs.fDy), fDz(rhs.fDz), fTol(rhs.fTol),
235  fMaxSize(rhs.fMaxSize)
236 {
237 }
238 
239 
241 //
242 // Assignment operator
243 
245 {
246  // Check assignment to self
247  //
248  if (this == &rhs) { return *this; }
249 
250  // Copy base class data
251  //
252  G4VSolid::operator=(rhs);
253 
254  // Copy data
255  //
257  fAnchor = rhs.fAnchor;
258  fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
262  fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
263  fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
264  fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
265  fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
266  fMaxSize = rhs.fMaxSize;
267  fRebuildPolyhedron = false;
268  delete fpPolyhedron; fpPolyhedron = 0;
269 
270  return *this;
271 }
272 
274 //
275 // CheckDegeneracy
276 
278  G4ThreeVector p2,
279  G4ThreeVector p3,
280  G4ThreeVector p4 )
281 {
282  G4bool result;
283  G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
284  delete object;
285  return result;
286 }
287 
289 //
290 // Dispatch to parameterisation for replication mechanism dimension
291 // computation & modification.
292 
294  const G4int ,
295  const G4VPhysicalVolume* )
296 {
297 }
298 
300 //
301 // Get bounding box
302 
304 {
305  pMin.set(fXMin,fYMin,fZMin);
306  pMax.set(fXMax,fYMax,fZMax);
307 
308  // Check correctness of the bounding box
309  //
310  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
311  {
312  std::ostringstream message;
313  message << "Bad bounding box (min >= max) for solid: "
314  << GetName() << " !"
315  << "\npMin = " << pMin
316  << "\npMax = " << pMax;
317  G4Exception("G4Tet::BoundingLimits()", "GeomMgt0001", JustWarning, message);
318  DumpInfo();
319  }
320 }
321 
323 //
324 // Calculate extent under transform and specified limit
325 
327  const G4VoxelLimits& pVoxelLimit,
328  const G4AffineTransform& pTransform,
329  G4double& pMin, G4double& pMax) const
330 {
331  G4ThreeVector bmin, bmax;
332  G4bool exist;
333 
334  // Check bounding box (bbox)
335  //
336  BoundingLimits(bmin,bmax);
337  G4BoundingEnvelope bbox(bmin,bmax);
338 #ifdef G4BBOX_EXTENT
339  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
340 #endif
341  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
342  {
343  return exist = (pMin < pMax) ? true : false;
344  }
345 
346  // Set bounding envelope (benv) and calculate extent
347  //
348  std::vector<G4ThreeVector> vec = GetVertices();
349 
350  G4ThreeVectorList anchor(1);
351  anchor[0].set(vec[0].x(),vec[0].y(),vec[0].z());
352 
354  base[0].set(vec[1].x(),vec[1].y(),vec[1].z());
355  base[1].set(vec[2].x(),vec[2].y(),vec[2].z());
356  base[2].set(vec[3].x(),vec[3].y(),vec[3].z());
357 
358  std::vector<const G4ThreeVectorList *> polygons(2);
359  polygons[0] = &anchor;
360  polygons[1] = &base;
361 
362  G4BoundingEnvelope benv(bmin,bmax,polygons);
363  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
364  return exist;
365 }
366 
368 //
369 // Return whether point inside/outside/on surface, using tolerance
370 
372 {
373  G4double r123, r134, r142, r234;
374 
375  // this is written to allow if-statement truncation so the outside test
376  // (where most of the world is) can fail very quickly and efficiently
377 
378  if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
379  (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
380  (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
381  (r234=p.dot(fNormal234)-fCdotN234) > fTol )
382  {
383  return kOutside; // at least one is out!
384  }
385  else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
386  {
387  return kInside; // all are definitively inside
388  }
389  else
390  {
391  return kSurface; // too close to tell
392  }
393 }
394 
396 //
397 // Calculate side nearest to p, and return normal
398 // If two sides are equidistant, normal of first side (x/y/z)
399 // encountered returned.
400 // This assumes that we are looking from the inside!
401 
403 {
404  G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
405  G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
406  G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
407  G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
408 
409  const G4double delta = 0.5*kCarTolerance;
410  G4ThreeVector sumnorm(0., 0., 0.);
411  G4int noSurfaces=0;
412 
413  if (r123 <= delta)
414  {
415  noSurfaces ++;
416  sumnorm= fNormal123;
417  }
418 
419  if (r134 <= delta)
420  {
421  noSurfaces ++;
422  sumnorm += fNormal134;
423  }
424 
425  if (r142 <= delta)
426  {
427  noSurfaces ++;
428  sumnorm += fNormal142;
429  }
430  if (r234 <= delta)
431  {
432  noSurfaces ++;
433  sumnorm += fNormal234;
434  }
435 
436  if( noSurfaces > 0 )
437  {
438  if( noSurfaces == 1 )
439  {
440  return sumnorm;
441  }
442  else
443  {
444  return sumnorm.unit();
445  }
446  }
447  else // Approximative Surface Normal
448  {
449 
450  if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
451  else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; }
452  else if (r142 <= r234) { return fNormal142; }
453  return fNormal234;
454  }
455 }
457 //
458 // Calculate distance to box from an outside point
459 // - return kInfinity if no intersection.
460 // All this is very unrolled, for speed.
461 
463  const G4ThreeVector& v) const
464 {
465  G4ThreeVector vu(v.unit()), hp;
466  G4double vdotn, t, tmin=kInfinity;
467 
468  G4double extraDistance=10.0*fTol; // a little ways into the solid
469 
470  vdotn=-vu.dot(fNormal123);
471  if(vdotn > 1e-12)
472  { // this is a candidate face, since it is pointing at us
473  t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection
474  if( (t>=-fTol) && (t<tmin) )
475  { // if not true, we're going away from this face or it's not close
476  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
477  if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
478  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
479  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
480  {
481  tmin=t;
482  }
483  }
484  }
485 
486  vdotn=-vu.dot(fNormal134);
487  if(vdotn > 1e-12)
488  { // # this is a candidate face, since it is pointing at us
489  t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection
490  if( (t>=-fTol) && (t<tmin) )
491  { // if not true, we're going away from this face
492  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
493  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
494  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
495  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
496  {
497  tmin=t;
498  }
499  }
500  }
501 
502  vdotn=-vu.dot(fNormal142);
503  if(vdotn > 1e-12)
504  { // # this is a candidate face, since it is pointing at us
505  t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection
506  if( (t>=-fTol) && (t<tmin) )
507  { // if not true, we're going away from this face
508  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
509  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
510  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
511  ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
512  {
513  tmin=t;
514  }
515  }
516  }
517 
518  vdotn=-vu.dot(fNormal234);
519  if(vdotn > 1e-12)
520  { // # this is a candidate face, since it is pointing at us
521  t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection
522  if( (t>=-fTol) && (t<tmin) )
523  { // if not true, we're going away from this face
524  hp=p+vu*(t+extraDistance); // a little beyond point of intersection
525  if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
526  ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
527  ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
528  {
529  tmin=t;
530  }
531  }
532  }
533 
534  return std::max(0.0,tmin);
535 }
536 
538 //
539 // Approximate distance to tet.
540 // returns distance to sphere centered on bounding box
541 // - If inside return 0
542 
544 {
545  G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
546  return std::max(0.0, dd);
547 }
548 
550 //
551 // Calcluate distance to surface of box from inside
552 // by calculating distances to box's x/y/z planes.
553 // Smallest distance is exact distance to exiting.
554 
556  const G4bool calcNorm,
557  G4bool *validNorm, G4ThreeVector *n) const
558 {
559  G4ThreeVector vu(v.unit());
561 
562  vdotn=vu.dot(fNormal123);
563  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
564  {
565  t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection
566  }
567 
568  vdotn=vu.dot(fNormal134);
569  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
570  {
571  t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection
572  }
573 
574  vdotn=vu.dot(fNormal142);
575  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
576  {
577  t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection
578  }
579 
580  vdotn=vu.dot(fNormal234);
581  if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
582  {
583  t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection
584  }
585 
586  tt=std::min(std::min(std::min(t1,t2),t3),t4);
587 
588  if (warningFlag && (tt == kInfinity || tt < -fTol))
589  {
590  DumpInfo();
591  std::ostringstream message;
592  message << "No good intersection found or already outside!?" << G4endl
593  << "p = " << p / mm << "mm" << G4endl
594  << "v = " << v << G4endl
595  << "t1, t2, t3, t4 (mm) "
596  << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
597  G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
598  JustWarning, message);
599  if(validNorm)
600  {
601  *validNorm=false; // flag normal as meaningless
602  }
603  }
604  else if(calcNorm && n)
605  {
607  if(tt==t1) { normal=fNormal123; }
608  else if (tt==t2) { normal=fNormal134; }
609  else if (tt==t3) { normal=fNormal142; }
610  else if (tt==t4) { normal=fNormal234; }
611  *n=normal;
612  if(validNorm) { *validNorm=true; }
613  }
614 
615  return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
616  // if we are right on a face
617 }
618 
620 //
621 // Calculate exact shortest distance to any boundary from inside
622 // - If outside return 0
623 
625 {
626  G4double t1,t2,t3,t4;
627  t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside
628  t2=fCdotN134-p.dot(fNormal134); // distance to plane
629  t3=fCdotN142-p.dot(fNormal142); // distance to plane
630  t4=fCdotN234-p.dot(fNormal234); // distance to plane
631 
632  // if any one of these is negative, we are outside,
633  // so return zero in that case
634 
635  G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
636  return (tmin < fTol)? 0:tmin;
637 }
638 
640 //
641 // GetEntityType
642 
644 {
645  return G4String("G4Tet");
646 }
647 
649 //
650 // Make a clone of the object
651 
653 {
654  return new G4Tet(*this);
655 }
656 
658 //
659 // Stream object contents to an output stream
660 
661 std::ostream& G4Tet::StreamInfo(std::ostream& os) const
662 {
663  G4int oldprc = os.precision(16);
664  os << "-----------------------------------------------------------\n"
665  << " *** Dump for solid - " << GetName() << " ***\n"
666  << " ===================================================\n"
667  << " Solid type: G4Tet\n"
668  << " Parameters: \n"
669  << " anchor: " << fAnchor/mm << " mm \n"
670  << " p2: " << fP2/mm << " mm \n"
671  << " p3: " << fP3/mm << " mm \n"
672  << " p4: " << fP4/mm << " mm \n"
673  << " normal123: " << fNormal123 << " \n"
674  << " normal134: " << fNormal134 << " \n"
675  << " normal142: " << fNormal142 << " \n"
676  << " normal234: " << fNormal234 << " \n"
677  << "-----------------------------------------------------------\n";
678  os.precision(oldprc);
679 
680  return os;
681 }
682 
683 
685 //
686 // GetPointOnFace
687 //
688 // Auxiliary method for get point on surface
689 
691  G4ThreeVector p3, G4double& area) const
692 {
693  G4double lambda1,lambda2;
694  G4ThreeVector v, w;
695 
696  v = p3 - p1;
697  w = p1 - p2;
698 
699  lambda1 = G4RandFlat::shoot(0.,1.);
700  lambda2 = G4RandFlat::shoot(0.,lambda1);
701 
702  area = 0.5*(v.cross(w)).mag();
703 
704  return (p2 + lambda1*w + lambda2*v);
705 }
706 
708 //
709 // GetPointOnSurface
710 
712 {
713  G4double chose,aOne,aTwo,aThree,aFour;
714  G4ThreeVector p1, p2, p3, p4;
715 
716  p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
717  p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
718  p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
719  p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
720 
721  chose = G4RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
722  if( (chose>=0.) && (chose <aOne) ) {return p1;}
723  else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
724  else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
725  return p4;
726 }
727 
729 //
730 // GetVertices
731 
732 std::vector<G4ThreeVector> G4Tet::GetVertices() const
733 {
734  std::vector<G4ThreeVector> vertices(4);
735  vertices[0] = fAnchor;
736  vertices[1] = fP2;
737  vertices[2] = fP3;
738  vertices[3] = fP4;
739 
740  return vertices;
741 }
742 
744 //
745 // GetCubicVolume
746 
748 {
749  return fCubicVolume;
750 }
751 
753 //
754 // GetSurfaceArea
755 
757 {
758  return fSurfaceArea;
759 }
760 
761 // Methods for visualisation
762 
764 //
765 // DescribeYourselfTo
766 
768 {
769  scene.AddSolid (*this);
770 }
771 
773 //
774 // GetExtent
775 
777 {
778  return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
779 }
780 
782 //
783 // CreatePolyhedron
784 
786 {
787  G4Polyhedron *ph=new G4Polyhedron;
788  G4double xyz[4][3];
789  const G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
790  xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
791  xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
792  xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
793  xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
794 
795  ph->createPolyhedron(4,4,xyz,faces);
796 
797  return ph;
798 }
799 
801 //
802 // GetPolyhedron
803 
805 {
806  if (!fpPolyhedron ||
810  {
811  G4AutoLock l(&polyhedronMutex);
812  delete fpPolyhedron;
814  fRebuildPolyhedron = false;
815  l.unlock();
816  }
817  return fpPolyhedron;
818 }
819 
820 //#endif
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
Float_t x
Definition: compare.C:6
G4double fCdotN134
Definition: G4Tet.hh:166
void set(double x, double y, double z)
G4ThreeVector fNormal142
Definition: G4Tet.hh:162
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ThreeVector fAnchor
Definition: G4Tet.hh:161
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tet.cc:661
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Definition: G4Tet.cc:303
CLHEP::Hep3Vector G4ThreeVector
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tet.cc:326
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
TTree * t1
Definition: plottest35.C:26
virtual void AddSolid(const G4Box &)=0
static const G4double kInfinity
Definition: geomdefs.hh:42
static const char CVSVers[]
Definition: G4Tet.hh:157
static G4bool CheckDegeneracy(G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4)
Definition: G4Tet.cc:277
G4ThreeVector GetPointOnSurface() const
Definition: G4Tet.cc:711
static constexpr double mm
Definition: G4SIunits.hh:115
G4double fDx
Definition: G4Tet.hh:168
G4double fZMax
Definition: G4Tet.hh:167
#define G4endl
Definition: G4ios.hh:61
Float_t y
Definition: compare.C:6
const char * p
Definition: xmltok.h:285
G4double fYMax
Definition: G4Tet.hh:167
Double_t z
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tet.cc:462
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
G4ThreeVector fNormal123
Definition: G4Tet.hh:162
virtual ~G4Tet()
Definition: G4Tet.cc:214
double z() const
G4double fCdotN123
Definition: G4Tet.hh:166
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
double dot(const Hep3Vector &) const
G4double fCubicVolume
Definition: G4Tet.hh:150
G4double GetSurfaceArea()
Definition: G4Tet.cc:756
G4double fXMin
Definition: G4Tet.hh:167
const XML_Char int const XML_Char int const XML_Char * base
Definition: expat.h:331
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tet.cc:371
G4double GetCubicVolume()
Definition: G4Tet.cc:747
G4double fYMin
Definition: G4Tet.hh:167
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:732
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double fDz
Definition: G4Tet.hh:168
G4bool BoundingBoxVsVoxelLimits(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
G4double fCdotN142
Definition: G4Tet.hh:166
G4ThreeVector GetPointOnFace(G4ThreeVector p1, G4ThreeVector p2, G4ThreeVector p3, G4double &area) const
Definition: G4Tet.cc:690
G4Polyhedron * fpPolyhedron
Definition: G4Tet.hh:153
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:88
G4Tet & operator=(const G4Tet &rhs)
Definition: G4Tet.cc:244
G4ThreeVector fMiddle
Definition: G4Tet.hh:161
G4double kCarTolerance
Definition: G4VSolid.hh:307
G4bool warningFlag
Definition: G4Tet.hh:164
G4ThreeVector fP3
Definition: G4Tet.hh:161
G4ThreeVector fNormal234
Definition: G4Tet.hh:162
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tet.cc:767
G4double fDy
Definition: G4Tet.hh:168
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double fZMin
Definition: G4Tet.hh:167
G4double fCdotN234
Definition: G4Tet.hh:166
Hep3Vector unit() const
Hep3Vector cross(const Hep3Vector &) const
static G4int GetNumberOfRotationSteps()
G4double G4ParticleHPJENDLHEData::G4double result
TTree * t2
Definition: plottest35.C:36
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:785
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
double mag() const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
int G4int
Definition: G4Types.hh:78
EInside
Definition: geomdefs.hh:58
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:98
EAxis
Definition: geomdefs.hh:54
G4String GetName() const
void DumpInfo() const
double x() const
G4ThreeVector fP2
Definition: G4Tet.hh:161
G4VSolid * Clone() const
Definition: G4Tet.cc:652
Char_t n[5]
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tet.cc:555
G4GeometryType GetEntityType() const
Definition: G4Tet.cc:643
double y() const
G4ThreeVector fP4
Definition: G4Tet.hh:161
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tet.cc:402
G4double fXMax
Definition: G4Tet.hh:167
G4VisExtent GetExtent() const
Definition: G4Tet.cc:776
G4ThreeVector fNormal134
Definition: G4Tet.hh:162
Definition: G4Tet.hh:66
G4Polyhedron * GetPolyhedron() const
Definition: G4Tet.cc:804
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tet.cc:293
G4double fMaxSize
Definition: G4Tet.hh:168
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double fSurfaceArea
Definition: G4Tet.hh:150
G4double fTol
Definition: G4Tet.hh:168
std::mutex G4Mutex
Definition: G4Threading.hh:84
G4bool fRebuildPolyhedron
Definition: G4Tet.hh:152