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