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