Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4UOrb.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 G4UOrb wrapper class
27//
28// 30.10.13 G.Cosmo, CERN
29// --------------------------------------------------------------------
30
31#include "G4Orb.hh"
32#include "G4UOrb.hh"
33
34#if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
35
36#include "G4TwoVector.hh"
37#include "G4AffineTransform.hh"
39#include "G4BoundingEnvelope.hh"
40
43
44using namespace CLHEP;
45
46////////////////////////////////////////////////////////////////////////
47//
48// constructor - check positive radius
49//
50
51G4UOrb::G4UOrb( const G4String& pName, G4double pRmax )
52 : Base_t(pName, pRmax)
53{
54}
55
56//////////////////////////////////////////////////////////////////////////
57//
58// Copy constructor
59
60G4UOrb::G4UOrb(const G4UOrb& rhs)
61 : Base_t(rhs)
62{
63}
64
65//////////////////////////////////////////////////////////////////////////
66//
67// Assignment operator
68
69G4UOrb& G4UOrb::operator = (const G4UOrb& rhs)
70{
71 // Check assignment to self
72 //
73 if (this == &rhs) { return *this; }
74
75 // Copy base class data
76 //
77 Base_t::operator=(rhs);
78
79 return *this;
80}
81
82//////////////////////////////////////////////////////////////////////////
83//
84// Accessors & modifiers
85
86G4double G4UOrb::GetRadius() const
87{
88 return Base_t::GetRadius();
89}
90
91void G4UOrb::SetRadius(G4double newRmax)
92{
93 Base_t::SetRadius(newRmax);
94 fRebuildPolyhedron = true;
95}
96
97G4double G4UOrb::GetRadialTolerance() const
98{
99 return Base_t::GetRadialTolerance();
100}
101
102//////////////////////////////////////////////////////////////////////////
103//
104// Dispatch to parameterisation for replication mechanism dimension
105// computation & modification.
106
107void G4UOrb::ComputeDimensions( G4VPVParameterisation* p,
108 const G4int n,
109 const G4VPhysicalVolume* pRep )
110{
111 p->ComputeDimensions(*(G4Orb*)this,n,pRep);
112}
113
114//////////////////////////////////////////////////////////////////////////
115//
116// Make a clone of the object
117
118G4VSolid* G4UOrb::Clone() const
119{
120 return new G4UOrb(*this);
121}
122
123//////////////////////////////////////////////////////////////////////////
124//
125// Get bounding box
126
127void G4UOrb::BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const
128{
129 G4double radius = GetRadius();
130 pMin.set(-radius,-radius,-radius);
131 pMax.set( radius, radius, radius);
132
133 // Check correctness of the bounding box
134 //
135 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
136 {
137 std::ostringstream message;
138 message << "Bad bounding box (min >= max) for solid: "
139 << GetName() << " !"
140 << "\npMin = " << pMin
141 << "\npMax = " << pMax;
142 G4Exception("G4UOrb::BoundingLimits()", "GeomMgt0001",
143 JustWarning, message);
144 StreamInfo(G4cout);
145 }
146}
147
148//////////////////////////////////////////////////////////////////////////
149//
150// Calculate extent under transform and specified limit
151
152G4bool
153G4UOrb::CalculateExtent(const EAxis pAxis,
154 const G4VoxelLimits& pVoxelLimit,
155 const G4AffineTransform& pTransform,
156 G4double& pMin, G4double& pMax) const
157{
158 G4ThreeVector bmin, bmax;
159 G4bool exist;
160
161 // Get bounding box
162 BoundingLimits(bmin,bmax);
163
164 // Check bounding box
165 G4BoundingEnvelope bbox(bmin,bmax);
166#ifdef G4BBOX_EXTENT
167 if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
168#endif
169 if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
170 {
171 return exist = pMin < pMax;
172 }
173
174 // Find bounding envelope and calculate extent
175 //
176 static const G4int NTHETA = 8; // number of steps along Theta
177 static const G4int NPHI = 16; // number of steps along Phi
178 static const G4double sinHalfTheta = std::sin(halfpi/NTHETA);
179 static const G4double cosHalfTheta = std::cos(halfpi/NTHETA);
180 static const G4double sinHalfPhi = std::sin(pi/NPHI);
181 static const G4double cosHalfPhi = std::cos(pi/NPHI);
182 static const G4double sinStepTheta = 2.*sinHalfTheta*cosHalfTheta;
183 static const G4double cosStepTheta = 1. - 2.*sinHalfTheta*sinHalfTheta;
184 static const G4double sinStepPhi = 2.*sinHalfPhi*cosHalfPhi;
185 static const G4double cosStepPhi = 1. - 2.*sinHalfPhi*sinHalfPhi;
186
187 G4double radius = GetRadius();
188 G4double rtheta = radius/cosHalfTheta;
189 G4double rphi = rtheta/cosHalfPhi;
190
191 // set reference circle
192 G4TwoVector xy[NPHI];
193 G4double sinCurPhi = sinHalfPhi;
194 G4double cosCurPhi = cosHalfPhi;
195 for (auto & k : xy)
196 {
197 k.set(cosCurPhi,sinCurPhi);
198 G4double sinTmpPhi = sinCurPhi;
199 sinCurPhi = sinCurPhi*cosStepPhi + cosCurPhi*sinStepPhi;
200 cosCurPhi = cosCurPhi*cosStepPhi - sinTmpPhi*sinStepPhi;
201 }
202
203 // set bounding circles
204 G4ThreeVectorList circles[NTHETA];
205 for (auto & circle : circles) circle.resize(NPHI);
206
207 G4double sinCurTheta = sinHalfTheta;
208 G4double cosCurTheta = cosHalfTheta;
209 for (auto & circle : circles)
210 {
211 G4double z = rtheta*cosCurTheta;
212 G4double rho = rphi*sinCurTheta;
213 for (G4int k=0; k<NPHI; ++k)
214 {
215 circle[k].set(rho*xy[k].x(),rho*xy[k].y(),z);
216 }
217 G4double sinTmpTheta = sinCurTheta;
218 sinCurTheta = sinCurTheta*cosStepTheta + cosCurTheta*sinStepTheta;
219 cosCurTheta = cosCurTheta*cosStepTheta - sinTmpTheta*sinStepTheta;
220 }
221
222 // set envelope and calculate extent
223 std::vector<const G4ThreeVectorList *> polygons;
224 polygons.resize(NTHETA);
225 for (G4int i=0; i<NTHETA; ++i) polygons[i] = &circles[i];
226
227 G4BoundingEnvelope benv(bmin,bmax,polygons);
228 exist = benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
229 return exist;
230}
231
232//////////////////////////////////////////////////////////////////////////
233//
234// Create polyhedron for visualization
235
236G4Polyhedron* G4UOrb::CreatePolyhedron() const
237{
238 return new G4PolyhedronSphere(0., GetRadius(), 0., twopi, 0., pi);
239}
240
241#endif // G4GEOM_USE_USOLIDS
std::vector< G4ThreeVector > G4ThreeVectorList
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
CLHEP::Hep2Vector G4TwoVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
void set(double x, double y, double z)
G4AffineTransform is a class for geometric affine transformations. It supports efficient arbitrary ro...
G4BoundingEnvelope is a helper class to facilitate calculation of the extent of a solid within the li...
G4Orb represents a full sphere.
Definition G4Orb.hh:59
G4VPVParameterisation ia an abstract base class for Parameterisation, able to compute the transformat...
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
G4VoxelLimits represents limitation/restrictions of space, where restrictions are only made perpendic...
EAxis
Definition geomdefs.hh:54