Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ReflectedSolid.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 result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id$
28 //
29 //
30 // Implementation for G4ReflectedSolid class for boolean
31 // operations between other solids
32 //
33 // Author: Vladimir Grichine, 23.07.01 (Vladimir.Grichine@cern.ch)
34 //
35 // --------------------------------------------------------------------
36 
37 #include "G4ReflectedSolid.hh"
38 
39 #include <sstream>
40 
41 #include "G4Point3D.hh"
42 #include "G4Normal3D.hh"
43 
44 #include "G4VoxelLimits.hh"
45 
46 #include "G4VPVParameterisation.hh"
47 
48 #include "G4VGraphicsScene.hh"
49 #include "G4Polyhedron.hh"
50 #include "G4NURBS.hh"
51 // #include "G4NURBSbox.hh"
52 
53 
55 //
56 // Constructor using HepTransform3D, in fact HepReflect3D
57 
59  G4VSolid* pSolid ,
60  const G4Transform3D& transform )
61  : G4VSolid(pName), fpPolyhedron(0)
62 {
63  fPtrSolid = pSolid ;
64  G4RotationMatrix rotMatrix ;
65 
67  new G4AffineTransform(rotMatrix, transform.getTranslation()) ;
69  new G4AffineTransform(rotMatrix, transform.getTranslation()) ;
70  fPtrTransform->Invert() ;
71 
72  fDirectTransform3D = new G4Transform3D(transform) ;
73  fPtrTransform3D = new G4Transform3D(transform.inverse()) ;
74 }
75 
77 //
78 
80 {
81  if(fPtrTransform)
82  {
83  delete fPtrTransform; fPtrTransform=0;
85  }
86  if(fPtrTransform3D)
87  {
90  }
91  delete fpPolyhedron;
92 }
93 
95 //
96 
98  : G4VSolid(rhs), fPtrSolid(rhs.fPtrSolid), fpPolyhedron(0)
99 {
104 }
105 
107 //
108 
110 {
111  // Check assignment to self
112  //
113  if (this == &rhs) { return *this; }
114 
115  // Copy base class data
116  //
117  G4VSolid::operator=(rhs);
118 
119  // Copy data
120  //
122  delete fPtrTransform;
124  delete fDirectTransform;
126  delete fPtrTransform3D;
128  delete fDirectTransform3D;
130 
131  return *this;
132 }
133 
135 //
136 
138 {
139  return G4String("G4ReflectedSolid");
140 }
141 
143 {
144  return this;
145 }
146 
148 {
149  return this;
150 }
151 
153 {
154  return fPtrSolid;
155 }
156 
158 
160 {
161  G4AffineTransform aTransform = *fPtrTransform;
162  return aTransform;
163 }
164 
166 {
167  fPtrTransform = &transform ;
168  fpPolyhedron = 0;
169 }
170 
172 
174 {
175  G4AffineTransform aTransform= *fDirectTransform;
176  return aTransform;
177 }
178 
180 {
181  fDirectTransform = &transform ;
182  fpPolyhedron = 0;
183 }
184 
186 
188 {
189  G4Transform3D aTransform = *fPtrTransform3D;
190  return aTransform;
191 }
192 
194 {
195  fPtrTransform3D = &transform ;
196  fpPolyhedron = 0;
197 }
198 
200 
202 {
203  G4Transform3D aTransform= *fDirectTransform3D;
204  return aTransform;
205 }
206 
208 {
209  fDirectTransform3D = &transform ;
210  fpPolyhedron = 0;
211 }
212 
214 
216 {
218  return InvRotation;
219 }
220 
222 {
224 }
225 
227 
229 {
230  return fPtrTransform->NetTranslation();
231 }
232 
234 {
236 }
237 
239 
241 {
243  return Rotation;
244 }
245 
247 {
248  fPtrTransform->SetNetRotation(matrix);
249 }
250 
252 
254 {
256 }
257 
259 {
261 }
262 
264 //
265 //
266 
267 G4bool
269  const G4VoxelLimits& pVoxelLimit,
270  const G4AffineTransform& pTransform,
271  G4double& pMin,
272  G4double& pMax ) const
273 {
274 
275  G4VoxelLimits unLimit;
276  G4AffineTransform unTransform;
277 
278  G4double x1 = -kInfinity, x2 = kInfinity,
279  y1 = -kInfinity, y2 = kInfinity,
280  z1 = -kInfinity, z2 = kInfinity;
281 
282  G4bool existsAfterClip = false ;
283  existsAfterClip =
284  fPtrSolid->CalculateExtent(kXAxis,unLimit,unTransform,x1,x2);
285  existsAfterClip =
286  fPtrSolid->CalculateExtent(kYAxis,unLimit,unTransform,y1,y2);
287  existsAfterClip =
288  fPtrSolid->CalculateExtent(kZAxis,unLimit,unTransform,z1,z2);
289 
290  existsAfterClip = false;
291  pMin = +kInfinity ;
292  pMax = -kInfinity ;
293 
294  G4Transform3D pTransform3D = G4Transform3D(pTransform.NetRotation().inverse(),
295  pTransform.NetTranslation());
296 
297  G4Transform3D transform3D = pTransform3D*(*fDirectTransform3D);
298 
299  G4Point3D tmpPoint;
300 
301  // Calculate rotated vertex coordinates
302 
303  G4ThreeVectorList* vertices = new G4ThreeVectorList();
304 
305  if (vertices)
306  {
307  vertices->reserve(8);
308 
309  G4ThreeVector vertex0(x1,y1,z1) ;
310  tmpPoint = transform3D*G4Point3D(vertex0);
311  vertex0 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
312  vertices->push_back(vertex0);
313 
314  G4ThreeVector vertex1(x2,y1,z1) ;
315  tmpPoint = transform3D*G4Point3D(vertex1);
316  vertex1 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
317  vertices->push_back(vertex1);
318 
319  G4ThreeVector vertex2(x2,y2,z1) ;
320  tmpPoint = transform3D*G4Point3D(vertex2);
321  vertex2 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
322  vertices->push_back(vertex2);
323 
324  G4ThreeVector vertex3(x1,y2,z1) ;
325  tmpPoint = transform3D*G4Point3D(vertex3);
326  vertex3 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
327  vertices->push_back(vertex3);
328 
329  G4ThreeVector vertex4(x1,y1,z2) ;
330  tmpPoint = transform3D*G4Point3D(vertex4);
331  vertex4 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
332  vertices->push_back(vertex4);
333 
334  G4ThreeVector vertex5(x2,y1,z2) ;
335  tmpPoint = transform3D*G4Point3D(vertex5);
336  vertex5 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
337  vertices->push_back(vertex5);
338 
339  G4ThreeVector vertex6(x2,y2,z2) ;
340  tmpPoint = transform3D*G4Point3D(vertex6);
341  vertex6 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
342  vertices->push_back(vertex6);
343 
344  G4ThreeVector vertex7(x1,y2,z2) ;
345  tmpPoint = transform3D*G4Point3D(vertex7);
346  vertex7 = G4ThreeVector(tmpPoint.x(),tmpPoint.y(),tmpPoint.z());
347  vertices->push_back(vertex7);
348  }
349  else
350  {
351  DumpInfo();
352  G4Exception("G4ReflectedSolid::CalculateExtent()",
353  "GeomMgt0003", FatalException,
354  "Error in allocation of vertices. Out of memory !");
355  }
356 
357  ClipCrossSection(vertices,0,pVoxelLimit,pAxis,pMin,pMax) ;
358  ClipCrossSection(vertices,4,pVoxelLimit,pAxis,pMin,pMax) ;
359  ClipBetweenSections(vertices,0,pVoxelLimit,pAxis,pMin,pMax) ;
360 
361  if (pVoxelLimit.IsLimited(pAxis) == false)
362  {
363  if ( pMin != kInfinity || pMax != -kInfinity )
364  {
365  existsAfterClip = true ;
366 
367  // Add 2*tolerance to avoid precision troubles
368 
369  pMin -= kCarTolerance;
370  pMax += kCarTolerance;
371  }
372  }
373  else
374  {
375  G4ThreeVector clipCentre(
376  ( pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
377  ( pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
378  ( pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5);
379 
380  if ( pMin != kInfinity || pMax != -kInfinity )
381  {
382  existsAfterClip = true ;
383 
384 
385  // Check to see if endpoints are in the solid
386 
387  clipCentre(pAxis) = pVoxelLimit.GetMinExtent(pAxis);
388 
389  if (Inside(transform3D.inverse()*G4Point3D(clipCentre)) != kOutside)
390  {
391  pMin = pVoxelLimit.GetMinExtent(pAxis);
392  }
393  else
394  {
395  pMin -= kCarTolerance;
396  }
397  clipCentre(pAxis) = pVoxelLimit.GetMaxExtent(pAxis);
398 
399  if (Inside(transform3D.inverse()*G4Point3D(clipCentre)) != kOutside)
400  {
401  pMax = pVoxelLimit.GetMaxExtent(pAxis);
402  }
403  else
404  {
405  pMax += kCarTolerance;
406  }
407  }
408  // Check for case where completely enveloping clipping volume
409  // If point inside then we are confident that the solid completely
410  // envelopes the clipping volume. Hence set min/max extents according
411  // to clipping volume extents along the specified axis.
412 
413  else if (Inside(transform3D.inverse()*G4Point3D(clipCentre)) != kOutside)
414  {
415  existsAfterClip = true ;
416  pMin = pVoxelLimit.GetMinExtent(pAxis) ;
417  pMax = pVoxelLimit.GetMaxExtent(pAxis) ;
418  }
419  }
420  delete vertices;
421  return existsAfterClip;
422 }
423 
425 //
426 //
427 
429 {
430 
431  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p) ;
432  // G4Point3D newPoint = (*fPtrTransform3D)*G4Point3D(p) ;
433 
434  return fPtrSolid->Inside(G4ThreeVector(newPoint.x(),
435  newPoint.y(),
436  newPoint.z())) ;
437 }
438 
440 //
441 //
442 
445 {
446  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p) ;
447  G4ThreeVector normal =
449  newPoint.y(),
450  newPoint.z() ) ) ;
451  G4Point3D newN = (*fDirectTransform3D)*G4Point3D(normal) ;
452  newN.unit() ;
453 
454  return G4ThreeVector(newN.x(),newN.y(),newN.z()) ;
455 }
456 
458 //
459 // The same algorithm as in DistanceToIn(p)
460 
461 G4double
463  const G4ThreeVector& v ) const
464 {
465  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p) ;
466  G4Point3D newDirection = (*fDirectTransform3D)*G4Point3D(v) ;
467  newDirection.unit() ;
468  return fPtrSolid->DistanceToIn(
469  G4ThreeVector(newPoint.x(),newPoint.y(),newPoint.z()),
470  G4ThreeVector(newDirection.x(),newDirection.y(),newDirection.z())) ;
471 }
472 
474 //
475 // Approximate nearest distance from the point p to the intersection of
476 // two solids
477 
478 G4double
480 {
481  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p) ;
482  return fPtrSolid->DistanceToIn(
483  G4ThreeVector(newPoint.x(),newPoint.y(),newPoint.z())) ;
484 }
485 
487 //
488 // The same algorithm as DistanceToOut(p)
489 
490 G4double
492  const G4ThreeVector& v,
493  const G4bool calcNorm,
494  G4bool *validNorm,
495  G4ThreeVector *n ) const
496 {
497  G4ThreeVector solNorm ;
498 
499  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p) ;
500  G4Point3D newDirection = (*fDirectTransform3D)*G4Point3D(v);
501  newDirection.unit() ;
502 
503  G4double dist =
505  G4ThreeVector(newPoint.x(),newPoint.y(),newPoint.z()),
506  G4ThreeVector(newDirection.x(),newDirection.y(),newDirection.z()),
507  calcNorm, validNorm, &solNorm) ;
508  if(calcNorm)
509  {
510  G4Point3D newN = (*fDirectTransform3D)*G4Point3D(solNorm);
511  newN.unit() ;
512  *n = G4ThreeVector(newN.x(),newN.y(),newN.z());
513  }
514  return dist ;
515 }
516 
518 //
519 // Inverted algorithm of DistanceToIn(p)
520 
521 G4double
523 {
524  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p);
525  return fPtrSolid->DistanceToOut(
526  G4ThreeVector(newPoint.x(),newPoint.y(),newPoint.z()));
527 }
528 
530 //
531 //
532 
533 void
535  const G4int,
536  const G4VPhysicalVolume* )
537 {
538  DumpInfo();
539  G4Exception("G4ReflectedSolid::ComputeDimensions()",
540  "GeomMgt0001", FatalException,
541  "Method not applicable in this context!");
542 }
543 
545 //
546 // Return a point (G4ThreeVector) randomly and uniformly selected
547 // on the solid surface
548 
550 {
552  G4Point3D newPoint = (*fDirectTransform3D)*G4Point3D(p);
553 
554  return G4ThreeVector(newPoint.x(),newPoint.y(),newPoint.z());
555 }
556 
558 //
559 // Make a clone of this object
560 
562 {
563  return new G4ReflectedSolid(*this);
564 }
565 
566 
568 //
569 // Stream object contents to an output stream
570 
571 std::ostream& G4ReflectedSolid::StreamInfo(std::ostream& os) const
572 {
573  os << "-----------------------------------------------------------\n"
574  << " *** Dump for Reflected solid - " << GetName() << " ***\n"
575  << " ===================================================\n"
576  << " Solid type: " << GetEntityType() << "\n"
577  << " Parameters of constituent solid: \n"
578  << "===========================================================\n";
579  fPtrSolid->StreamInfo(os);
580  os << "===========================================================\n"
581  << " Transformations: \n"
582  << " Direct transformation - translation : \n"
583  << " " << fDirectTransform->NetTranslation() << "\n"
584  << " - rotation : \n"
585  << " ";
587  os << "\n"
588  << "===========================================================\n";
589 
590  return os;
591 }
592 
594 //
595 //
596 
597 void
599 {
600  scene.AddSolid (*this);
601 }
602 
604 //
605 //
606 
607 G4Polyhedron*
609 {
610  G4Polyhedron* polyhedron = fPtrSolid->CreatePolyhedron();
611  if (polyhedron)
612  {
613  polyhedron->Transform(*fDirectTransform3D);
614  return polyhedron;
615  }
616  else
617  {
618  std::ostringstream message;
619  message << "Solid - " << GetName()
620  << " - original solid has no" << G4endl
621  << "corresponding polyhedron. Returning NULL!";
622  G4Exception("G4ReflectedSolid::CreatePolyhedron()",
623  "GeomMgt1001", JustWarning, message);
624  return 0;
625  }
626 }
627 
629 //
630 //
631 
632 G4NURBS*
634 {
635  // Take into account local transformation - see CreatePolyhedron.
636  // return fPtrSolid->CreateNURBS() ;
637  return 0;
638 }
639 
641 //
642 //
643 
646 {
647  if (!fpPolyhedron ||
650  {
651  delete fpPolyhedron;
653  }
654  return fpPolyhedron;
655 }