Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UTubs.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 G4UTubs wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4Tubs.hh"
34 #include "G4UTubs.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 G4UTubs::G4UTubs( const G4String& pName,
51  G4double pRMin, G4double pRMax,
52  G4double pDz,
53  G4double pSPhi, G4double pDPhi )
54  : Base_t(pName, pRMin, pRMax, pDz, pSPhi, pDPhi)
55 {
56 }
57 
59 //
60 // Fake default constructor - sets only member data and allocates memory
61 // for usage restricted to object persistency.
62 //
63 G4UTubs::G4UTubs( __void__& a )
64  : Base_t(a)
65 {
66 }
67 
69 //
70 // Destructor
71 
72 G4UTubs::~G4UTubs()
73 {
74 }
75 
77 //
78 // Copy constructor
79 
80 G4UTubs::G4UTubs(const G4UTubs& rhs)
81  : Base_t(rhs)
82 {
83 }
84 
86 //
87 // Assignment operator
88 
89 G4UTubs& G4UTubs::operator = (const G4UTubs& rhs)
90 {
91  // Check assignment to self
92  //
93  if (this == &rhs) { return *this; }
94 
95  // Copy base class data
96  //
97  Base_t::operator=(rhs);
98 
99  return *this;
100 }
101 
103 //
104 // Accessors and modifiers
105 
106 G4double G4UTubs::GetInnerRadius() const
107 {
108  return rmin();
109 }
110 G4double G4UTubs::GetOuterRadius() const
111 {
112  return rmax();
113 }
114 G4double G4UTubs::GetZHalfLength() const
115 {
116  return z();
117 }
118 G4double G4UTubs::GetStartPhiAngle() const
119 {
120  return sphi();
121 }
122 G4double G4UTubs::GetDeltaPhiAngle() const
123 {
124  return dphi();
125 }
126 G4double G4UTubs::GetSinStartPhi() const
127 {
128  return std::sin(GetStartPhiAngle());
129 }
130 G4double G4UTubs::GetCosStartPhi() const
131 {
132  return std::cos(GetStartPhiAngle());
133 }
134 G4double G4UTubs::GetSinEndPhi() const
135 {
136  return std::sin(GetStartPhiAngle()+GetDeltaPhiAngle());
137 }
138 G4double G4UTubs::GetCosEndPhi() const
139 {
140  return std::cos(GetStartPhiAngle()+GetDeltaPhiAngle());
141 }
142 
143 void G4UTubs::SetInnerRadius(G4double newRMin)
144 {
145  SetRMin(newRMin);
146  fRebuildPolyhedron = true;
147 }
148 void G4UTubs::SetOuterRadius(G4double newRMax)
149 {
150  SetRMax(newRMax);
151  fRebuildPolyhedron = true;
152 }
153 void G4UTubs::SetZHalfLength(G4double newDz)
154 {
155  SetDz(newDz);
156  fRebuildPolyhedron = true;
157 }
158 void G4UTubs::SetStartPhiAngle(G4double newSPhi, G4bool)
159 {
160  SetSPhi(newSPhi);
161  fRebuildPolyhedron = true;
162 }
163 void G4UTubs::SetDeltaPhiAngle(G4double newDPhi)
164 {
165  SetDPhi(newDPhi);
166  fRebuildPolyhedron = true;
167 }
168 
170 //
171 // Dispatch to parameterisation for replication mechanism dimension
172 // computation & modification.
173 
174 void G4UTubs::ComputeDimensions( G4VPVParameterisation* p,
175  const G4int n,
176  const G4VPhysicalVolume* pRep )
177 {
178  p->ComputeDimensions(*(G4Tubs*)this,n,pRep) ;
179 }
180 
182 //
183 // Make a clone of the object
184 
185 G4VSolid* G4UTubs::Clone() const
186 {
187  return new G4UTubs(*this);
188 }
189 
191 //
192 // Get bounding box
193 
194 void G4UTubs::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 
202  // Find bounding box
203  //
204  if (GetDeltaPhiAngle() < twopi)
205  {
206  G4TwoVector vmin,vmax;
207  G4GeomTools::DiskExtent(rmin,rmax,
208  GetSinStartPhi(),GetCosStartPhi(),
209  GetSinEndPhi(),GetCosEndPhi(),
210  vmin,vmax);
211  pMin.set(vmin.x(),vmin.y(),-dz);
212  pMax.set(vmax.x(),vmax.y(), dz);
213  }
214  else
215  {
216  pMin.set(-rmax,-rmax,-dz);
217  pMax.set( rmax, rmax, dz);
218  }
219 
220  // Check correctness of the bounding box
221  //
222  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
223  {
224  std::ostringstream message;
225  message << "Bad bounding box (min >= max) for solid: "
226  << GetName() << " !"
227  << "\npMin = " << pMin
228  << "\npMax = " << pMax;
229  G4Exception("G4UTubs::BoundingLimits()", "GeomMgt0001",
230  JustWarning, message);
231  StreamInfo(G4cout);
232  }
233 
234  // Check consistency of bounding boxes
235  //
236  if (checkBBox)
237  {
238  U3Vector vmin, vmax;
239  Extent(vmin,vmax);
240  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
241  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
242  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
243  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
244  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
245  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
246  {
247  std::ostringstream message;
248  message << "Inconsistency in bounding boxes for solid: "
249  << GetName() << " !"
250  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
251  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
252  G4Exception("G4UTubs::BoundingLimits()", "GeomMgt0001",
253  JustWarning, message);
254  checkBBox = false;
255  }
256  }
257 }
258 
260 //
261 // Calculate extent under transform and specified limit
262 
263 G4bool
264 G4UTubs::CalculateExtent(const EAxis pAxis,
265  const G4VoxelLimits& pVoxelLimit,
266  const G4AffineTransform& pTransform,
267  G4double& pMin, G4double& pMax) const
268 {
269  G4ThreeVector bmin, bmax;
270  G4bool exist;
271 
272  // Get bounding box
273  BoundingLimits(bmin,bmax);
274 
275  // Check bounding box
276  G4BoundingEnvelope bbox(bmin,bmax);
277 #ifdef G4BBOX_EXTENT
278  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
279 #endif
280  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
281  {
282  return exist = (pMin < pMax) ? true : false;
283  }
284 
285  // Get parameters of the solid
286  G4double rmin = GetInnerRadius();
287  G4double rmax = GetOuterRadius();
288  G4double dz = GetZHalfLength();
289  G4double dphi = GetDeltaPhiAngle();
290 
291  // Find bounding envelope and calculate extent
292  //
293  const G4int NSTEPS = 24; // number of steps for whole circle
294  G4double astep = twopi/NSTEPS; // max angle for one step
295  G4int ksteps = (dphi <= astep) ? 1 : (G4int)((dphi-deg)/astep) + 1;
296  G4double ang = dphi/ksteps;
297 
298  G4double sinHalf = std::sin(0.5*ang);
299  G4double cosHalf = std::cos(0.5*ang);
300  G4double sinStep = 2.*sinHalf*cosHalf;
301  G4double cosStep = 1. - 2.*sinHalf*sinHalf;
302  G4double rext = rmax/cosHalf;
303 
304  // bounding envelope for full cylinder consists of two polygons,
305  // in other cases it is a sequence of quadrilaterals
306  if (rmin == 0 && dphi == twopi)
307  {
308  G4double sinCur = sinHalf;
309  G4double cosCur = cosHalf;
310 
311  G4ThreeVectorList baseA(NSTEPS),baseB(NSTEPS);
312  for (G4int k=0; k<NSTEPS; ++k)
313  {
314  baseA[k].set(rext*cosCur,rext*sinCur,-dz);
315  baseB[k].set(rext*cosCur,rext*sinCur, dz);
316 
317  G4double sinTmp = sinCur;
318  sinCur = sinCur*cosStep + cosCur*sinStep;
319  cosCur = cosCur*cosStep - sinTmp*sinStep;
320  }
321  std::vector<const G4ThreeVectorList *> polygons(2);
322  polygons[0] = &baseA;
323  polygons[1] = &baseB;
324  G4BoundingEnvelope benv(bmin,bmax,polygons);
325  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
326  }
327  else
328  {
329  G4double sinStart = GetSinStartPhi();
330  G4double cosStart = GetCosStartPhi();
331  G4double sinEnd = GetSinEndPhi();
332  G4double cosEnd = GetCosEndPhi();
333  G4double sinCur = sinStart*cosHalf + cosStart*sinHalf;
334  G4double cosCur = cosStart*cosHalf - sinStart*sinHalf;
335 
336  // set quadrilaterals
337  G4ThreeVectorList pols[NSTEPS+2];
338  for (G4int k=0; k<ksteps+2; ++k) pols[k].resize(4);
339  pols[0][0].set(rmin*cosStart,rmin*sinStart, dz);
340  pols[0][1].set(rmin*cosStart,rmin*sinStart,-dz);
341  pols[0][2].set(rmax*cosStart,rmax*sinStart,-dz);
342  pols[0][3].set(rmax*cosStart,rmax*sinStart, dz);
343  for (G4int k=1; k<ksteps+1; ++k)
344  {
345  pols[k][0].set(rmin*cosCur,rmin*sinCur, dz);
346  pols[k][1].set(rmin*cosCur,rmin*sinCur,-dz);
347  pols[k][2].set(rext*cosCur,rext*sinCur,-dz);
348  pols[k][3].set(rext*cosCur,rext*sinCur, dz);
349 
350  G4double sinTmp = sinCur;
351  sinCur = sinCur*cosStep + cosCur*sinStep;
352  cosCur = cosCur*cosStep - sinTmp*sinStep;
353  }
354  pols[ksteps+1][0].set(rmin*cosEnd,rmin*sinEnd, dz);
355  pols[ksteps+1][1].set(rmin*cosEnd,rmin*sinEnd,-dz);
356  pols[ksteps+1][2].set(rmax*cosEnd,rmax*sinEnd,-dz);
357  pols[ksteps+1][3].set(rmax*cosEnd,rmax*sinEnd, dz);
358 
359  // set envelope and calculate extent
360  std::vector<const G4ThreeVectorList *> polygons;
361  polygons.resize(ksteps+2);
362  for (G4int k=0; k<ksteps+2; ++k) polygons[k] = &pols[k];
363  G4BoundingEnvelope benv(bmin,bmax,polygons);
364  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
365  }
366  return exist;
367 }
368 
370 //
371 // Create polyhedron for visualization
372 //
373 G4Polyhedron* G4UTubs::CreatePolyhedron() const
374 {
375  return new G4PolyhedronTubs(GetInnerRadius(),
376  GetOuterRadius(),
377  GetZHalfLength(),
378  GetStartPhiAngle(),
379  GetDeltaPhiAngle());
380 }
381 
382 #endif // G4GEOM_USE_USOLIDS
void set(double x, double y, double z)
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
Definition: G4Tubs.hh:85
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
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
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
int G4int
Definition: G4Types.hh:78
EAxis
Definition: geomdefs.hh:54
double y() const
G4GLOB_DLL std::ostream G4cout
double x() const
Char_t n[5]
double y() const