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