source: trunk/source/geometry/solids/BREPS/src/G4BREPSolidSphere.cc@ 1202

Last change on this file since 1202 was 1058, checked in by garnier, 17 years ago

file release beta

File size: 4.8 KB
Line 
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: G4BREPSolidSphere.cc,v 1.11 2006/06/29 18:41:32 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// ----------------------------------------------------------------------
31// GEANT 4 class source file
32//
33// G4BREPSolidSphere.cc
34//
35// ----------------------------------------------------------------------
36
37#include "G4BREPSolidSphere.hh"
38#include "G4SphericalSurface.hh"
39
40G4BREPSolidSphere::G4BREPSolidSphere(const G4String& name,
41 const G4Vector3D& origin,
42 const G4Vector3D& xhat,
43 const G4Vector3D& zhat,
44 G4double radius)
45 : G4BREPSolid(name)
46{
47 SurfaceVec = new G4Surface*[1];
48 G4double ph1 = 0;
49 G4double ph2 = 2*pi;
50 G4double th1 = 0;
51 G4double th2 = pi;
52 SurfaceVec[0] = new G4SphericalSurface(origin, xhat, zhat, radius, ph1, ph2, th1, th2);
53 nb_of_surfaces = 1;
54
55 constructorParams.origin = origin;
56 constructorParams.xhat = xhat;
57 constructorParams.zhat = zhat;
58 constructorParams.radius = radius;
59
60 active=1;
61 Initialize();
62}
63
64G4BREPSolidSphere::G4BREPSolidSphere( __void__& a )
65 : G4BREPSolid(a)
66{
67}
68
69G4BREPSolidSphere::~G4BREPSolidSphere()
70{
71}
72
73EInside G4BREPSolidSphere::Inside(register const G4ThreeVector& Pt) const
74{
75 G4double Dist = SurfaceVec[0]->HowNear(Pt);
76 if(Dist > 0+kCarTolerance) return kInside;
77 if(Dist < 0-kCarTolerance) return kOutside;
78 return kSurface;
79}
80
81
82G4ThreeVector G4BREPSolidSphere::SurfaceNormal(const G4ThreeVector& Pt) const
83{
84 G4Vector3D n = SurfaceVec[0]->Normal(Pt);
85 G4ThreeVector norm(n.x(), n.y(), n.z());
86 return norm;
87}
88
89
90G4double G4BREPSolidSphere::DistanceToIn(const G4ThreeVector& Pt) const
91{
92 return std::fabs(SurfaceVec[0]->HowNear(Pt));
93}
94
95
96G4double G4BREPSolidSphere::DistanceToIn(register const G4ThreeVector& Pt,
97 register const G4ThreeVector& V) const
98{
99 // SphReset();
100 G4Vector3D Pttmp(Pt);
101 G4Vector3D Vtmp(V);
102 G4Ray r(Pttmp, Vtmp);
103 G4int Result = SurfaceVec[0]->Intersect( r );
104
105 if(Result>0)
106 {
107 ShortestDistance = SurfaceVec[0]->GetDistance();
108 return std::sqrt(ShortestDistance);
109 }
110 return kInfinity;
111}
112
113
114G4double G4BREPSolidSphere::DistanceToOut(register const G4ThreeVector& Pt,
115 register const G4ThreeVector& V,
116 const G4bool calcNorm,
117 G4bool *validNorm,
118 G4ThreeVector *n) const
119{
120 if(validNorm) *validNorm = false;
121 // SphReset();
122 G4Vector3D Pttmp(Pt);
123 G4Vector3D Vtmp(V);
124 G4Ray r(Pttmp, Vtmp);
125
126 if(SurfaceVec[0]->Intersect( r ))
127 {
128 if(calcNorm)
129 {
130 *validNorm = true;
131 *n = SurfaceNormal(Pt);
132 }
133
134 ShortestDistance = SurfaceVec[0]->GetDistance();
135 return std::sqrt(ShortestDistance);
136 }
137 return kInfinity;
138}
139
140
141G4double G4BREPSolidSphere::DistanceToOut(const G4ThreeVector& Pt) const
142{
143 return std::fabs(SurfaceVec[0]->HowNear(Pt));
144}
145
146// Streams solid contents to output stream.
147std::ostream& G4BREPSolidSphere::StreamInfo(std::ostream& os) const
148{
149 G4BREPSolid::StreamInfo( os )
150 << "\n origin: " << constructorParams.origin
151 << "\n xhat: " << constructorParams.xhat
152 << "\n zhat: " << constructorParams.zhat
153 << "\n radius: " << constructorParams.radius
154 << "\n-----------------------------------------------------------\n";
155
156 return os;
157}
158
Note: See TracBrowser for help on using the repository browser.