Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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
59const char G4Tet::CVSVers[]="$Id: G4Tet.cc,v 1.16 2010-10-20 08:54:18 gcosmo Exp $";
60
61#include "G4VoxelLimits.hh"
62#include "G4AffineTransform.hh"
63
65
66#include "Randomize.hh"
67
68#include "G4VGraphicsScene.hh"
69#include "G4Polyhedron.hh"
70#include "G4NURBS.hh"
71#include "G4NURBSbox.hh"
72#include "G4VisExtent.hh"
73
74#include "G4ThreeVector.hh"
75
76#include <cmath>
77
78using namespace CLHEP;
79
80////////////////////////////////////////////////////////////////////////
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
90 G4ThreeVector anchor,
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
184//////////////////////////////////////////////////////////////////////////
185//
186// Fake default constructor - sets only member data and allocates memory
187// for usage restricted to object persistency.
188//
189G4Tet::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
200//////////////////////////////////////////////////////////////////////////
201//
202// Destructor
203
205{
206 delete fpPolyhedron;
207}
208
209///////////////////////////////////////////////////////////////////////////////
210//
211// Copy constructor
212
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{
227}
228
229
230///////////////////////////////////////////////////////////////////////////////
231//
232// Assignment operator
233
235{
236 // Check assignment to self
237 //
238 if (this == &rhs) { return *this; }
239
240 // Copy base class data
241 //
243
244 // Copy data
245 //
246 fCubicVolume = rhs.fCubicVolume; fSurfaceArea = rhs.fSurfaceArea;
247 fpPolyhedron = 0; fAnchor = rhs.fAnchor;
248 fP2 = rhs.fP2; fP3 = rhs.fP3; fP4 = rhs.fP4; fMiddle = rhs.fMiddle;
249 fNormal123 = rhs.fNormal123; fNormal142 = rhs.fNormal142;
250 fNormal134 = rhs.fNormal134; fNormal234 = rhs.fNormal234;
251 warningFlag = rhs.warningFlag; fCdotN123 = rhs.fCdotN123;
252 fCdotN142 = rhs.fCdotN142; fCdotN134 = rhs.fCdotN134;
253 fCdotN234 = rhs.fCdotN234; fXMin = rhs.fXMin; fXMax = rhs.fXMax;
254 fYMin = rhs.fYMin; fYMax = rhs.fYMax; fZMin = rhs.fZMin; fZMax = rhs.fZMax;
255 fDx = rhs.fDx; fDy = rhs.fDy; fDz = rhs.fDz; fTol = rhs.fTol;
256 fMaxSize = rhs.fMaxSize;
257
258 return *this;
259}
260
261//////////////////////////////////////////////////////////////////////////
262//
263// CheckDegeneracy
264
266 G4ThreeVector p2,
267 G4ThreeVector p3,
268 G4ThreeVector p4 )
269{
270 G4bool result;
271 G4Tet *object=new G4Tet("temp",anchor,p2,p3,p4,&result);
272 delete object;
273 return result;
274}
275
276//////////////////////////////////////////////////////////////////////////
277//
278// Dispatch to parameterisation for replication mechanism dimension
279// computation & modification.
280
282 const G4int ,
283 const G4VPhysicalVolume* )
284{
285}
286
287//////////////////////////////////////////////////////////////////////////
288//
289// Calculate extent under transform and specified limit
290
292 const G4VoxelLimits& pVoxelLimit,
293 const G4AffineTransform& pTransform,
294 G4double& pMin, G4double& pMax) const
295{
296 G4double xMin,xMax;
297 G4double yMin,yMax;
298 G4double zMin,zMax;
299
300 if (pTransform.IsRotated())
301 {
302 G4ThreeVector pp0=pTransform.TransformPoint(fAnchor);
303 G4ThreeVector pp1=pTransform.TransformPoint(fP2);
304 G4ThreeVector pp2=pTransform.TransformPoint(fP3);
305 G4ThreeVector pp3=pTransform.TransformPoint(fP4);
306
307 xMin = std::min(std::min(std::min(pp0.x(), pp1.x()),pp2.x()),pp3.x());
308 xMax = std::max(std::max(std::max(pp0.x(), pp1.x()),pp2.x()),pp3.x());
309 yMin = std::min(std::min(std::min(pp0.y(), pp1.y()),pp2.y()),pp3.y());
310 yMax = std::max(std::max(std::max(pp0.y(), pp1.y()),pp2.y()),pp3.y());
311 zMin = std::min(std::min(std::min(pp0.z(), pp1.z()),pp2.z()),pp3.z());
312 zMax = std::max(std::max(std::max(pp0.z(), pp1.z()),pp2.z()),pp3.z());
313
314 }
315 else
316 {
317 G4double xoffset = pTransform.NetTranslation().x() ;
318 xMin = xoffset + fXMin;
319 xMax = xoffset + fXMax;
320 G4double yoffset = pTransform.NetTranslation().y() ;
321 yMin = yoffset + fYMin;
322 yMax = yoffset + fYMax;
323 G4double zoffset = pTransform.NetTranslation().z() ;
324 zMin = zoffset + fZMin;
325 zMax = zoffset + fZMax;
326 }
327
328 if (pVoxelLimit.IsXLimited())
329 {
330 if ( (xMin > pVoxelLimit.GetMaxXExtent()+fTol) ||
331 (xMax < pVoxelLimit.GetMinXExtent()-fTol) ) { return false; }
332 else
333 {
334 xMin = std::max(xMin, pVoxelLimit.GetMinXExtent());
335 xMax = std::min(xMax, pVoxelLimit.GetMaxXExtent());
336 }
337 }
338
339 if (pVoxelLimit.IsYLimited())
340 {
341 if ( (yMin > pVoxelLimit.GetMaxYExtent()+fTol) ||
342 (yMax < pVoxelLimit.GetMinYExtent()-fTol) ) { return false; }
343 else
344 {
345 yMin = std::max(yMin, pVoxelLimit.GetMinYExtent());
346 yMax = std::min(yMax, pVoxelLimit.GetMaxYExtent());
347 }
348 }
349
350 if (pVoxelLimit.IsZLimited())
351 {
352 if ( (zMin > pVoxelLimit.GetMaxZExtent()+fTol) ||
353 (zMax < pVoxelLimit.GetMinZExtent()-fTol) ) { return false; }
354 else
355 {
356 zMin = std::max(zMin, pVoxelLimit.GetMinZExtent());
357 zMax = std::min(zMax, pVoxelLimit.GetMaxZExtent());
358 }
359 }
360
361 switch (pAxis)
362 {
363 case kXAxis:
364 pMin=xMin;
365 pMax=xMax;
366 break;
367 case kYAxis:
368 pMin=yMin;
369 pMax=yMax;
370 break;
371 case kZAxis:
372 pMin=zMin;
373 pMax=zMax;
374 break;
375 default:
376 break;
377 }
378
379 return true;
380}
381
382/////////////////////////////////////////////////////////////////////////
383//
384// Return whether point inside/outside/on surface, using tolerance
385
387{
388 G4double r123, r134, r142, r234;
389
390 // this is written to allow if-statement truncation so the outside test
391 // (where most of the world is) can fail very quickly and efficiently
392
393 if ( (r123=p.dot(fNormal123)-fCdotN123) > fTol ||
394 (r134=p.dot(fNormal134)-fCdotN134) > fTol ||
395 (r142=p.dot(fNormal142)-fCdotN142) > fTol ||
396 (r234=p.dot(fNormal234)-fCdotN234) > fTol )
397 {
398 return kOutside; // at least one is out!
399 }
400 else if( (r123 < -fTol)&&(r134 < -fTol)&&(r142 < -fTol)&&(r234 < -fTol) )
401 {
402 return kInside; // all are definitively inside
403 }
404 else
405 {
406 return kSurface; // too close to tell
407 }
408}
409
410///////////////////////////////////////////////////////////////////////
411//
412// Calculate side nearest to p, and return normal
413// If two sides are equidistant, normal of first side (x/y/z)
414// encountered returned.
415// This assumes that we are looking from the inside!
416
418{
419 G4double r123=std::fabs(p.dot(fNormal123)-fCdotN123);
420 G4double r134=std::fabs(p.dot(fNormal134)-fCdotN134);
421 G4double r142=std::fabs(p.dot(fNormal142)-fCdotN142);
422 G4double r234=std::fabs(p.dot(fNormal234)-fCdotN234);
423
424 if( (r123<=r134) && (r123<=r142) && (r123<=r234) ) { return fNormal123; }
425 else if ( (r134<=r142) && (r134<=r234) ) { return fNormal134; }
426 else if (r142 <= r234) { return fNormal142; }
427 return fNormal234;
428}
429
430///////////////////////////////////////////////////////////////////////////
431//
432// Calculate distance to box from an outside point
433// - return kInfinity if no intersection.
434// All this is very unrolled, for speed.
435
437 const G4ThreeVector& v) const
438{
439 G4ThreeVector vu(v.unit()), hp;
440 G4double vdotn, t, tmin=kInfinity;
441
442 G4double extraDistance=10.0*fTol; // a little ways into the solid
443
444 vdotn=-vu.dot(fNormal123);
445 if(vdotn > 1e-12)
446 { // this is a candidate face, since it is pointing at us
447 t=(p.dot(fNormal123)-fCdotN123)/vdotn; // # distance to intersection
448 if( (t>=-fTol) && (t<tmin) )
449 { // if not true, we're going away from this face or it's not close
450 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
451 if ( ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
452 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
453 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
454 {
455 tmin=t;
456 }
457 }
458 }
459
460 vdotn=-vu.dot(fNormal134);
461 if(vdotn > 1e-12)
462 { // # this is a candidate face, since it is pointing at us
463 t=(p.dot(fNormal134)-fCdotN134)/vdotn; // # distance to intersection
464 if( (t>=-fTol) && (t<tmin) )
465 { // if not true, we're going away from this face
466 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
467 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
468 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) &&
469 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
470 {
471 tmin=t;
472 }
473 }
474 }
475
476 vdotn=-vu.dot(fNormal142);
477 if(vdotn > 1e-12)
478 { // # this is a candidate face, since it is pointing at us
479 t=(p.dot(fNormal142)-fCdotN142)/vdotn; // # distance to intersection
480 if( (t>=-fTol) && (t<tmin) )
481 { // if not true, we're going away from this face
482 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
483 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
484 ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
485 ( hp.dot(fNormal234)-fCdotN234 < 0.0 ) )
486 {
487 tmin=t;
488 }
489 }
490 }
491
492 vdotn=-vu.dot(fNormal234);
493 if(vdotn > 1e-12)
494 { // # this is a candidate face, since it is pointing at us
495 t=(p.dot(fNormal234)-fCdotN234)/vdotn; // # distance to intersection
496 if( (t>=-fTol) && (t<tmin) )
497 { // if not true, we're going away from this face
498 hp=p+vu*(t+extraDistance); // a little beyond point of intersection
499 if ( ( hp.dot(fNormal123)-fCdotN123 < 0.0 ) &&
500 ( hp.dot(fNormal134)-fCdotN134 < 0.0 ) &&
501 ( hp.dot(fNormal142)-fCdotN142 < 0.0 ) )
502 {
503 tmin=t;
504 }
505 }
506 }
507
508 return std::max(0.0,tmin);
509}
510
511//////////////////////////////////////////////////////////////////////////
512//
513// Approximate distance to tet.
514// returns distance to sphere centered on bounding box
515// - If inside return 0
516
518{
519 G4double dd=(p-fMiddle).mag() - fMaxSize - fTol;
520 return std::max(0.0, dd);
521}
522
523/////////////////////////////////////////////////////////////////////////
524//
525// Calcluate distance to surface of box from inside
526// by calculating distances to box's x/y/z planes.
527// Smallest distance is exact distance to exiting.
528
530 const G4bool calcNorm,
531 G4bool *validNorm, G4ThreeVector *n) const
532{
533 G4ThreeVector vu(v.unit());
534 G4double t1=kInfinity,t2=kInfinity,t3=kInfinity,t4=kInfinity, vdotn, tt;
535
536 vdotn=vu.dot(fNormal123);
537 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
538 {
539 t1=(fCdotN123-p.dot(fNormal123))/vdotn; // # distance to intersection
540 }
541
542 vdotn=vu.dot(fNormal134);
543 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
544 {
545 t2=(fCdotN134-p.dot(fNormal134))/vdotn; // # distance to intersection
546 }
547
548 vdotn=vu.dot(fNormal142);
549 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
550 {
551 t3=(fCdotN142-p.dot(fNormal142))/vdotn; // # distance to intersection
552 }
553
554 vdotn=vu.dot(fNormal234);
555 if(vdotn > 1e-12) // #we're heading towards this face, so it is a candidate
556 {
557 t4=(fCdotN234-p.dot(fNormal234))/vdotn; // # distance to intersection
558 }
559
560 tt=std::min(std::min(std::min(t1,t2),t3),t4);
561
562 if (warningFlag && (tt == kInfinity || tt < -fTol))
563 {
564 DumpInfo();
565 std::ostringstream message;
566 message << "No good intersection found or already outside!?" << G4endl
567 << "p = " << p / mm << "mm" << G4endl
568 << "v = " << v << G4endl
569 << "t1, t2, t3, t4 (mm) "
570 << t1/mm << ", " << t2/mm << ", " << t3/mm << ", " << t4/mm;
571 G4Exception("G4Tet::DistanceToOut(p,v,...)", "GeomSolids1002",
572 JustWarning, message);
573 if(validNorm)
574 {
575 *validNorm=false; // flag normal as meaningless
576 }
577 }
578 else if(calcNorm && n)
579 {
580 static G4ThreeVector normal;
581 if(tt==t1) { normal=fNormal123; }
582 else if (tt==t2) { normal=fNormal134; }
583 else if (tt==t3) { normal=fNormal142; }
584 else if (tt==t4) { normal=fNormal234; }
585 n=&normal;
586 if(validNorm) { *validNorm=true; }
587 }
588
589 return std::max(tt,0.0); // avoid tt<0.0 by a tiny bit
590 // if we are right on a face
591}
592
593////////////////////////////////////////////////////////////////////////////
594//
595// Calculate exact shortest distance to any boundary from inside
596// - If outside return 0
597
599{
600 G4double t1,t2,t3,t4;
601 t1=fCdotN123-p.dot(fNormal123); // distance to plane, positive if inside
602 t2=fCdotN134-p.dot(fNormal134); // distance to plane
603 t3=fCdotN142-p.dot(fNormal142); // distance to plane
604 t4=fCdotN234-p.dot(fNormal234); // distance to plane
605
606 // if any one of these is negative, we are outside,
607 // so return zero in that case
608
609 G4double tmin=std::min(std::min(std::min(t1,t2),t3),t4);
610 return (tmin < fTol)? 0:tmin;
611}
612
613////////////////////////////////////////////////////////////////////////
614//
615// Create a List containing the transformed vertices
616// Note: Caller has deletion responsibility
617
620{
621 G4ThreeVectorList* vertices = new G4ThreeVectorList();
622
623 if (vertices)
624 {
625 vertices->reserve(4);
626 G4ThreeVector vertex0(fAnchor);
627 G4ThreeVector vertex1(fP2);
628 G4ThreeVector vertex2(fP3);
629 G4ThreeVector vertex3(fP4);
630
631 vertices->push_back(pTransform.TransformPoint(vertex0));
632 vertices->push_back(pTransform.TransformPoint(vertex1));
633 vertices->push_back(pTransform.TransformPoint(vertex2));
634 vertices->push_back(pTransform.TransformPoint(vertex3));
635 }
636 else
637 {
638 DumpInfo();
639 G4Exception("G4Tet::CreateRotatedVertices()",
640 "GeomSolids0003", FatalException,
641 "Error in allocation of vertices. Out of memory !");
642 }
643 return vertices;
644}
645
646//////////////////////////////////////////////////////////////////////////
647//
648// GetEntityType
649
651{
652 return G4String("G4Tet");
653}
654
655//////////////////////////////////////////////////////////////////////////
656//
657// Make a clone of the object
658
660{
661 return new G4Tet(*this);
662}
663
664//////////////////////////////////////////////////////////////////////////
665//
666// Stream object contents to an output stream
667
668std::ostream& G4Tet::StreamInfo(std::ostream& os) const
669{
670 G4int oldprc = os.precision(16);
671 os << "-----------------------------------------------------------\n"
672 << " *** Dump for solid - " << GetName() << " ***\n"
673 << " ===================================================\n"
674 << " Solid type: G4Tet\n"
675 << " Parameters: \n"
676 << " anchor: " << fAnchor/mm << " mm \n"
677 << " p2: " << fP2/mm << " mm \n"
678 << " p3: " << fP3/mm << " mm \n"
679 << " p4: " << fP4/mm << " mm \n"
680 << " normal123: " << fNormal123 << " \n"
681 << " normal134: " << fNormal134 << " \n"
682 << " normal142: " << fNormal142 << " \n"
683 << " normal234: " << fNormal234 << " \n"
684 << "-----------------------------------------------------------\n";
685 os.precision(oldprc);
686
687 return os;
688}
689
690
691////////////////////////////////////////////////////////////////////////
692//
693// GetPointOnFace
694//
695// Auxiliary method for get point on surface
696
697G4ThreeVector G4Tet::GetPointOnFace(G4ThreeVector p1, G4ThreeVector p2,
698 G4ThreeVector p3, G4double& area) const
699{
700 G4double lambda1,lambda2;
701 G4ThreeVector v, w;
702
703 v = p3 - p1;
704 w = p1 - p2;
705
706 lambda1 = RandFlat::shoot(0.,1.);
707 lambda2 = RandFlat::shoot(0.,lambda1);
708
709 area = 0.5*(v.cross(w)).mag();
710
711 return (p2 + lambda1*w + lambda2*v);
712}
713
714////////////////////////////////////////////////////////////////////////////
715//
716// GetPointOnSurface
717
719{
720 G4double chose,aOne,aTwo,aThree,aFour;
721 G4ThreeVector p1, p2, p3, p4;
722
723 p1 = GetPointOnFace(fAnchor,fP2,fP3,aOne);
724 p2 = GetPointOnFace(fAnchor,fP4,fP3,aTwo);
725 p3 = GetPointOnFace(fAnchor,fP4,fP2,aThree);
726 p4 = GetPointOnFace(fP4,fP3,fP2,aFour);
727
728 chose = RandFlat::shoot(0.,aOne+aTwo+aThree+aFour);
729 if( (chose>=0.) && (chose <aOne) ) {return p1;}
730 else if( (chose>=aOne) && (chose < aOne+aTwo) ) {return p2;}
731 else if( (chose>=aOne+aTwo) && (chose<aOne+aTwo+aThree) ) {return p3;}
732 return p4;
733}
734
735////////////////////////////////////////////////////////////////////////
736//
737// GetVertices
738
739std::vector<G4ThreeVector> G4Tet::GetVertices() const
740{
741 std::vector<G4ThreeVector> vertices(4);
742 vertices[0] = fAnchor;
743 vertices[1] = fP2;
744 vertices[2] = fP3;
745 vertices[3] = fP4;
746
747 return vertices;
748}
749
750////////////////////////////////////////////////////////////////////////
751//
752// GetCubicVolume
753
755{
756 return fCubicVolume;
757}
758
759////////////////////////////////////////////////////////////////////////
760//
761// GetSurfaceArea
762
764{
765 return fSurfaceArea;
766}
767
768// Methods for visualisation
769
770////////////////////////////////////////////////////////////////////////
771//
772// DescribeYourselfTo
773
775{
776 scene.AddSolid (*this);
777}
778
779////////////////////////////////////////////////////////////////////////
780//
781// GetExtent
782
784{
785 return G4VisExtent (fXMin, fXMax, fYMin, fYMax, fZMin, fZMax);
786}
787
788////////////////////////////////////////////////////////////////////////
789//
790// CreatePolyhedron
791
793{
795 G4double xyz[4][3];
796 static G4int faces[4][4]={{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
797 xyz[0][0]=fAnchor.x(); xyz[0][1]=fAnchor.y(); xyz[0][2]=fAnchor.z();
798 xyz[1][0]=fP2.x(); xyz[1][1]=fP2.y(); xyz[1][2]=fP2.z();
799 xyz[2][0]=fP3.x(); xyz[2][1]=fP3.y(); xyz[2][2]=fP3.z();
800 xyz[3][0]=fP4.x(); xyz[3][1]=fP4.y(); xyz[3][2]=fP4.z();
801
802 ph->createPolyhedron(4,4,xyz,faces);
803
804 return ph;
805}
806
807////////////////////////////////////////////////////////////////////////
808//
809// CreateNURBS
810
812{
813 return new G4NURBSbox (fDx, fDy, fDz);
814}
815
816////////////////////////////////////////////////////////////////////////
817//
818// GetPolyhedron
819
821{
822 if (!fpPolyhedron ||
824 fpPolyhedron->GetNumberOfRotationSteps())
825 {
826 delete fpPolyhedron;
827 fpPolyhedron = CreatePolyhedron();
828 }
829 return fpPolyhedron;
830}
@ JustWarning
@ FatalException
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:85
#define G4endl
Definition: G4ios.hh:52
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
static double shoot()
Definition: RandFlat.cc:59
G4bool IsRotated() const
G4ThreeVector NetTranslation() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
Definition: G4Tet.hh:57
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Tet.cc:436
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Tet.cc:281
G4Tet(const G4String &pName, G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4, G4bool *degeneracyFlag=0)
Definition: G4Tet.cc:89
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pmin, G4double &pmax) const
Definition: G4Tet.cc:291
G4Tet & operator=(const G4Tet &rhs)
Definition: G4Tet.cc:234
virtual ~G4Tet()
Definition: G4Tet.cc:204
G4Polyhedron * GetPolyhedron() const
Definition: G4Tet.cc:820
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Definition: G4Tet.cc:417
G4GeometryType GetEntityType() const
Definition: G4Tet.cc:650
G4ThreeVector GetPointOnSurface() const
Definition: G4Tet.cc:718
G4double GetSurfaceArea()
Definition: G4Tet.cc:763
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Tet.cc:668
static G4bool CheckDegeneracy(G4ThreeVector anchor, G4ThreeVector p2, G4ThreeVector p3, G4ThreeVector p4)
Definition: G4Tet.cc:265
G4VSolid * Clone() const
Definition: G4Tet.cc:659
G4double GetCubicVolume()
Definition: G4Tet.cc:754
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Definition: G4Tet.cc:529
G4Polyhedron * CreatePolyhedron() const
Definition: G4Tet.cc:792
std::vector< G4ThreeVector > GetVertices() const
Definition: G4Tet.cc:739
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
Definition: G4Tet.cc:619
G4VisExtent GetExtent() const
Definition: G4Tet.cc:783
EInside Inside(const G4ThreeVector &p) const
Definition: G4Tet.cc:386
G4NURBS * CreateNURBS() const
Definition: G4Tet.cc:811
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Definition: G4Tet.cc:774
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
void DumpInfo() const
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4bool IsYLimited() const
G4double GetMinZExtent() const
G4bool IsXLimited() const
G4double GetMaxYExtent() const
G4double GetMaxZExtent() const
G4double GetMinYExtent() const
G4double GetMinXExtent() const
G4bool IsZLimited() const
G4double GetMaxXExtent() const
static G4int GetNumberOfRotationSteps()
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
EAxis
Definition: geomdefs.hh:54
@ kYAxis
Definition: geomdefs.hh:54
@ kXAxis
Definition: geomdefs.hh:54
@ kZAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:58
@ kInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
@ kSurface
Definition: geomdefs.hh:58
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Definition: DoubConv.h:17