Geant4  10.01.p01
G4UGenericTrap.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 //
27 // $Id:$
28 //
29 //
30 // Implementation of G4UGenericTrap wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4GenericTrap.hh"
34 #include "G4UGenericTrap.hh"
35 #include "G4Polyhedron.hh"
36 #include "G4PolyhedronArbitrary.hh"
37 
39 //
40 // Constructor (generic parameters)
41 //
43  const std::vector<G4TwoVector>& vertices)
44  : G4USolid(name, new UGenericTrap())
45 {
46  SetZHalfLength(halfZ);
47  std::vector<UVector2> v;
48  for (size_t n=0; n<vertices.size(); ++n)
49  {
50  v.push_back(UVector2(vertices[n].x(),vertices[n].y()));
51  }
52  GetShape()->SetName(name);
53  GetShape()->Initialise(v);
54 }
55 
56 
58 //
59 // Fake default constructor - sets only member data and allocates memory
60 // for usage restricted to object persistency.
61 //
63  : G4USolid(a)
64 {
65 }
66 
67 
69 //
70 // Destructor
71 //
73 {
74 }
75 
76 
78 //
79 // Copy constructor
80 //
82  : G4USolid(source)
83 {
84 }
85 
86 
88 //
89 // Assignment operator
90 //
93 {
94  if (this == &source) return *this;
95 
96  G4USolid::operator=( source );
97 
98  return *this;
99 }
100 
102 //
103 // CreatePolyhedron()
104 //
106 {
107  // Approximation of Twisted Side
108  // Construct extra Points, if Twisted Side
109  //
110  G4PolyhedronArbitrary* polyhedron;
111  size_t nVertices, nFacets;
112  G4double fDz = GetZHalfLength();
113 
114  G4int subdivisions=0;
115  G4int i;
116  if(IsTwisted())
117  {
118  if ( GetVisSubdivisions()!= 0 )
119  {
120  subdivisions=GetVisSubdivisions();
121  }
122  else
123  {
124  // Estimation of Number of Subdivisions for smooth visualisation
125  //
126  G4double maxTwist=0.;
127  for(i=0; i<4; i++)
128  {
129  if(GetTwistAngle(i)>maxTwist) { maxTwist=GetTwistAngle(i); }
130  }
131 
132  // Computes bounding vectors for the shape
133  //
134  G4double Dx,Dy;
135  UVector3 minBox = GetShape()->GetMinimumBBox();
136  UVector3 maxBox = GetShape()->GetMaximumBBox();
137  G4ThreeVector minVec(minBox.x, minBox.y, minBox.z);
138  G4ThreeVector maxVec(maxBox.x, maxBox.y, maxBox.z);
139  Dx = 0.5*(maxVec.x()- minVec.y());
140  Dy = 0.5*(maxVec.y()- minVec.y());
141  if (Dy > Dx) { Dx=Dy; }
142 
143  subdivisions=8*G4int(maxTwist/(Dx*Dx*Dx)*fDz);
144  if (subdivisions<4) { subdivisions=4; }
145  if (subdivisions>30) { subdivisions=30; }
146  }
147  }
148  G4int sub4=4*subdivisions;
149  nVertices = 8+subdivisions*4;
150  nFacets = 6+subdivisions*4;
151  G4double cf=1./(subdivisions+1);
152  polyhedron = new G4PolyhedronArbitrary (nVertices, nFacets);
153 
154  // Add Vertex
155  //
156  for (i=0;i<4;i++)
157  {
158  polyhedron->AddVertex(G4ThreeVector(GetVertex(i).x(),
159  GetVertex(i).y(),-fDz));
160  }
161  for( i=0;i<subdivisions;i++)
162  {
163  for(G4int j=0;j<4;j++)
164  {
165  G4TwoVector u=GetVertex(j)+cf*(i+1)*( GetVertex(j+4)-GetVertex(j));
166  polyhedron->AddVertex(G4ThreeVector(u.x(),u.y(),-fDz+cf*2*fDz*(i+1)));
167  }
168  }
169  for (i=4;i<8;i++)
170  {
171  polyhedron->AddVertex(G4ThreeVector(GetVertex(i).x(),
172  GetVertex(i).y(),fDz));
173  }
174 
175  // Add Facets
176  //
177  polyhedron->AddFacet(1,4,3,2); //Z-plane
178  for (i=0;i<subdivisions+1;i++)
179  {
180  G4int is=i*4;
181  polyhedron->AddFacet(5+is,8+is,4+is,1+is);
182  polyhedron->AddFacet(8+is,7+is,3+is,4+is);
183  polyhedron->AddFacet(7+is,6+is,2+is,3+is);
184  polyhedron->AddFacet(6+is,5+is,1+is,2+is);
185  }
186  polyhedron->AddFacet(5+sub4,6+sub4,7+sub4,8+sub4); //Z-plane
187 
188  polyhedron->SetReferences();
189  polyhedron->InvertFacets();
190 
191  return (G4Polyhedron*) polyhedron;
192 }
G4Polyhedron * CreatePolyhedron() const
CLHEP::Hep3Vector G4ThreeVector
UGenericTrap * GetShape() const
G4double GetZHalfLength() const
G4TwoVector GetVertex(G4int index) const
G4int GetVisSubdivisions() const
G4String name
Definition: TRTMaterials.hh:40
G4UGenericTrap(const G4String &name, G4double halfZ, const std::vector< G4TwoVector > &vertices)
G4double a
Definition: TRTMaterials.hh:39
int G4int
Definition: G4Types.hh:78
double x
Definition: UVector3.hh:136
void Initialise(const std::vector< UVector2 > &vertices)
void SetName(const std::string &aName)
Definition: VUSolid.hh:104
const G4int n
void AddVertex(const G4ThreeVector &v)
UVector3 GetMaximumBBox() const
void SetZHalfLength(G4double)
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
G4USolid & operator=(const G4USolid &rhs)
Definition: G4USolid.cc:372
G4UGenericTrap & operator=(const G4UGenericTrap &source)
G4double GetTwistAngle(G4int index) const
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
double z
Definition: UVector3.hh:138
double G4double
Definition: G4Types.hh:76
G4bool IsTwisted() const
UVector3 GetMinimumBBox() const
double y
Definition: UVector3.hh:137