ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UCutTubs.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4UCutTubs.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 G4UCutTubs wrapper class
27 //
28 // 07.07.17 G.Cosmo, CERN/PH
29 // --------------------------------------------------------------------
30 
31 #include "G4CutTubs.hh"
32 #include "G4UCutTubs.hh"
33 
34 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35 
36 #include "G4GeomTools.hh"
37 #include "G4AffineTransform.hh"
38 #include "G4VPVParameterisation.hh"
39 #include "G4BoundingEnvelope.hh"
40 
41 using namespace CLHEP;
42 
44 //
45 // Constructor - check parameters, convert angles so 0<sphi+dpshi<=2_PI
46 // - note if pdphi>2PI then reset to 2PI
47 
48 G4UCutTubs::G4UCutTubs( const G4String& pName,
49  G4double pRMin, G4double pRMax,
50  G4double pDz,
51  G4double pSPhi, G4double pDPhi,
52  G4ThreeVector pLowNorm,
53  G4ThreeVector pHighNorm )
54  : Base_t(pName, pRMin, pRMax, pDz, pSPhi, pDPhi,
55  pLowNorm.x(), pLowNorm.y(), pLowNorm.z(),
56  pHighNorm.x(), pHighNorm.y(), pHighNorm.z())
57 {
58 }
59 
61 //
62 // Fake default constructor - sets only member data and allocates memory
63 // for usage restricted to object persistency.
64 //
65 G4UCutTubs::G4UCutTubs( __void__& a )
66  : Base_t(a)
67 {
68 }
69 
71 //
72 // Destructor
73 
74 G4UCutTubs::~G4UCutTubs()
75 {
76 }
77 
79 //
80 // Copy constructor
81 
82 G4UCutTubs::G4UCutTubs(const G4UCutTubs& rhs)
83  : Base_t(rhs)
84 {
85 }
86 
88 //
89 // Assignment operator
90 
91 G4UCutTubs& G4UCutTubs::operator = (const G4UCutTubs& rhs)
92 {
93  // Check assignment to self
94  //
95  if (this == &rhs) { return *this; }
96 
97  // Copy base class data
98  //
99  Base_t::operator=(rhs);
100 
101  return *this;
102 }
103 
105 //
106 // Accessors and modifiers
107 
108 G4double G4UCutTubs::GetInnerRadius() const
109 {
110  return rmin();
111 }
112 G4double G4UCutTubs::GetOuterRadius() const
113 {
114  return rmax();
115 }
116 G4double G4UCutTubs::GetZHalfLength() const
117 {
118  return z();
119 }
120 G4double G4UCutTubs::GetStartPhiAngle() const
121 {
122  return sphi();
123 }
124 G4double G4UCutTubs::GetDeltaPhiAngle() const
125 {
126  return dphi();
127 }
128 G4double G4UCutTubs::GetSinStartPhi() const
129 {
130  return std::sin(GetStartPhiAngle());
131 }
132 G4double G4UCutTubs::GetCosStartPhi() const
133 {
134  return std::cos(GetStartPhiAngle());
135 }
136 G4double G4UCutTubs::GetSinEndPhi() const
137 {
138  return std::sin(GetStartPhiAngle()+GetDeltaPhiAngle());
139 }
140 G4double G4UCutTubs::GetCosEndPhi() const
141 {
142  return std::cos(GetStartPhiAngle()+GetDeltaPhiAngle());
143 }
144 G4ThreeVector G4UCutTubs::GetLowNorm () const
145 {
146  U3Vector lc = BottomNormal();
147  return G4ThreeVector(lc.x(), lc.y(), lc.z());
148 }
149 G4ThreeVector G4UCutTubs::GetHighNorm () const
150 {
151  U3Vector hc = TopNormal();
152  return G4ThreeVector(hc.x(), hc.y(), hc.z());
153 }
154 
155 void G4UCutTubs::SetInnerRadius(G4double newRMin)
156 {
157  SetRMin(newRMin);
158  fRebuildPolyhedron = true;
159 }
160 void G4UCutTubs::SetOuterRadius(G4double newRMax)
161 {
162  SetRMax(newRMax);
163  fRebuildPolyhedron = true;
164 }
165 void G4UCutTubs::SetZHalfLength(G4double newDz)
166 {
167  SetDz(newDz);
168  fRebuildPolyhedron = true;
169 }
170 void G4UCutTubs::SetStartPhiAngle(G4double newSPhi, G4bool)
171 {
172  SetSPhi(newSPhi);
173  fRebuildPolyhedron = true;
174 }
175 void G4UCutTubs::SetDeltaPhiAngle(G4double newDPhi)
176 {
177  SetDPhi(newDPhi);
178  fRebuildPolyhedron = true;
179 }
180 
182 //
183 // Make a clone of the object
184 
185 G4VSolid* G4UCutTubs::Clone() const
186 {
187  return new G4UCutTubs(*this);
188 }
189 
191 //
192 // Get bounding box
193 
194 void G4UCutTubs::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
195 {
196  static G4bool checkBBox = true;
197 
198  G4double rmin = GetInnerRadius();
199  G4double rmax = GetOuterRadius();
200  G4double dz = GetZHalfLength();
201  G4double dphi = GetDeltaPhiAngle();
202 
203  G4double sinSphi = GetSinStartPhi();
204  G4double cosSphi = GetCosStartPhi();
205  G4double sinEphi = GetSinEndPhi();
206  G4double cosEphi = GetCosEndPhi();
207 
209  G4double mag, topx, topy, dists, diste;
210  G4bool iftop;
211 
212  // Find Zmin
213  //
214  G4double zmin;
215  norm = GetLowNorm();
216  mag = std::sqrt(norm.x()*norm.x() + norm.y()*norm.y());
217  topx = (mag == 0) ? 0 : -rmax*norm.x()/mag;
218  topy = (mag == 0) ? 0 : -rmax*norm.y()/mag;
219  dists = sinSphi*topx - cosSphi*topy;
220  diste = -sinEphi*topx + cosEphi*topy;
221  if (dphi > pi)
222  {
223  iftop = true;
224  if (dists > 0 && diste > 0)iftop = false;
225  }
226  else
227  {
228  iftop = false;
229  if (dists <= 0 && diste <= 0) iftop = true;
230  }
231  if (iftop)
232  {
233  zmin = -(norm.x()*topx + norm.y()*topy)/norm.z() - dz;
234  }
235  else
236  {
237  G4double z1 = -rmin*(norm.x()*cosSphi + norm.y()*sinSphi)/norm.z() - dz;
238  G4double z2 = -rmin*(norm.x()*cosEphi + norm.y()*sinEphi)/norm.z() - dz;
239  G4double z3 = -rmax*(norm.x()*cosSphi + norm.y()*sinSphi)/norm.z() - dz;
240  G4double z4 = -rmax*(norm.x()*cosEphi + norm.y()*sinEphi)/norm.z() - dz;
241  zmin = std::min(std::min(std::min(z1,z2),z3),z4);
242  }
243 
244  // Find Zmax
245  //
246  G4double zmax;
247  norm = GetHighNorm();
248  mag = std::sqrt(norm.x()*norm.x() + norm.y()*norm.y());
249  topx = (mag == 0) ? 0 : -rmax*norm.x()/mag;
250  topy = (mag == 0) ? 0 : -rmax*norm.y()/mag;
251  dists = sinSphi*topx - cosSphi*topy;
252  diste = -sinEphi*topx + cosEphi*topy;
253  if (dphi > pi)
254  {
255  iftop = true;
256  if (dists > 0 && diste > 0) iftop = false;
257  }
258  else
259  {
260  iftop = false;
261  if (dists <= 0 && diste <= 0) iftop = true;
262  }
263  if (iftop)
264  {
265  zmax = -(norm.x()*topx + norm.y()*topy)/norm.z() + dz;
266  }
267  else
268  {
269  G4double z1 = -rmin*(norm.x()*cosSphi + norm.y()*sinSphi)/norm.z() + dz;
270  G4double z2 = -rmin*(norm.x()*cosEphi + norm.y()*sinEphi)/norm.z() + dz;
271  G4double z3 = -rmax*(norm.x()*cosSphi + norm.y()*sinSphi)/norm.z() + dz;
272  G4double z4 = -rmax*(norm.x()*cosEphi + norm.y()*sinEphi)/norm.z() + dz;
273  zmax = std::max(std::max(std::max(z1,z2),z3),z4);
274  }
275 
276  // Find bounding box
277  //
278  if (GetDeltaPhiAngle() < twopi)
279  {
280  G4TwoVector vmin,vmax;
281  G4GeomTools::DiskExtent(rmin,rmax,
282  GetSinStartPhi(),GetCosStartPhi(),
283  GetSinEndPhi(),GetCosEndPhi(),
284  vmin,vmax);
285  pMin.set(vmin.x(),vmin.y(), zmin);
286  pMax.set(vmax.x(),vmax.y(), zmax);
287  }
288  else
289  {
290  pMin.set(-rmax,-rmax, zmin);
291  pMax.set( rmax, rmax, zmax);
292  }
293 
294  // Check correctness of the bounding box
295  //
296  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
297  {
298  std::ostringstream message;
299  message << "Bad bounding box (min >= max) for solid: "
300  << GetName() << " !"
301  << "\npMin = " << pMin
302  << "\npMax = " << pMax;
303  G4Exception("G4CUutTubs::BoundingLimits()", "GeomMgt0001",
304  JustWarning, message);
305  StreamInfo(G4cout);
306  }
307 
308  // Check consistency of bounding boxes
309  //
310  if (checkBBox)
311  {
312  U3Vector vmin, vmax;
313  Extent(vmin,vmax);
314  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
315  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
316  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
317  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
318  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
319  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
320  {
321  std::ostringstream message;
322  message << "Inconsistency in bounding boxes for solid: "
323  << GetName() << " !"
324  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
325  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
326  G4Exception("G4UCutTubs::BoundingLimits()", "GeomMgt0001",
327  JustWarning, message);
328  checkBBox = false;
329  }
330  }
331 }
332 
334 //
335 // Calculate extent under transform and specified limit
336 
337 G4bool
338 G4UCutTubs::CalculateExtent(const EAxis pAxis,
339  const G4VoxelLimits& pVoxelLimit,
340  const G4AffineTransform& pTransform,
341  G4double& pMin, G4double& pMax) const
342 {
343  G4ThreeVector bmin, bmax;
344  G4bool exist;
345 
346  // Get bounding box
347  BoundingLimits(bmin,bmax);
348 
349  // Check bounding box
350  G4BoundingEnvelope bbox(bmin,bmax);
351 #ifdef G4BBOX_EXTENT
352  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
353 #endif
354  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
355  {
356  return exist = (pMin < pMax) ? true : false;
357  }
358 
359  // Get parameters of the solid
360  G4double rmin = GetInnerRadius();
361  G4double rmax = GetOuterRadius();
362  G4double dphi = GetDeltaPhiAngle();
363  G4double zmin = bmin.z();
364  G4double zmax = bmax.z();
365 
366  // Find bounding envelope and calculate extent
367  //
368  const G4int NSTEPS = 24; // number of steps for whole circle
369  G4double astep = twopi/NSTEPS; // max angle for one step
370  G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
371  G4double ang = dphi/ksteps;
372 
373  G4double sinHalf = std::sin(0.5*ang);
374  G4double cosHalf = std::cos(0.5*ang);
375  G4double sinStep = 2.*sinHalf*cosHalf;
376  G4double cosStep = 1. - 2.*sinHalf*sinHalf;
377  G4double rext = rmax/cosHalf;
378 
379  // bounding envelope for full cylinder consists of two polygons,
380  // in other cases it is a sequence of quadrilaterals
381  if (rmin == 0 && dphi == twopi)
382  {
383  G4double sinCur = sinHalf;
384  G4double cosCur = cosHalf;
385 
386  G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
387  for (G4int k=0; k<NSTEPS; ++k)
388  {
389  baseA[k].set(rext*cosCur,rext*sinCur,zmin);
390  baseB[k].set(rext*cosCur,rext*sinCur,zmax);
391 
392  G4double sinTmp = sinCur;
393  sinCur = sinCur*cosStep + cosCur*sinStep;
394  cosCur = cosCur*cosStep - sinTmp*sinStep;
395  }
396  std::vector<const G4ThreeVectorList *> polygons(2);
397  polygons[0] = &baseA;
398  polygons[1] = &baseB;
399  G4BoundingEnvelope benv(bmin,bmax,polygons);
400  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
401  }
402  else
403  {
404  G4double sinStart = GetSinStartPhi();
405  G4double cosStart = GetCosStartPhi();
406  G4double sinEnd = GetSinEndPhi();
407  G4double cosEnd = GetCosEndPhi();
408  G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
409  G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
410 
411  // set quadrilaterals
412  G4ThreeVectorList pols[NSTEPS+2];
413  for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(4);
414  pols[0][0].set(rmin*cosStart,rmin*sinStart,zmax);
415  pols[0][1].set(rmin*cosStart,rmin*sinStart,zmin);
416  pols[0][2].set(rmax*cosStart,rmax*sinStart,zmin);
417  pols[0][3].set(rmax*cosStart,rmax*sinStart,zmax);
418  for (G4int k=1; k<ksteps+1; ++k)
419  {
420  pols[k][0].set(rmin*cosCur,rmin*sinCur,zmax);
421  pols[k][1].set(rmin*cosCur,rmin*sinCur,zmin);
422  pols[k][2].set(rext*cosCur,rext*sinCur,zmin);
423  pols[k][3].set(rext*cosCur,rext*sinCur,zmax);
424 
425  G4double sinTmp = sinCur;
426  sinCur = sinCur*cosStep + cosCur*sinStep;
427  cosCur = cosCur*cosStep - sinTmp*sinStep;
428  }
429  pols[ksteps+1][0].set(rmin*cosEnd,rmin*sinEnd,zmax);
430  pols[ksteps+1][1].set(rmin*cosEnd,rmin*sinEnd,zmin);
431  pols[ksteps+1][2].set(rmax*cosEnd,rmax*sinEnd,zmin);
432  pols[ksteps+1][3].set(rmax*cosEnd,rmax*sinEnd,zmax);
433 
434  // set envelope and calculate extent
435  std::vector<const G4ThreeVectorList *> polygons;
436  polygons.resize(ksteps+2);
437  for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
438  G4BoundingEnvelope benv(bmin,bmax,polygons);
439  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
440  }
441  return exist;
442 }
443 
445 //
446 // Return real Z coordinate of point on Cutted +/- fDZ plane
447 
448 G4double G4UCutTubs::GetCutZ(const G4ThreeVector& p) const
449 {
450  G4double newz = p.z(); // p.z() should be either +fDz or -fDz
451  G4ThreeVector fLowNorm = GetLowNorm();
452  G4ThreeVector fHighNorm = GetHighNorm();
453 
454  if (p.z()<0)
455  {
456  if(fLowNorm.z()!=0.)
457  {
458  newz = -GetZHalfLength()
459  - (p.x()*fLowNorm.x()+p.y()*fLowNorm.y())/fLowNorm.z();
460  }
461  }
462  else
463  {
464  if(fHighNorm.z()!=0.)
465  {
466  newz = GetZHalfLength()
467  - (p.x()*fHighNorm.x()+p.y()*fHighNorm.y())/fHighNorm.z();
468  }
469  }
470  return newz;
471 }
472 
474 //
475 // Create polyhedron for visualization
476 //
477 G4Polyhedron* G4UCutTubs::CreatePolyhedron() const
478 {
479  typedef G4double G4double3[3];
480  typedef G4int G4int4[4];
481 
482  G4Polyhedron *ph = new G4Polyhedron;
483  G4Polyhedron *ph1 = new G4PolyhedronTubs(GetInnerRadius(),
484  GetOuterRadius(),
485  GetZHalfLength(),
486  GetStartPhiAngle(),
487  GetDeltaPhiAngle());
488  G4int nn=ph1->GetNoVertices();
489  G4int nf=ph1->GetNoFacets();
490  G4double3* xyz = new G4double3[nn]; // number of nodes
491  G4int4* faces = new G4int4[nf] ; // number of faces
492  G4double fDz = GetZHalfLength();
493 
494  for(G4int i=0; i<nn; ++i)
495  {
496  xyz[i][0]=ph1->GetVertex(i+1).x();
497  xyz[i][1]=ph1->GetVertex(i+1).y();
498  G4double tmpZ=ph1->GetVertex(i+1).z();
499  if (tmpZ>=fDz-kCarTolerance)
500  {
501  xyz[i][2]=GetCutZ(G4ThreeVector(xyz[i][0],xyz[i][1],fDz));
502  }
503  else if(tmpZ<=-fDz+kCarTolerance)
504  {
505  xyz[i][2]=GetCutZ(G4ThreeVector(xyz[i][0],xyz[i][1],-fDz));
506  }
507  else
508  {
509  xyz[i][2]=tmpZ;
510  }
511  }
512  G4int iNodes[4];
513  G4int *iEdge=0;
514  G4int n;
515  for(G4int i=0; i<nf; ++i)
516  {
517  ph1->GetFacet(i+1,n,iNodes,iEdge);
518  for(G4int k=0; k<n; ++k)
519  {
520  faces[i][k]=iNodes[k];
521  }
522  for(G4int k=n; k<4; ++k)
523  {
524  faces[i][k]=0;
525  }
526  }
527  ph->createPolyhedron(nn,nf,xyz,faces);
528 
529  delete [] xyz;
530  delete [] faces;
531  delete ph1;
532 
533  return ph;
534 }
535 
536 #endif // G4GEOM_USE_USOLIDS