source: trunk/source/geometry/solids/CSG/include/G4Cons.hh@ 1139

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

file release beta

File size: 7.2 KB
RevLine 
[831]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//
[921]27// $Id: G4Cons.hh,v 1.21 2008/11/06 11:04:00 gcosmo Exp $
[1058]28// GEANT4 tag $Name: geant4-09-02-ref-02 $
[831]29//
30//
31// --------------------------------------------------------------------
32// GEANT 4 class header file
33//
34// G4Cons
35//
36// Class description:
37//
38// A G4Cons is, in the general case, a Phi segment of a cone, with
39// half-length fDz, inner and outer radii specified at -fDz and +fDz.
40// The Phi segment is described by a starting fSPhi angle, and the
41// +fDPhi delta angle for the shape.
42// If the delta angle is >=2*pi, the shape is treated as continuous
43// in Phi
44//
45// Member Data:
46//
47// fRmin1 inside radius at -fDz
48// fRmin2 inside radius at +fDz
49// fRmax1 outside radius at -fDz
50// fRmax2 outside radius at +fDz
51// fDz half length in z
52//
53// fSPhi starting angle of the segment in radians
54// fDPhi delta angle of the segment in radians
55//
[921]56// fPhiFullCone Boolean variable used for indicate the Phi Section
57//
[831]58// Note:
59// Internally fSPhi & fDPhi are adjusted so that fDPhi<=2PI,
60// and fDPhi+fSPhi<=2PI. This enables simpler comparisons to be
61// made with (say) Phi of a point.
62
63// History:
64// 19.3.94 P.Kent: Old C++ code converted to tolerant geometry
65// 13.9.96 V.Grichine: Final modifications to commit
66// --------------------------------------------------------------------
67
68#ifndef G4Cons_HH
69#define G4Cons_HH
70
71#include "G4CSGSolid.hh"
72
73class G4Cons : public G4CSGSolid
74{
75 public: // with description
76
[921]77 G4Cons(const G4String& pName,
78 G4double pRmin1, G4double pRmax1,
79 G4double pRmin2, G4double pRmax2,
80 G4double pDz,
81 G4double pSPhi, G4double pDPhi);
[831]82
[921]83 virtual ~G4Cons() ;
[831]84
[921]85 // Accessors
[831]86
[921]87 inline G4double GetInnerRadiusMinusZ() const;
88 inline G4double GetOuterRadiusMinusZ() const;
89 inline G4double GetInnerRadiusPlusZ() const;
90 inline G4double GetOuterRadiusPlusZ() const;
[831]91
[921]92 inline G4double GetZHalfLength() const;
[831]93
[921]94 inline G4double GetStartPhiAngle () const;
95 inline G4double GetDeltaPhiAngle () const;
[831]96
[921]97 // Modifiers
[831]98
[921]99 inline void SetInnerRadiusMinusZ( G4double Rmin1 );
100 inline void SetOuterRadiusMinusZ( G4double Rmax1 );
101 inline void SetInnerRadiusPlusZ ( G4double Rmin2 );
102 inline void SetOuterRadiusPlusZ ( G4double Rmax2 );
[831]103
[921]104 inline void SetZHalfLength ( G4double newDz );
105 inline void SetStartPhiAngle ( G4double newSPhi);
106 inline void SetDeltaPhiAngle ( G4double newDPhi);
[831]107
[921]108 // Other methods for solid
[831]109
[921]110 inline G4double GetCubicVolume();
111 inline G4double GetSurfaceArea();
[831]112
[921]113 void ComputeDimensions(G4VPVParameterisation* p,
114 const G4int n,
115 const G4VPhysicalVolume* pRep);
[831]116
[921]117 G4bool CalculateExtent(const EAxis pAxis,
118 const G4VoxelLimits& pVoxelLimit,
119 const G4AffineTransform& pTransform,
120 G4double& pmin, G4double& pmax) const;
[831]121
[921]122 EInside Inside(const G4ThreeVector& p) const;
[831]123
[921]124 G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
[831]125
[921]126 G4double DistanceToIn (const G4ThreeVector& p,
127 const G4ThreeVector& v) const;
128 G4double DistanceToIn (const G4ThreeVector& p) const;
129 G4double DistanceToOut(const G4ThreeVector& p,
130 const G4ThreeVector& v,
131 const G4bool calcNorm=G4bool(false),
132 G4bool *validNorm=0,
133 G4ThreeVector *n=0) const;
134 G4double DistanceToOut(const G4ThreeVector& p) const;
[831]135
[921]136 G4GeometryType GetEntityType() const;
[831]137
[921]138 G4ThreeVector GetPointOnSurface() const;
[831]139
[921]140 std::ostream& StreamInfo(std::ostream& os) const;
[831]141
[921]142 // Visualisation functions
[831]143
[921]144 void DescribeYourselfTo( G4VGraphicsScene& scene ) const;
145 G4Polyhedron* CreatePolyhedron() const;
146 G4NURBS* CreateNURBS() const;
[831]147
148 public: // without description
149
[921]150 G4Cons(__void__&);
151 //
152 // Fake default constructor for usage restricted to direct object
153 // persistency for clients requiring preallocation of memory for
154 // persistifiable objects.
[831]155
[921]156 // Old access functions
[831]157
[921]158 inline G4double GetRmin1() const;
159 inline G4double GetRmax1() const;
160 inline G4double GetRmin2() const;
161 inline G4double GetRmax2() const;
[831]162
[921]163 inline G4double GetDz() const;
[831]164
[921]165 inline G4double GetSPhi() const;
166 inline G4double GetDPhi() const;
[831]167
168 protected:
169
[921]170 G4ThreeVectorList*
171 CreateRotatedVertices(const G4AffineTransform& pTransform) const;
[831]172
[921]173 G4double fRmin1, fRmin2, fRmax1, fRmax2, fDz, fSPhi, fDPhi;
174 G4bool fPhiFullCone;
175
176 // Used by distanceToOut
[831]177
[921]178 enum ESide {kNull,kRMin,kRMax,kSPhi,kEPhi,kPZ,kMZ};
[831]179
[921]180 // used by normal
[831]181
[921]182 enum ENorm {kNRMin,kNRMax,kNSPhi,kNEPhi,kNZ};
[831]183
184 private:
185
[921]186 inline void Initialise();
187 // Reset relevant values to zero
[831]188
[921]189 inline void InitializeTrigonometry();
190 //
191 // Recompute relevant trigonometric values and cache them
192
193 G4ThreeVector ApproxSurfaceNormal(const G4ThreeVector& p) const;
194 //
195 // Algorithm for SurfaceNormal() following the original
196 // specification for points not on the surface
197
[831]198 private:
199
[921]200 G4double kRadTolerance, kAngTolerance;
201 //
202 // Radial and angular tolerances
[831]203
[921]204 G4double sinCPhi, cosCPhi, cosHDPhiOT, cosHDPhiIT,
205 sinSPhi, cosSPhi, sinEPhi, cosEPhi;
206 //
207 // Cached trigonometric values
[831]208};
209
210#include "G4Cons.icc"
211
212#endif
Note: See TracBrowser for help on using the repository browser.