Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UPara.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 G4UPara wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4Para.hh"
34 #include "G4UPara.hh"
35 
36 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
37 
38 #include "G4AffineTransform.hh"
39 #include "G4VPVParameterisation.hh"
40 #include "G4BoundingEnvelope.hh"
41 
42 using namespace CLHEP;
43 
45 //
46 // Constructor - set & check half widths
47 
48 G4UPara::G4UPara(const G4String& pName,
49  G4double pDx, G4double pDy, G4double pDz,
50  G4double pAlpha, G4double pTheta, G4double pPhi)
51  : Base_t(pName, pDx, pDy, pDz, pAlpha, pTheta, pPhi)
52 {
53  fTalpha = std::tan(pAlpha);
54  fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
55  fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
56  CheckParameters();
57  MakePlanes();
58 }
59 
61 //
62 // Constructor - design of trapezoid based on 8 vertices
63 
64 G4UPara::G4UPara( const G4String& pName,
65  const G4ThreeVector pt[8] )
66  : Base_t(pName)
67 {
68  // Find dimensions and trigonometric values
69  //
70  G4double fDx = (pt[3].x() - pt[2].x())*0.5;
71  G4double fDy = (pt[2].y() - pt[1].y())*0.5;
72  G4double fDz = pt[7].z();
73  SetDimensions(fDx, fDy, fDz);
74  CheckParameters(); // check dimensions
75 
76  fTalpha = (pt[2].x() + pt[3].x() - pt[1].x() - pt[0].x())*0.25/fDy;
77  fTthetaCphi = (pt[4].x() + fDy*fTalpha + fDx)/fDz;
78  fTthetaSphi = (pt[4].y() + fDy)/fDz;
79  SetAlpha(std::atan(fTalpha));
80  SetTheta(std::atan(std::sqrt(fTthetaSphi*fTthetaSphi
81  + fTthetaCphi*fTthetaCphi)));
82  SetPhi (std::atan2(fTthetaSphi, fTthetaCphi));
83  MakePlanes();
84 
85  // Recompute vertices
86  //
87  G4ThreeVector v[8];
88  G4double DyTalpha = fDy*fTalpha;
89  G4double DzTthetaSphi = fDz*fTthetaSphi;
90  G4double DzTthetaCphi = fDz*fTthetaCphi;
91  v[0].set(-DzTthetaCphi-DyTalpha-fDx, -DzTthetaSphi-fDy, -fDz);
92  v[1].set(-DzTthetaCphi-DyTalpha+fDx, -DzTthetaSphi-fDy, -fDz);
93  v[2].set(-DzTthetaCphi+DyTalpha-fDx, -DzTthetaSphi+fDy, -fDz);
94  v[3].set(-DzTthetaCphi+DyTalpha+fDx, -DzTthetaSphi+fDy, -fDz);
95  v[4].set( DzTthetaCphi-DyTalpha-fDx, DzTthetaSphi-fDy, fDz);
96  v[5].set( DzTthetaCphi-DyTalpha+fDx, DzTthetaSphi-fDy, fDz);
97  v[6].set( DzTthetaCphi+DyTalpha-fDx, DzTthetaSphi+fDy, fDz);
98  v[7].set( DzTthetaCphi+DyTalpha+fDx, DzTthetaSphi+fDy, fDz);
99 
100  // Compare with original vertices
101  //
102  for (G4int i=0; i<8; ++i)
103  {
104  G4double delx = std::abs(pt[i].x() - v[i].x());
105  G4double dely = std::abs(pt[i].y() - v[i].y());
106  G4double delz = std::abs(pt[i].z() - v[i].z());
107  G4double discrepancy = std::max(std::max(delx,dely),delz);
108  if (discrepancy > 0.1*kCarTolerance)
109  {
110  std::ostringstream message;
111  G4int oldprc = message.precision(16);
112  message << "Invalid vertice coordinates for Solid: " << GetName()
113  << "\nVertix #" << i << ", discrepancy = " << discrepancy
114  << "\n original : " << pt[i]
115  << "\n recomputed : " << v[i];
116  G4cout.precision(oldprc);
117  G4Exception("G4UPara::G4UPara()", "GeomSolids0002",
118  FatalException, message);
119 
120  }
121  }
122 }
123 
125 //
126 // Fake default constructor - sets only member data and allocates memory
127 // for usage restricted to object persistency
128 
129 G4UPara::G4UPara( __void__& a )
130  : Base_t(a)
131 {
132  SetAllParameters(1., 1., 1., 0., 0., 0.);
133  fRebuildPolyhedron = false;
134 }
135 
137 //
138 // Destructor
139 
140 G4UPara::~G4UPara()
141 {
142 }
143 
145 //
146 // Copy constructor
147 
148 G4UPara::G4UPara(const G4UPara& rhs)
149  : Base_t(rhs), fTalpha(rhs.fTalpha),
150  fTthetaCphi(rhs.fTthetaCphi),fTthetaSphi(rhs.fTthetaSphi)
151 {
152  for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
153 }
154 
156 //
157 // Assignment operator
158 
159 G4UPara& G4UPara::operator = (const G4UPara& rhs)
160 {
161  // Check assignment to self
162  //
163  if (this == &rhs) { return *this; }
164 
165  // Copy base class data
166  //
167  Base_t::operator=(rhs);
168 
169  // Copy data
170  //
171  fTalpha = rhs.fTalpha;
172  fTthetaCphi = rhs.fTthetaCphi;
173  fTthetaSphi = rhs.fTthetaSphi;
174  for (G4int i=0; i<4; ++i) { fPlanes[i] = rhs.fPlanes[i]; }
175 
176  return *this;
177 }
178 
180 //
181 // Accessors & modifiers
182 
183 G4double G4UPara::GetZHalfLength() const
184 {
185  return GetZ();
186 }
187 G4double G4UPara::GetYHalfLength() const
188 {
189  return GetY();
190 }
191 G4double G4UPara::GetXHalfLength() const
192 {
193  return GetZ();
194 }
195 G4ThreeVector G4UPara::GetSymAxis() const
196 {
197  return G4ThreeVector(fTthetaCphi,fTthetaSphi,1.).unit();
198 }
199 G4double G4UPara::GetTanAlpha() const
200 {
201  return fTalpha;
202 }
203 
204 void G4UPara::SetXHalfLength(G4double val)
205 {
206  SetDimensions(val, GetY(), GetZ());
207  fRebuildPolyhedron = true;
208 
209  CheckParameters();
210  MakePlanes();
211 }
212 void G4UPara::SetYHalfLength(G4double val)
213 {
214  SetDimensions(GetX(), val, GetZ());
215  fRebuildPolyhedron = true;
216 
217  CheckParameters();
218  MakePlanes();
219 }
220 void G4UPara::SetZHalfLength(G4double val)
221 {
222  SetDimensions(GetX(), GetY(), val);
223  fRebuildPolyhedron = true;
224 
225  CheckParameters();
226  MakePlanes();
227 }
228 void G4UPara::SetAlpha(G4double alpha)
229 {
230  Base_t::SetAlpha(alpha);
231  fTalpha = std::tan(alpha);
232  fRebuildPolyhedron = true;
233 
234  MakePlanes();
235 }
236 void G4UPara::SetTanAlpha(G4double val)
237 {
238  fTalpha = val;
239  fRebuildPolyhedron = true;
240 
241  MakePlanes();
242 }
243 void G4UPara::SetThetaAndPhi(double pTheta, double pPhi)
244 {
245  Base_t::SetThetaAndPhi(pTheta, pPhi);
246  G4double tanTheta = std::tan(pTheta);
247  fTthetaCphi = tanTheta*std::cos(pPhi);
248  fTthetaSphi = tanTheta*std::sin(pPhi);
249  fRebuildPolyhedron = true;
250 
251  MakePlanes();
252 }
253 
255 //
256 // Set all parameters, as for constructor - set and check half-widths
257 
258 void G4UPara::SetAllParameters(G4double pDx, G4double pDy, G4double pDz,
259  G4double pAlpha, G4double pTheta, G4double pPhi)
260 {
261  // Reset data of the base class
262  fRebuildPolyhedron = true;
263 
264  // Set parameters
265  SetDimensions(pDx, pDy, pDz);
266  Base_t::SetAlpha(pAlpha);
267  Base_t::SetThetaAndPhi(pTheta, pPhi);
268  fTalpha = std::tan(pAlpha);
269  fTthetaCphi = std::tan(pTheta)*std::cos(pPhi);
270  fTthetaSphi = std::tan(pTheta)*std::sin(pPhi);
271 
272  CheckParameters();
273  MakePlanes();
274 }
275 
277 //
278 // Check dimensions
279 
280 void G4UPara::CheckParameters()
281 {
282  if (GetX() < 2*kCarTolerance ||
283  GetY() < 2*kCarTolerance ||
284  GetZ() < 2*kCarTolerance)
285  {
286  std::ostringstream message;
287  message << "Invalid (too small or negative) dimensions for Solid: "
288  << GetName()
289  << "\n X - " << GetX()
290  << "\n Y - " << GetY()
291  << "\n Z - " << GetZ();
292  G4Exception("G4UPara::CheckParameters()", "GeomSolids0002",
293  FatalException, message);
294  }
295 }
296 
298 //
299 // Set side planes
300 
301 void G4UPara::MakePlanes()
302 {
303  G4ThreeVector vx(1, 0, 0);
304  G4ThreeVector vy(fTalpha, 1, 0);
305  G4ThreeVector vz(fTthetaCphi, fTthetaSphi, 1);
306 
307  // Set -Y & +Y planes
308  //
309  G4ThreeVector ynorm = (vx.cross(vz)).unit();
310 
311  fPlanes[0].a = 0.;
312  fPlanes[0].b = ynorm.y();
313  fPlanes[0].c = ynorm.z();
314  fPlanes[0].d = fPlanes[0].b*GetY(); // point (0,fDy,0) is on plane
315 
316  fPlanes[1].a = 0.;
317  fPlanes[1].b = -fPlanes[0].b;
318  fPlanes[1].c = -fPlanes[0].c;
319  fPlanes[1].d = fPlanes[0].d;
320 
321  // Set -X & +X planes
322  //
323  G4ThreeVector xnorm = (vz.cross(vy)).unit();
324 
325  fPlanes[2].a = xnorm.x();
326  fPlanes[2].b = xnorm.y();
327  fPlanes[2].c = xnorm.z();
328  fPlanes[2].d = fPlanes[2].a*GetZ(); // point (fDx,0,0) is on plane
329 
330  fPlanes[3].a = -fPlanes[2].a;
331  fPlanes[3].b = -fPlanes[2].b;
332  fPlanes[3].c = -fPlanes[2].c;
333  fPlanes[3].d = fPlanes[2].d;
334 }
335 
337 //
338 // Dispatch to parameterisation for replication mechanism dimension
339 // computation & modification
340 
341 void G4UPara::ComputeDimensions( G4VPVParameterisation* p,
342  const G4int n,
343  const G4VPhysicalVolume* pRep )
344 {
345  p->ComputeDimensions(*(G4Para*)this,n,pRep);
346 }
347 
349 //
350 // Get bounding box
351 
352 void G4UPara::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
353 {
354  G4double dz = GetZHalfLength();
355  G4double dx = GetXHalfLength();
356  G4double dy = GetYHalfLength();
357 
358  G4double x0 = dz*fTthetaCphi;
359  G4double x1 = dy*GetTanAlpha();
360  G4double xmin =
361  std::min(
362  std::min(
363  std::min(-x0-x1-dx,-x0+x1-dx),x0-x1-dx),x0+x1-dx);
364  G4double xmax =
365  std::max(
366  std::max(
367  std::max(-x0-x1+dx,-x0+x1+dx),x0-x1+dx),x0+x1+dx);
368 
369  G4double y0 = dz*fTthetaSphi;
370  G4double ymin = std::min(-y0-dy,y0-dy);
371  G4double ymax = std::max(-y0+dy,y0+dy);
372 
373  pMin.set(xmin,ymin,-dz);
374  pMax.set(xmax,ymax, dz);
375 
376  // Check correctness of the bounding box
377  //
378  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
379  {
380  std::ostringstream message;
381  message << "Bad bounding box (min >= max) for solid: "
382  << GetName() << " !"
383  << "\npMin = " << pMin
384  << "\npMax = " << pMax;
385  G4Exception("G4UPara::BoundingLimits()", "GeomMgt0001",
386  JustWarning, message);
387  StreamInfo(G4cout);
388  }
389 }
390 
392 //
393 // Calculate extent under transform and specified limit
394 
395 G4bool G4UPara::CalculateExtent( const EAxis pAxis,
396  const G4VoxelLimits& pVoxelLimit,
397  const G4AffineTransform& pTransform,
398  G4double& pMin, G4double& pMax ) const
399 {
400  G4ThreeVector bmin, bmax;
401  G4bool exist;
402 
403  // Check bounding box (bbox)
404  //
405  BoundingLimits(bmin,bmax);
406  G4BoundingEnvelope bbox(bmin,bmax);
407 #ifdef G4BBOX_EXTENT
408  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
409 #endif
410  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
411  {
412  return exist = (pMin < pMax) ? true : false;
413  }
414 
415  // Set bounding envelope (benv) and calculate extent
416  //
417  G4double dz = GetZHalfLength();
418  G4double dx = GetXHalfLength();
419  G4double dy = GetYHalfLength();
420 
421  G4double x0 = dz*fTthetaCphi;
422  G4double x1 = dy*GetTanAlpha();
423  G4double y0 = dz*fTthetaSphi;
424 
425  G4ThreeVectorList baseA(4), baseB(4);
426  baseA[0].set(-x0-x1-dx,-y0-dy,-dz);
427  baseA[1].set(-x0-x1+dx,-y0-dy,-dz);
428  baseA[2].set(-x0+x1+dx,-y0+dy,-dz);
429  baseA[3].set(-x0+x1-dx,-y0+dy,-dz);
430 
431  baseB[0].set(+x0-x1-dx, y0-dy, dz);
432  baseB[1].set(+x0-x1+dx, y0-dy, dz);
433  baseB[2].set(+x0+x1+dx, y0+dy, dz);
434  baseB[3].set(+x0+x1-dx, y0+dy, dz);
435 
436  std::vector<const G4ThreeVectorList *> polygons(2);
437  polygons[0] = &baseA;
438  polygons[1] = &baseB;
439 
440  G4BoundingEnvelope benv(bmin,bmax,polygons);
441  exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
442  return exist;
443 }
444 
446 //
447 // Make a clone of the object
448 //
449 G4VSolid* G4UPara::Clone() const
450 {
451  return new G4UPara(*this);
452 }
453 
455 //
456 // Methods for visualisation
457 
458 G4Polyhedron* G4UPara::CreatePolyhedron () const
459 {
460  return new G4PolyhedronPara(GetX(), GetY(), GetZ(),
461  GetAlpha(), GetTheta(), GetPhi());
462 }
463 
464 #endif // G4GEOM_USE_USOLIDS
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
CLHEP::Hep3Vector G4ThreeVector
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
Float_t x1[n_points_granero]
Definition: compare.C:5
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
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
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
std::vector< G4ThreeVector > G4ThreeVectorList
static const G4double alpha
TMarker * pt
Definition: egs.C:25
Hep3Vector unit() const
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
G4GLOB_DLL std::ostream G4cout
double x() const
Char_t n[5]
double y() const
Definition: G4Para.hh:86
T min(const T t1, const T t2)
brief Return the smallest of the two arguments