Geant4-11
G4IntersectionSolid.hh
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// G4IntersectionSolid
27//
28// Class description:
29//
30// Class for description of intersection of two solids.
31
32// 12.09.98 V.Grichine: initial design and implementation
33// --------------------------------------------------------------------
34#ifndef G4INTERSECTIONSOLID_HH
35#define G4INTERSECTIONSOLID_HH
36
37#include "G4BooleanSolid.hh"
38#include "G4VSolid.hh"
39
40#include "G4RotationMatrix.hh"
41#include "G4ThreeVector.hh"
42#include "G4Transform3D.hh"
43#include "G4AffineTransform.hh"
44
46{
47 public: // with description
48
49 G4IntersectionSolid( const G4String& pName,
50 G4VSolid* pSolidA ,
51 G4VSolid* pSolidB ) ;
52
53 G4IntersectionSolid( const G4String& pName,
54 G4VSolid* pSolidA ,
55 G4VSolid* pSolidB,
56 G4RotationMatrix* rotMatrix,
57 const G4ThreeVector& transVector ) ;
58
59 G4IntersectionSolid( const G4String& pName,
60 G4VSolid* pSolidA ,
61 G4VSolid* pSolidB,
62 const G4Transform3D& transform ) ;
63
64 virtual ~G4IntersectionSolid() ;
65
67
68 G4VSolid* Clone() const;
69
70 public: // without description
71
72 G4IntersectionSolid(__void__&);
73 // Fake default constructor for usage restricted to direct object
74 // persistency for clients requiring preallocation of memory for
75 // persistifiable objects.
76
79 // Copy constructor and assignment operator.
80
82
83 G4bool CalculateExtent( const EAxis pAxis,
84 const G4VoxelLimits& pVoxelLimit,
85 const G4AffineTransform& pTransform,
86 G4double& pMin, G4double& pMax) const ;
87
88 EInside Inside( const G4ThreeVector& p ) const ;
89
90 G4ThreeVector SurfaceNormal( const G4ThreeVector& p ) const ;
91
93 const G4ThreeVector& v ) const ;
94
95 G4double DistanceToIn( const G4ThreeVector& p) const ;
96
98 const G4ThreeVector& v,
99 const G4bool calcNorm=false,
100 G4bool *validNorm=0,
101 G4ThreeVector *n=0 ) const ;
102
103 G4double DistanceToOut( const G4ThreeVector& p ) const ;
104
105
107 const G4int n,
108 const G4VPhysicalVolume* pRep ) ;
109
110 void DescribeYourselfTo ( G4VGraphicsScene& scene ) const ;
112};
113
114#endif
115
static const G4double pMax
static const G4double pMin
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
G4VSolid * Clone() const
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4IntersectionSolid & operator=(const G4IntersectionSolid &rhs)
G4IntersectionSolid(const G4String &pName, G4VSolid *pSolidA, G4VSolid *pSolidB)
G4GeometryType GetEntityType() const
G4Polyhedron * CreatePolyhedron() const
EInside Inside(const G4ThreeVector &p) const
void DescribeYourselfTo(G4VGraphicsScene &scene) const
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
EAxis
Definition: geomdefs.hh:54
EInside
Definition: geomdefs.hh:67
G4bool transform(G4String &input, const G4String &type)