Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4UExtrudedSolid.cc
이 파일의 문서화 페이지로 가기
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 G4UExtrudedSolid wrapper class
31 // --------------------------------------------------------------------
32 
33 #include "G4ExtrudedSolid.hh"
34 #include "G4UExtrudedSolid.hh"
35 
36 #if ( defined(G4GEOM_USE_USOLIDS) || defined(G4GEOM_USE_PARTIAL_USOLIDS) )
37 
38 #include "G4GeomTools.hh"
39 #include "G4AffineTransform.hh"
40 #include "G4BoundingEnvelope.hh"
41 
42 #include "G4PolyhedronArbitrary.hh"
43 
45 //
46 // Constructors
47 //
48 G4UExtrudedSolid::G4UExtrudedSolid(const G4String& name,
49  std::vector<G4TwoVector> polygon,
50  std::vector<ZSection> zsections)
51  : Base_t(name) // General constructor
52 {
53  unsigned int nVertices = polygon.size();
54  unsigned int nSections = zsections.size();
55 
56  vecgeom::XtruVertex2* vertices = new vecgeom::XtruVertex2[nVertices];
57  vecgeom::XtruSection* sections = new vecgeom::XtruSection[nSections];
58 
59  for (unsigned int i = 0; i < nVertices; ++i)
60  {
61  vertices[i].x = polygon[i].x();
62  vertices[i].y = polygon[i].y();
63  }
64  for (unsigned int i = 0; i < nSections; ++i)
65  {
66  sections[i].fOrigin.Set(zsections[i].fOffset.x(),
67  zsections[i].fOffset.y(),
68  zsections[i].fZ);
69  sections[i].fScale = zsections[i].fScale;
70  }
71  Base_t::Initialize(nVertices, vertices, nSections, sections);
72  delete[] vertices;
73  delete[] sections;
74 }
75 
76 
77 G4UExtrudedSolid::G4UExtrudedSolid(const G4String& name,
78  std::vector<G4TwoVector> polygon,
79  G4double halfZ,
80  G4TwoVector off1, G4double scale1,
81  G4TwoVector off2, G4double scale2)
82  : Base_t(name) // Special constructor for 2 sections
83 {
84  unsigned int nVertices = polygon.size();
85  unsigned int nSections = 2;
86 
87  vecgeom::XtruVertex2* vertices = new vecgeom::XtruVertex2[nVertices];
88  vecgeom::XtruSection* sections = new vecgeom::XtruSection[nSections];
89 
90  for (unsigned int i = 0; i < nVertices; ++i)
91  {
92  vertices[i].x = polygon[i].x();
93  vertices[i].y = polygon[i].y();
94  }
95  sections[0].fOrigin.Set(off1.x(), off1.y(), -halfZ);
96  sections[0].fScale = scale1;
97  sections[1].fOrigin.Set(off2.x(), off2.y(), halfZ);
98  sections[1].fScale = scale2;
99  Base_t::Initialize(nVertices, vertices, nSections, sections);
100  delete[] vertices;
101  delete[] sections;
102 }
103 
105 //
106 // Fake default constructor - sets only member data and allocates memory
107 // for usage restricted to object persistency.
108 //
109 G4UExtrudedSolid::G4UExtrudedSolid(__void__& a)
110  : Base_t(a)
111 {
112 }
113 
114 
116 //
117 // Destructor
118 //
119 G4UExtrudedSolid::~G4UExtrudedSolid()
120 {
121 }
122 
123 
125 //
126 // Copy constructor
127 //
128 G4UExtrudedSolid::G4UExtrudedSolid(const G4UExtrudedSolid &source)
129  : Base_t(source)
130 {
131 }
132 
133 
135 //
136 // Assignment operator
137 //
138 G4UExtrudedSolid&
139 G4UExtrudedSolid::operator=(const G4UExtrudedSolid &source)
140 {
141  if (this == &source) return *this;
142 
143  Base_t::operator=( source );
144 
145  return *this;
146 }
147 
148 
150 //
151 // Accessors
152 
153 G4int G4UExtrudedSolid::GetNofVertices() const
154 {
155  return Base_t::GetNVertices();
156 }
157 G4TwoVector G4UExtrudedSolid::GetVertex(G4int i) const
158 {
159  G4double xx, yy;
160  Base_t::GetVertex(i, xx, yy);
161  return G4TwoVector(xx, yy);
162 }
163 std::vector<G4TwoVector> G4UExtrudedSolid::GetPolygon() const
164 {
165  std::vector<G4TwoVector> pol;
166  for (unsigned int i = 0; i < Base_t::GetNVertices(); ++i)
167  {
168  pol.push_back(GetVertex(i));
169  }
170  return pol;
171 }
172 G4int G4UExtrudedSolid::GetNofZSections() const
173 {
174  return Base_t::GetNSections();
175 }
176 G4UExtrudedSolid::ZSection G4UExtrudedSolid::GetZSection(G4int i) const
177 {
178  vecgeom::XtruSection sect = Base_t::GetSection(i);
179  return ZSection(sect.fOrigin[2],
180  G4TwoVector(sect.fOrigin[0], sect.fOrigin[1]),
181  sect.fScale);
182 }
183 std::vector<G4UExtrudedSolid::ZSection> G4UExtrudedSolid::GetZSections() const
184 {
185  std::vector<G4UExtrudedSolid::ZSection> sections;
186  for (unsigned int i = 0; i < Base_t::GetNSections(); ++i)
187  {
188  vecgeom::XtruSection sect = Base_t::GetSection(i);
189  sections.push_back(ZSection(sect.fOrigin[2],
190  G4TwoVector(sect.fOrigin[0], sect.fOrigin[1]),
191  sect.fScale));
192  }
193  return sections;
194 }
195 
196 
198 //
199 // Get bounding box
200 
201 void G4UExtrudedSolid::BoundingLimits(G4ThreeVector& pMin,
202  G4ThreeVector& pMax) const
203 {
204  static G4bool checkBBox = true;
205 
206  G4double xmin0 = kInfinity, xmax0 = -kInfinity;
207  G4double ymin0 = kInfinity, ymax0 = -kInfinity;
208 
209  for (G4int i=0; i<GetNofVertices(); ++i)
210  {
211  G4TwoVector vertex = GetVertex(i);
212  G4double x = vertex.x();
213  if (x < xmin0) xmin0 = x;
214  if (x > xmax0) xmax0 = x;
215  G4double y = vertex.y();
216  if (y < ymin0) ymin0 = y;
217  if (y > ymax0) ymax0 = y;
218  }
219 
220  G4double xmin = kInfinity, xmax = -kInfinity;
222 
223  G4int nsect = GetNofZSections();
224  for (G4int i=0; i<nsect; ++i)
225  {
226  ZSection zsect = GetZSection(i);
227  G4double dx = zsect.fOffset.x();
228  G4double dy = zsect.fOffset.y();
229  G4double scale = zsect.fScale;
230  xmin = std::min(xmin,xmin0*scale+dx);
231  xmax = std::max(xmax,xmax0*scale+dx);
232  ymin = std::min(ymin,ymin0*scale+dy);
233  ymax = std::max(ymax,ymax0*scale+dy);
234  }
235 
236  G4double zmin = GetZSection(0).fZ;
237  G4double zmax = GetZSection(nsect-1).fZ;
238 
239  pMin.set(xmin,ymin,zmin);
240  pMax.set(xmax,ymax,zmax);
241 
242  // Check correctness of the bounding box
243  //
244  if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
245  {
246  std::ostringstream message;
247  message << "Bad bounding box (min >= max) for solid: "
248  << GetName() << " !"
249  << "\npMin = " << pMin
250  << "\npMax = " << pMax;
251  G4Exception("G4UExtrudedSolid::BoundingLimits()", "GeomMgt0001",
252  JustWarning, message);
253  StreamInfo(G4cout);
254  }
255 
256  // Check consistency of bounding boxes
257  //
258  if (checkBBox)
259  {
260  UVector3 vmin, vmax;
261  Base_t::Extent(vmin,vmax);
262  if (std::abs(pMin.x()-vmin.x()) > kCarTolerance ||
263  std::abs(pMin.y()-vmin.y()) > kCarTolerance ||
264  std::abs(pMin.z()-vmin.z()) > kCarTolerance ||
265  std::abs(pMax.x()-vmax.x()) > kCarTolerance ||
266  std::abs(pMax.y()-vmax.y()) > kCarTolerance ||
267  std::abs(pMax.z()-vmax.z()) > kCarTolerance)
268  {
269  std::ostringstream message;
270  message << "Inconsistency in bounding boxes for solid: "
271  << GetName() << " !"
272  << "\nBBox min: wrapper = " << pMin << " solid = " << vmin
273  << "\nBBox max: wrapper = " << pMax << " solid = " << vmax;
274  G4Exception("G4UExtrudedSolid::BoundingLimits()", "GeomMgt0001",
275  JustWarning, message);
276  checkBBox = false;
277  }
278  }
279 }
280 
281 
283 //
284 // Calculate extent under transform and specified limit
285 
286 G4bool
287 G4UExtrudedSolid::CalculateExtent(const EAxis pAxis,
288  const G4VoxelLimits& pVoxelLimit,
289  const G4AffineTransform& pTransform,
290  G4double& pMin, G4double& pMax) const
291 {
292  G4ThreeVector bmin, bmax;
293  G4bool exist;
294 
295  // Check bounding box (bbox)
296  //
297  BoundingLimits(bmin,bmax);
298  G4BoundingEnvelope bbox(bmin,bmax);
299 #ifdef G4BBOX_EXTENT
300  if (true) return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
301 #endif
302  if (bbox.BoundingBoxVsVoxelLimits(pAxis,pVoxelLimit,pTransform,pMin,pMax))
303  {
304  return exist = (pMin < pMax) ? true : false;
305  }
306 
307  // To find the extent, the base polygon is subdivided in triangles.
308  // The extent is calculated as cumulative extent of the parts
309  // formed by extrusion of the triangles
310  //
311  G4TwoVectorList basePolygon = GetPolygon();
312  G4TwoVectorList triangles;
313  G4double eminlim = pVoxelLimit.GetMinExtent(pAxis);
314  G4double emaxlim = pVoxelLimit.GetMaxExtent(pAxis);
315 
316  // triangulate the base polygon
317  if (!G4GeomTools::TriangulatePolygon(basePolygon,triangles))
318  {
319  std::ostringstream message;
320  message << "Triangulation of the base polygon has failed for solid: "
321  << GetName() << " !"
322  << "\nExtent has been calculated using boundary box";
323  G4Exception("G4UExtrudedSolid::CalculateExtent()",
324  "GeomMgt1002",JustWarning,message);
325  return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
326  }
327 
328  // allocate vector lists
329  G4int nsect = GetNofZSections();
330  std::vector<const G4ThreeVectorList *> polygons;
331  polygons.resize(nsect);
332  for (G4int k=0; k<nsect; ++k) { polygons[k] = new G4ThreeVectorList(3); }
333 
334  // main loop along triangles
335  pMin = kInfinity;
336  pMax = -kInfinity;
337  G4int ntria = triangles.size()/3;
338  for (G4int i=0; i<ntria; ++i)
339  {
340  G4int i3 = i*3;
341  for (G4int k=0; k<nsect; ++k) // extrude triangle
342  {
343  ZSection zsect = GetZSection(k);
344  G4double z = zsect.fZ;
345  G4double dx = zsect.fOffset.x();
346  G4double dy = zsect.fOffset.y();
347  G4double scale = zsect.fScale;
348 
349  G4ThreeVectorList* ptr = const_cast<G4ThreeVectorList*>(polygons[k]);
350  G4ThreeVectorList::iterator iter = ptr->begin();
351  G4double x0 = triangles[i3+0].x()*scale+dx;
352  G4double y0 = triangles[i3+0].y()*scale+dy;
353  iter->set(x0,y0,z);
354  iter++;
355  G4double x1 = triangles[i3+1].x()*scale+dx;
356  G4double y1 = triangles[i3+1].y()*scale+dy;
357  iter->set(x1,y1,z);
358  iter++;
359  G4double x2 = triangles[i3+2].x()*scale+dx;
360  G4double y2 = triangles[i3+2].y()*scale+dy;
361  iter->set(x2,y2,z);
362  }
363 
364  // set sub-envelope and adjust extent
365  G4double emin,emax;
366  G4BoundingEnvelope benv(polygons);
367  if (!benv.CalculateExtent(pAxis,pVoxelLimit,pTransform,emin,emax)) continue;
368  if (emin < pMin) pMin = emin;
369  if (emax > pMax) pMax = emax;
370  if (eminlim > pMin && emaxlim < pMax) break; // max possible extent
371  }
372  // free memory
373  for (G4int k=0; k<nsect; ++k) { delete polygons[k]; polygons[k]=0;}
374  return (pMin < pMax);
375 }
376 
377 
379 //
380 // CreatePolyhedron()
381 //
382 G4Polyhedron* G4UExtrudedSolid::CreatePolyhedron () const
383 {
384  unsigned int nFacets = Base_t::GetStruct().fTslHelper.fFacets.size();
385  unsigned int fVertices = 3*nFacets;
386 
387  G4PolyhedronArbitrary *polyhedron =
388  new G4PolyhedronArbitrary (fVertices,nFacets);
389 
390  for (unsigned int i = 0; i < nFacets; ++i)
391  {
392  G4int v[3]; // Facets are only triangular in VecGeom
393  for (unsigned int j = 0; j < 3; ++j)
394  {
395  UVector3 vtx = Base_t::GetStruct().fTslHelper.fFacets[i]->fVertices[j];
396  polyhedron->AddVertex(G4ThreeVector(vtx.x(), vtx.y(), vtx.z()));
397  v[j] = Base_t::GetStruct().fTslHelper.fFacets[i]->fIndices[j]+1;
398  }
399  polyhedron->AddFacet(v[0],v[1],v[2],0);
400  }
401  polyhedron->SetReferences();
402 
403  return (G4Polyhedron*) polyhedron;
404 }
405 
406 #endif // G4GEOM_USE_USOLIDS
Float_t x
Definition: compare.C:6
void set(double x, double y, double z)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const XML_Char * name
Definition: expat.h:151
static G4bool TriangulatePolygon(const G4TwoVectorList &polygon, G4TwoVectorList &result)
Definition: G4GeomTools.cc:199
CLHEP::Hep3Vector G4ThreeVector
Double_t xx
G4double GetMaxExtent(const EAxis pAxis) const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
static const G4double kInfinity
Definition: geomdefs.hh:42
Float_t y1[n_points_granero]
Definition: compare.C:5
void AddFacet(const G4int iv1, const G4int iv2, const G4int iv3, const G4int iv4=0)
void AddVertex(const G4ThreeVector &v)
Float_t x1[n_points_granero]
Definition: compare.C:5
Float_t y
Definition: compare.C:6
void Initialize()
Definition: errprop.cc:101
CLHEP::Hep2Vector G4TwoVector
Definition: G4TwoVector.hh:42
Double_t z
void message(RunManager *runmanager)
Definition: ts_scorers.cc:72
double z() const
static const G4double emax
const G4double kCarTolerance
Float_t y2[n_points_geant4]
Definition: compare.C:26
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
Double_t scale
std::vector< G4ThreeVector > G4ThreeVectorList
G4double GetMinExtent(const EAxis pAxis) const
double x() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
int G4int
Definition: G4Types.hh:78
std::vector< G4TwoVector > G4TwoVectorList
Definition: G4GeomTools.hh:50
EAxis
Definition: geomdefs.hh:54
double y() const
G4GLOB_DLL std::ostream G4cout
double x() const
double y() const
Float_t x2[n_points_geant4]
Definition: compare.C:26
T min(const T t1, const T t2)
brief Return the smallest of the two arguments