source: trunk/source/geometry/solids/BREPS/src/G4BREPSolidCone.cc @ 1246

Last change on this file since 1246 was 1228, checked in by garnier, 15 years ago

update geant4.9.3 tag

File size: 7.6 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// $Id: G4BREPSolidCone.cc,v 1.15 2006/06/29 18:41:16 gunter Exp $
27// GEANT4 tag $Name: geant4-09-03 $
28//
29// ----------------------------------------------------------------------
30// GEANT 4 class source file
31//
32// G4BREPSolidCone.cc
33//
34// ----------------------------------------------------------------------
35
36#include "G4BREPSolidCone.hh"
37#include "G4FPlane.hh"
38#include "G4FConicalSurface.hh"
39#include "G4FCylindricalSurface.hh"
40#include "G4CircularCurve.hh"
41
42G4BREPSolidCone::G4BREPSolidCone(const G4String& name,
43                                 const G4ThreeVector& origin,
44                                 const G4ThreeVector& axis,
45                                 const G4ThreeVector& direction,
46                                 G4double length,
47                                 G4double radius,
48                                 G4double large_radius)
49  : G4BREPSolid(name)
50{
51  SurfaceVec              = new G4Surface*[3];
52  G4Point3D    ArcStart1  = G4Point3D(origin + (radius*direction));
53  G4Vector3D   tmpaxis(axis);
54  G4Vector3D   tmporigin(origin); 
55  G4Point3D    tmppoint;
56
57  tmppoint= G4Point3D(origin) + (length*tmpaxis);
58  G4Point3D origin2(tmppoint.x(), tmppoint.y(), tmppoint.z());
59
60  tmppoint=  origin2 + (large_radius*tmpaxis);
61  G4Point3D ArcStart2(tmppoint.x(), tmppoint.y(), tmppoint.z());
62
63  G4Ray::Vcross(tmpaxis, axis, direction);
64  G4ThreeVector axis2(tmpaxis.x(),tmpaxis.y(), tmpaxis.z());
65
66  G4CurveVector CVec;
67  G4CircularCurve* tmp;
68
69  tmp = new G4CircularCurve();
70  tmp->Init(G4Axis2Placement3D(direction, axis2, origin) , large_radius);
71  tmp->SetBounds(ArcStart1, ArcStart1);
72  CVec.push_back(tmp);
73
74  tmp = new G4CircularCurve();
75  tmp->Init(G4Axis2Placement3D(direction, axis2, origin2), large_radius);
76  tmp->SetBounds(ArcStart2, ArcStart2);
77  CVec.push_back(tmp);
78
79  SurfaceVec[0] = new G4FConicalSurface(tmporigin, axis, 
80                                        length, radius, large_radius);
81  SurfaceVec[0]->SetBoundaries(&CVec);
82
83  // new G4AdvancedFace("G4FConicalSurface", tmporigin, direction,
84  // axis, CVec, 1, 0,0,length, radius, large_radius);
85
86  // Create end planes & boundaries for cone solid
87  G4CurveVector CVec2;
88  tmp = new G4CircularCurve();
89  tmp->Init(G4Axis2Placement3D(direction, axis2, origin), radius);
90  tmp->SetBounds(ArcStart1, ArcStart1);
91  CVec2.push_back(tmp);
92
93  SurfaceVec[1] = new G4FPlane(tmpaxis, direction, origin2);
94  //new G4AdvancedFace("G4FPlane" , origin2, direction, tmpaxis, CVec2, 1);
95  SurfaceVec[1]->SetBoundaries(&CVec2);
96
97  CVec2[0] = tmp = new G4CircularCurve();
98  tmp->Init(G4Axis2Placement3D(direction, axis2, origin2), large_radius);
99  tmp->SetBounds(ArcStart2, ArcStart2);
100
101  SurfaceVec[2] = new G4FPlane(tmpaxis, direction, origin);
102  //new G4AdvancedFace("G4FPlane", origin, direction, tmpaxis, CVec2, 1); 
103  SurfaceVec[2]->SetBoundaries(&CVec2);
104
105  nb_of_surfaces = 3;
106  active=1;
107 
108  // Save constructor parameters
109  constructorParams.origin       = origin;
110  constructorParams.axis         = axis;
111  constructorParams.direction    = direction;
112  constructorParams.length       = length;
113  constructorParams.radius       = radius;
114  constructorParams.large_radius = large_radius;
115 
116  Initialize();
117}
118
119G4BREPSolidCone::G4BREPSolidCone( __void__& a )
120  : G4BREPSolid(a)
121{
122}
123
124G4BREPSolidCone::~G4BREPSolidCone()
125{
126}
127
128void G4BREPSolidCone::Initialize()
129{
130  // Calc bounding box for solids and surfaces
131  // Convert concave planes to convex     
132  ShortestDistance=1000000;
133  CheckSurfaceNormals();
134  if(!Box || !AxisBox)
135    IsConvex();
136  CalcBBoxes();
137}
138
139EInside G4BREPSolidCone::Inside(register const G4ThreeVector& Pt) const
140{
141  G4double dist1 = SurfaceVec[0]->HowNear(Pt);
142  G4double dist2 = SurfaceVec[1]->ClosestDistanceToPoint(Pt);
143  G4double dist3 = SurfaceVec[2]->ClosestDistanceToPoint(Pt); 
144  if(dist1 > dist2) dist1 = dist2;
145  if(dist1 > dist3) dist1 = dist3; 
146  if(dist1 > 0) return kInside;
147  if(dist1 < 0) return kOutside;
148  return kSurface;
149}
150
151G4ThreeVector G4BREPSolidCone::SurfaceNormal(const G4ThreeVector& Pt) const
152{
153  G4Vector3D n =  SurfaceVec[0]->Normal(Pt);
154  G4ThreeVector norm(n.x(), n.y(), n.z());
155  return norm;
156}
157
158G4double G4BREPSolidCone::DistanceToIn(const G4ThreeVector& Pt) const
159{
160  G4double dist1 = std::fabs(SurfaceVec[0]->HowNear(Pt));
161  G4double dist2 = std::fabs(SurfaceVec[1]->ClosestDistanceToPoint(Pt));
162  G4double dist3 = std::fabs(SurfaceVec[2]->ClosestDistanceToPoint(Pt)); 
163  if(dist1 > dist2) dist1 = dist2;
164  if(dist1 > dist3) dist1 = dist3; 
165  return dist1;
166 
167}
168
169G4double G4BREPSolidCone::DistanceToIn(register const G4ThreeVector& Pt, 
170                                       register const G4ThreeVector& V) const
171{
172  Reset(); 
173  G4Vector3D Pttmp(Pt);
174  G4Vector3D Vtmp(V);   
175  //  G4double kInfinity = 10e20;
176  G4Ray r(Pttmp, Vtmp);
177
178  if(SurfaceVec[0]->Intersect( r ))
179  {
180    ShortestDistance = SurfaceVec[0]->GetDistance();
181    return ShortestDistance;
182  }
183  return kInfinity; 
184}
185
186G4double G4BREPSolidCone::DistanceToOut(register const G4ThreeVector& Pt, 
187                                        register const G4ThreeVector& V, 
188                                        const G4bool, 
189                                        G4bool *validNorm, 
190                                        G4ThreeVector *) const
191{
192  if(validNorm)
193    *validNorm = false;
194  Reset(); 
195
196  G4Vector3D Pttmp(Pt);
197  G4Vector3D Vtmp(V);   
198  //  G4double kInfinity = 10e20;
199
200  G4Ray r(Pttmp, Vtmp);
201  if(SurfaceVec[0]->Intersect( r ))
202  {
203    ShortestDistance = SurfaceVec[0]->GetDistance();
204    return ShortestDistance;
205  }
206  return kInfinity; 
207}
208
209G4double G4BREPSolidCone::DistanceToOut(const G4ThreeVector& Pt) const
210{
211  G4double dist1 = std::fabs(SurfaceVec[0]->HowNear(Pt));
212  G4double dist2 = std::fabs(SurfaceVec[1]->ClosestDistanceToPoint(Pt));
213  G4double dist3 = std::fabs(SurfaceVec[2]->ClosestDistanceToPoint(Pt)); 
214  if(dist1 > dist2) dist1 = dist2;
215  if(dist1 > dist3) dist1 = dist3; 
216  return dist1;
217}
218
219// Streams solid contents to output stream.
220std::ostream& G4BREPSolidCone::StreamInfo(std::ostream& os) const
221{
222  G4BREPSolid::StreamInfo( os )
223  << "\n origin:       " << constructorParams.origin
224  << "\n axis:         " << constructorParams.axis
225  << "\n direction:    " << constructorParams.direction
226  << "\n length:       " << constructorParams.length
227  << "\n radius:       " << constructorParams.radius
228  << "\n large_radius: " << constructorParams.large_radius
229  << "\n-----------------------------------------------------------\n";
230
231  return os;
232}
233
Note: See TracBrowser for help on using the repository browser.