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