ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UTet.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4UTet.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 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 // Implementation for G4UTet wrapper class
27 //
28 // 1.11.13 G.Cosmo, CERN
29 // --------------------------------------------------------------------
30 
31 #include "G4Tet.hh"
32 #include "G4UTet.hh"
33 
34 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35 
36 #include "G4AffineTransform.hh"
37 #include "G4VPVParameterisation.hh"
38 #include "G4BoundingEnvelope.hh"
39 
40 using namespace CLHEP;
41 
43 //
44 // Constructor - create a tetrahedron
45 // This class is implemented separately from general polyhedra,
46 // because the simplex geometry can be computed very quickly,
47 // which may become important in situations imported from mesh generators,
48 // in which a very large number of G4Tets are created.
49 // A Tet has all of its geometrical information precomputed
50 //
51 G4UTet::G4UTet(const G4String& pName,
52  G4ThreeVector anchor,
53  G4ThreeVector p2,
54  G4ThreeVector p3,
55  G4ThreeVector p4, G4bool* degeneracyFlag)
56  : Base_t(pName, U3Vector(anchor.x(),anchor.y(),anchor.z()),
57  U3Vector(p2.x(), p2.y(), p2.z()),
58  U3Vector(p3.x(), p3.y(), p3.z()),
59  U3Vector(p4.x(), p4.y(), p4.z()))
60 {
61  G4double fXMin=std::min(std::min(std::min(anchor.x(), p2.x()),p3.x()),p4.x());
62  G4double fXMax=std::max(std::max(std::max(anchor.x(), p2.x()),p3.x()),p4.x());
63  G4double fYMin=std::min(std::min(std::min(anchor.y(), p2.y()),p3.y()),p4.y());
64  G4double fYMax=std::max(std::max(std::max(anchor.y(), p2.y()),p3.y()),p4.y());
65  G4double fZMin=std::min(std::min(std::min(anchor.z(), p2.z()),p3.z()),p4.z());
66  G4double fZMax=std::max(std::max(std::max(anchor.z(), p2.z()),p3.z()),p4.z());
67 
68  G4ThreeVector fMiddle=G4ThreeVector(fXMax+fXMin,fYMax+fYMin,fZMax+fZMin)*0.5;
69  G4double fMaxSize=std::max(std::max(std::max((anchor-fMiddle).mag(),
70  (p2-fMiddle).mag()),
71  (p3-fMiddle).mag()),
72  (p4-fMiddle).mag());
73  // fV<x><y> is vector from vertex <y> to vertex <x>
74  //
75  G4ThreeVector fV21=p2-anchor;
76  G4ThreeVector fV31=p3-anchor;
77  G4ThreeVector fV41=p4-anchor;
78 
79  // make sure this is a correctly oriented set of points for the tetrahedron
80  //
81  G4double signed_vol=fV21.cross(fV31).dot(fV41);
82  G4bool degenerate=std::fabs(signed_vol) < 1e-9*fMaxSize*fMaxSize*fMaxSize;
83 
84  if(degeneracyFlag) *degeneracyFlag = degenerate;
85  else if (degenerate)
86  {
87  G4Exception("G4UTet::G4UTet()", "GeomSolids0002", FatalException,
88  "Degenerate tetrahedron not allowed.");
89  }
90 }
91 
93 //
94 // Fake default constructor - sets only member data and allocates memory
95 // for usage restricted to object persistency.
96 //
97 G4UTet::G4UTet( __void__& a )
98  : Base_t(a)
99 {
100 }
101 
103 //
104 // Destructor
105 //
106 G4UTet::~G4UTet()
107 {
108 }
109 
111 //
112 // Copy constructor
113 //
114 G4UTet::G4UTet(const G4UTet& rhs)
115  : Base_t(rhs)
116 {
117 }
118 
119 
121 //
122 // Assignment operator
123 //
124 G4UTet& G4UTet::operator = (const G4UTet& rhs)
125 {
126  // Check assignment to self
127  //
128  if (this == &rhs) { return *this; }
129 
130  // Copy base class data
131  //
132  Base_t::operator=(rhs);
133 
134  return *this;
135 }
136 
138 //
139 // Dispatch to parameterisation for replication mechanism dimension
140 // computation & modification.
141 //
142 void G4UTet::ComputeDimensions(G4VPVParameterisation*,
143  const G4int,
144  const G4VPhysicalVolume*)
145 {
146 }
147 
149 //
150 // Make a clone of the object
151 //
152 G4VSolid* G4UTet::Clone() const
153 {
154  return new G4UTet(*this);
155 }
156 
158 //
159 // Accessors
160 //
161 std::vector<G4ThreeVector> G4UTet::GetVertices() const
162 {
163  std::vector<U3Vector> vec(4);
164  Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
165  std::vector<G4ThreeVector> vertices;
166  for (unsigned int i=0; i<4; ++i)
167  {
168  G4ThreeVector v(vec[i].x(), vec[i].y(), vec[i].z());
169  vertices.push_back(v);
170  }
171  return vertices;
172 }
173 
175 //
176 // Get bounding box
177 
178 void G4UTet::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
179 {
180  U3Vector vmin, vmax;
181  Base_t::Extent(vmin,vmax);
182  pMin.set(vmin.x(),vmin.y(),vmin.z());
183  pMax.set(vmax.x(),vmax.y(),vmax.z());
184 
185  // Check correctness of the bounding box
186  //
187  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
188  {
189  std::ostringstream message;
190  message << "Bad bounding box (min >= max) for solid: "
191  << GetName() << " !"
192  << "\npMin = " << pMin
193  << "\npMax = " << pMax;
194  G4Exception("G4UTet::BoundingLimits()", "GeomMgt0001",
195  JustWarning, message);
196  StreamInfo(G4cout);
197  }
198 }
199 
201 //
202 // Calculate extent under transform and specified limit
203 
204 G4bool
205 G4UTet::CalculateExtent(const EAxis pAxis,
206  const G4VoxelLimits& pVoxelLimit,
207  const G4AffineTransform& pTransform,
208  G4double& pMin, G4double& pMax) const
209 {
210  G4ThreeVector bmin, bmax;
211 
212  // Check bounding box (bbox)
213  //
214  BoundingLimits(bmin,bmax);
215  G4BoundingEnvelope bbox(bmin,bmax);
216 
217  // Use simple bounding-box to help in the case of complex 3D meshes
218  //
219  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
220 
221 #if 0
222  // Precise extent computation (disabled by default for this shape)
223  //
224  G4bool exist;
225  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
226  {
227  return exist = (pMin < pMax) ? true : false;
228  }
229 
230  // Set bounding envelope (benv) and calculate extent
231  //
232  std::vector<G4ThreeVector> vec = GetVertices();
233 
234  G4ThreeVectorList anchor(1);
235  anchor[0] = vec[0];
236 
238  base[0] = vec[1];
239  base[1] = vec[2];
240  base[2] = vec[3];
241 
242  std::vector<const G4ThreeVectorList *> polygons(2);
243  polygons[0] = &anchor;
244  polygons[1] = &base;
245 
246  G4BoundingEnvelope benv(bmin,bmax,polygons);
247  return exists = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
248 #endif
249 }
250 
252 //
253 // CreatePolyhedron
254 //
255 G4Polyhedron* G4UTet::CreatePolyhedron() const
256 {
257  std::vector<U3Vector> vec(4);
258  Base_t::GetVertices(vec[0], vec[1], vec[2], vec[3]);
259 
260  G4double xyz[4][3];
261  const G4int faces[4][4] = {{1,3,2,0},{1,4,3,0},{1,2,4,0},{2,3,4,0}};
262  for (unsigned int i=0; i<4; ++i)
263  {
264  xyz[i][0] = vec[i].x();
265  xyz[i][1] = vec[i].y();
266  xyz[i][2] = vec[i].z();
267  }
268 
269  G4Polyhedron* ph = new G4Polyhedron;
270  ph->createPolyhedron(4,4,xyz,faces);
271  return ph;
272 }
273 
274 #endif // G4GEOM_USE_USOLIDS