g4tools  5.4.0
polyhedron
Go to the documentation of this file.
1 // Copyright (C) 2010, Guy Barrand. All rights reserved.
2 // See the file tools.license for terms.
3 
4 #ifndef tools_hep_polyhedron
5 #define tools_hep_polyhedron
6 
7 // see (lengthy) doc and disclaimer at end.
8 
9 #include "../lina/vec3d"
10 #include "../lina/rotd"
11 
12 #ifdef TOOLS_MEM
13 #include "../mem"
14 #include "../S_STRING"
15 #endif
16 
17 #include <string>
18 #include <vector>
19 
20 //#define TOOLS_HEP_PH_OUT_ERR
21 //#define TOOLS_HEP_PH_OUT_ERR_TRD2
22 //#define TOOLS_HEP_PH_NOT_OPT
23 
24 #ifdef TOOLS_HEP_PH_OUT_ERR
25 #include <iostream>
26 #endif
27 
28 #ifdef TOOLS_HEP_PH_OUT_ERR_TRD2
29 #include <iostream>
30 #endif
31 
32 #include <cmath>
33 
34 namespace tools {
35 namespace hep {
36 
37 typedef vec3d HVPoint3D;
40 
41 //WARNING : with SbFacet, take care of exlib::geant4::polyhedron
42 // that attempts to copy the private content of G4Facet.
43 // (see WARNING here).
44 
45 class SbFacet {
46 #ifdef TOOLS_MEM
48 #endif
49 
50  friend class polyhedron;
51 #ifndef SWIG
52  friend std::ostream& operator<<(std::ostream&, const SbFacet &facet);
53  //G.Barrand
54  friend int operator == (const SbFacet & v1, const SbFacet & v2);
55  friend int operator != (const SbFacet & v1, const SbFacet & v2);
56 #endif
57 
58  private:
59  typedef struct { int v,f; } edge_t; //G.Barrand
60  edge_t edge[4];
61 
62  public:
63  SbFacet(int v1=0, int f1=0, int v2=0, int f2=0,
64  int v3=0, int f3=0, int v4=0, int f4=0)
65  {
66 #ifdef TOOLS_MEM
67  mem::increment(s_class().c_str());
68 #endif
69  edge[0].v=v1; edge[0].f=f1; edge[1].v=v2; edge[1].f=f2;
70  edge[2].v=v3; edge[2].f=f3; edge[3].v=v4; edge[3].f=f4; }
71 
72  virtual ~SbFacet() {
73 #ifdef TOOLS_MEM
74  mem::decrement(s_class().c_str());
75 #endif
76  }
77 //public:
78 protected:
79  SbFacet(const SbFacet & aFrom) {
80 #ifdef TOOLS_MEM
81  mem::increment(s_class().c_str());
82 #endif
83  edge[0].v = aFrom.edge[0].v;
84  edge[0].f = aFrom.edge[0].f;
85  edge[1].v = aFrom.edge[1].v;
86  edge[1].f = aFrom.edge[1].f;
87  edge[2].v = aFrom.edge[2].v;
88  edge[2].f = aFrom.edge[2].f;
89  edge[3].v = aFrom.edge[3].v;
90  edge[3].f = aFrom.edge[3].f;
91  }
92  SbFacet& operator=(const SbFacet& aFrom) {
93  edge[0].v = aFrom.edge[0].v;
94  edge[0].f = aFrom.edge[0].f;
95  edge[1].v = aFrom.edge[1].v;
96  edge[1].f = aFrom.edge[1].f;
97  edge[2].v = aFrom.edge[2].v;
98  edge[2].f = aFrom.edge[2].f;
99  edge[3].v = aFrom.edge[3].v;
100  edge[3].f = aFrom.edge[3].f;
101  return *this;
102  }
103 public:
104  bool isEqual(const SbFacet& aFrom) const { //G.Barrand
105  if(edge[0].v!=aFrom.edge[0].v) return false;
106  if(edge[0].f!=aFrom.edge[0].f) return false;
107 
108  if(edge[1].v!=aFrom.edge[1].v) return false;
109  if(edge[1].f!=aFrom.edge[1].f) return false;
110 
111  if(edge[2].v!=aFrom.edge[2].v) return false;
112  if(edge[2].f!=aFrom.edge[2].f) return false;
113 
114  if(edge[3].v!=aFrom.edge[3].v) return false;
115  if(edge[3].f!=aFrom.edge[3].f) return false;
116 
117  return true;
118  }
119  void GetEdge(int i,int& v,int& f) const { //G.Barrand
120  v = edge[i].v;
121  f = edge[i].f;
122  }
123 
124  void set(int v1, int f1, int v2, int f2, //G.Barrand
125  int v3, int f3, int v4, int f4) {
126  edge[0].v=v1; edge[0].f=f1; edge[1].v=v2; edge[1].f=f2;
127  edge[2].v=v3; edge[2].f=f3; edge[3].v=v4; edge[3].f=f4;
128  }
129 
130 
131  void Set(int v[8]) //G.Barrand
132  { edge[0].v = v[0]; edge[0].f = v[1];
133  edge[1].v = v[2]; edge[1].f = v[3];
134  edge[2].v = v[4]; edge[2].f = v[5];
135  edge[3].v = v[6]; edge[3].f = v[7]; }
136 };
137 
138 //G.Barrand :
139 //int operator == (const SbFacet & v1, const SbFacet & v2);
140 //int operator != (const SbFacet & v1, const SbFacet & v2);
141 
142 class polyhedron {
143 #ifdef TOOLS_MEM
145 #endif
146 
147 #ifndef SWIG
148  friend std::ostream& operator<<(std::ostream&, const polyhedron &ph);
149  //G.Barrand
150  friend int operator == (const polyhedron & v1, const polyhedron & v2);
151  friend int operator != (const polyhedron & v1, const polyhedron & v2);
152 #endif
153 
154 private: //G.Barrand
155  //std::string* m_name; //have a pointer to optimize memory.
156 protected:
157  int nvert, nface;
160 private:
161  int fNumberOfRotationSteps;
162 
163 protected:
164  static double _M_PI() {return 3.1415926535897931160E0;}
165  //static double _M_PI_2() {return 1.5707963267948965580E0;}
166 
167  // Allocate memory for polyhedron
168  void AllocateMemory(int Nvert, int Nface);
169 
170  // Find neighbouring facet
171  int FindNeighbour(int iFace, int iNode, int iOrder) const;
172 
173  // Find normal at node
174  HVNormal3D FindNodeNormal(int iFace, int iNode) const;
175 
176  // Create polyhedron for prism with quadrilateral base
177  void CreatePrism();
178 
179  // Generate facets by revolving an edge around Z-axis
180  void RotateEdge(int k1, int k2, double r1, double r2,
181  int v1, int v2, int vEdge,
182  bool ifWholeCircle, int ns, int &kface);
183 
184  // Set side facets for the case of incomplete rotation
185  void SetSideFacets(int ii[4], int vv[4],
186  int *kk, double *r,
187  double dphi, int ns, int &kface);
188 
189  // Create polyhedron for body of revolution around Z-axis
190  void RotateAroundZ(int nstep, double phi, double dphi,
191  int np1, int np2,
192  const double *z, double *r,
193  int nodeVis, int edgeVis);
194 
195  // For each edge set reference to neighbouring facet
197 
198  // Invert the order on nodes in facets
199  void InvertFacets();
200 public: //public for iv2sg
201  static int NUMBER_OF_STEPS() {return 24;}
202 public: //for iv2sg
203  static void do_not_set_NUMBER_OF_STEPS(int){}
204 public:
205  polyhedron(int Nvert=0, int Nface=0)
206  : /*m_name(0) //G.Barrand
207  ,*/nvert(Nvert),nface(Nface)
208  ,pV(Nvert ? new HVPoint3D[Nvert+1] : 0)
209  ,pF(Nface ? new SbFacet[Nface+1] : 0)
210  ,fNumberOfRotationSteps(NUMBER_OF_STEPS())
211  {
212 #ifdef TOOLS_MEM
213  mem::increment(s_class().c_str());
214 #endif
215  }
216 public:
217  virtual ~polyhedron() {
218  //delete m_name; //G.Barrand.
219  delete [] pV; delete [] pF;
220 #ifdef TOOLS_MEM
221  mem::decrement(s_class().c_str());
222 #endif
223  }
224 public:
225  polyhedron(const polyhedron & from);
227 public:
228  //G.Barrand : handle a name to help debugging.
229 /*
230  void setName(const std::string& aName) {
231  delete m_name;
232  m_name = new std::string(aName);
233  }
234  const std::string& getName() const {
235  if(!m_name) return s_empty();
236  return *m_name;
237  }
238 */
239  //G.Barrand :end
240 
241  void Set(int Nvert, HVPoint3D* aV,
242  int Nface, SbFacet* aF) //G.Barrand
243  { delete [] pV; delete [] pF;
244  nvert = Nvert; nface = Nface; pV = aV; pF = aF;}
245 
246  void Empty() //G.Barrand
247  { nvert = 0; nface = 0; pV = 0;pF = 0;}
248 
249  // Get number of vertices
250  int GetNoVertices() const { return nvert; }
251 
252  // Get number of facets
253  int GetNoFacets() const { return nface; }
254 
255  // Transform the polyhedron
256  polyhedron& Translate(double,double,double);
257  polyhedron& Transform(const rotd& rot,double,double,double);
258  polyhedron& Transform(const rotd& rot,const vec3d& trans);
259 
260  // Get next vertex index of the quadrilateral
261  //G.Barrand
262  bool GetNextVertexIndex(int & index, int & edgeFlag) const;
263 
264  // Get vertex by index
265  HVPoint3D GetVertex(int index) const;
266  const HVPoint3D& GetVertexFast(int index) const; //G.Barrand
267 
268  //G.Barrand : to optimize SoPolyhedron.
269  HVPoint3D* GetPV() const {return pV;} //G.Barrand
270  SbFacet* GetPF() const {return pF;} //G.Barrand
271 
272  // Get next vertex + edge visibility of the quadrilateral
273  bool GetNextVertex(HVPoint3D & vertex, int & edgeFlag) const;
274 
275  // Get next vertex + edge visibility + normal of the quadrilateral
276  bool GetNextVertex(HVPoint3D & vertex, int & edgeFlag,
277  HVNormal3D & normal) const;
278 
279  // Get indeces of the next edge with indeces of the faces
280  bool GetNextEdgeIndeces(int & i1, int & i2, int & edgeFlag,
281  int & iface1, int & iface2) const;
282 
283  // Get indeces of the next edge
284  bool GetNextEdgeIndeces(int & i1, int & i2, int & edgeFlag) const;
285 
286  // Get next edge
287  bool GetNextEdge(HVPoint3D &p1, HVPoint3D &p2, int &edgeFlag) const;
288 
289  // Get next edge
290  bool GetNextEdge(HVPoint3D &p1, HVPoint3D &p2, int &edgeFlag,
291  int &iface1, int &iface2) const;
292 
293  // Get face by index
294  void GetFacet(int iFace, int &n, int *iNodes,
295  int *edgeFlags = 0, int *iFaces = 0) const;
296 
297  // Get face by index
298  void GetFacet(int iFace, int &n, HVPoint3D *nodes,
299  int *edgeFlags = 0, HVNormal3D *normals = 0) const;
300 
301  // Get next face with normals at the nodes
302  bool GetNextFacet(int &n, HVPoint3D *nodes,
303  int *edgeFlags=0, HVNormal3D *normals=0) const;
304 
305  // Get normal of the face given by index
306  HVNormal3D GetNormal(int iFace) const;
307 
308  // Get unit normal of the face given by index
309  HVNormal3D GetUnitNormal(int iFace) const;
310 
311  // Get normal of the next face
312  bool GetNextNormal(HVNormal3D &normal) const;
313 
314  // Get normal of unit length of the next face
315  bool GetNextUnitNormal(HVNormal3D &normal) const;
316 
317  // Boolean operations
318  polyhedron add(const polyhedron &p) const;
319  polyhedron subtract(const polyhedron &p) const;
320  polyhedron intersect(const polyhedron &p) const;
321 
322  // Get area of the surface of the polyhedron
323  double GetSurfaceArea() const;
324 
325  // Get volume of the polyhedron
326  double GetVolume() const;
327 
328  bool isEqual(const polyhedron &p) const; //G.Barrand
329  bool isConsistent(const char* = 0) const; //G.Barrand
330  void dump(std::ostream&) const;
331 
332  // Get number of steps for whole circle
333  int GetNumberOfRotationSteps(); //G.Barrand : no more static.
334 
335  // Set number of steps for whole circle
337 
338  // Reset number of steps for whole circle to default value
339  void ResetNumberOfRotationSteps(); //G.Barrand : have code in .cxx.
340 
341 public:
342  //G.Barrand : have the below set_ to optimize exlib/sg/polyhedron setup.
343  bool set_polyhedron_cons(double Rmn1, double Rmx1,
344  double Rmn2, double Rmx2, double Dz,
345  double Phi1, double Dphi,
346  int nstep = 0); //G.Barrand
347  bool set_polyhedron_tube(double Rmin, double Rmax, double Dz,
348  int nstep = 0){
349  return set_polyhedron_cons(Rmin, Rmax, Rmin, Rmax, Dz, 0, 2*_M_PI(), nstep);
350  }
351 
352  double vxy(const double* xy,int i,int j) {return xy[i*2+j];}
353  bool set_polyhedron_arb8(double Dz,const double* xy);
354 
355  bool set_polyhedron_trd2(double Dx1, double Dx2,
356  double Dy1, double Dy2, double Dz);
357  bool set_polyhedron_box(double Dx, double Dy, double Dz){
358  return set_polyhedron_trd2(Dx, Dx, Dy, Dy, Dz);
359  }
360  bool set_polyhedron_trd1(double Dx1, double Dx2,
361  double Dy, double Dz){
362  return set_polyhedron_trd2(Dx1, Dx2, Dy, Dy, Dz);
363  }
364  bool set_polyhedron_trap(double Dz, double Theta, double Phi,
365  double Dy1,
366  double Dx1, double Dx2, double Alp1,
367  double Dy2,
368  double Dx3, double Dx4, double Alp2);
369 
370  bool set_polyhedron_para(double Dx, double Dy, double Dz,
371  double Alpha, double Theta,
372  double Phi){
373  return set_polyhedron_trap(Dz,Theta,Phi,Dy,Dx,Dx,Alpha,Dy,Dx,Dx,Alpha);
374  }
375 
376  bool set_polyhedron_pgon(double phi, double dphi, int npdv, int nz,
377  const double *z,
378  const double *rmin,
379  const double *rmax);
380  bool set_polyhedron_pcon(double phi, double dphi, int nz,
381  const double *z,
382  const double *rmin,
383  const double *rmax) {
384  return set_polyhedron_pgon(phi, dphi, 0, nz, z, rmin, rmax);
385  }
386 
387  bool set_polyhedron_tubs(double Rmin, double Rmax,
388  double Dz,
389  double Phi1, double Dphi,
390  int nstep) {//G.Barrand
391  return set_polyhedron_cons(Rmin, Rmax, Rmin, Rmax, Dz, Phi1, Dphi, nstep);
392  }
393 
394  bool set_polyhedron_cone(double Rmn1, double Rmx1,
395  double Rmn2, double Rmx2,
396  double Dz,
397  int nstep) {
398  return set_polyhedron_cons(Rmn1, Rmx1, Rmn2, Rmx2, Dz, 0, 2*_M_PI(), nstep);
399  }
400 
401  bool set_polyhedron_torus(double rmin,double rmax,double rtor,
402  double phi,double dphi,
403  int nphi, //G.Barrand
404  int nthe); //G.Barrand
405 
406  bool set_polyhedron_xtru(int a_npts,int a_nz,
407  double* a_xs,double* a_ys,double* a_zs,
408  bool a_acw = true,
409  bool a_zfb = true);
410 
411  bool set_polyhedron_sphere(double rmin, double rmax,
412  double phi, double dphi,
413  double the, double dthe,
414  int nphi = 0,
415  int nthe = 0); //G.Barrand
416 
417  bool set_polyhedron_hype(double a_st_in,double a_st_out,
418  double a_rmin,double a_rmax,double a_dz,
419  int a_nz = 10,int a_nphi = 24);
420  bool set_polyhedron_eltu(double a_dx,double a_dy,double a_dz,
421  int a_nz = 10,int a_nphi = 24);
422 private: //G.Barrand
423  int _ixy(int,int,int,int,bool,bool);
424  void _clear(){
425  //used in set_polyhedronXxx()
426  delete [] pV;
427  pV = 0;
428  delete [] pF;
429  pF = 0;
430  nvert = 0;
431  nface = 0;
432  }
433  bool CHECK_INDEX(const char* a_method,int a_index) const;
434 };
435 
436 //G.Barrand :
437 //int operator == (const polyhedron & v1, const polyhedron & v2);
438 //int operator != (const polyhedron & v1, const polyhedron & v2);
439 
440 // G.Barrand : introduce iabs to avoid a mess with cmath and some compiler.
441 inline int Sb_iabs(int a) {
442  return a < 0 ? -a : a;
443 }
444 
445 inline //G.Barrand
446 bool polyhedron::GetNextVertexIndex(int &index, int &edgeFlag) const
447 /***********************************************************************
448  * *
449  * Name: polyhedron::GetNextVertexIndex Date: 03.09.96 *
450  * Author: Yasuhide Sawada Revised: *
451  * *
452  * Function: *
453  * *
454  ***********************************************************************/
455 {
456  static int iFace = 1;
457  static int iQVertex = 0;
458  //G.Barrand : int vIndex = pF[iFace].edge[iQVertex].v;
459  SbFacet::edge_t* edge = pF[iFace].edge; //G.Barrand : optimize.
460  int vIndex = edge[iQVertex].v;
461 
462  edgeFlag = (vIndex > 0) ? 1 : 0;
463  index = Sb_iabs(vIndex);
464 
465  if(index>nvert) {
466 #ifdef TOOLS_HEP_PH_OUT_ERR
467  std::cerr << "polyhedron::GetNextVertexIndex: pV index problem "
468  << index << " exceed " << nvert << std::endl;
469 #endif
470  index = 0;
471  }
472 
473  //G.Barrand : if (iQVertex >= 3 || pF[iFace].edge[iQVertex+1].v == 0) {
474  if (iQVertex >= 3 || edge[iQVertex+1].v == 0) {
475  iQVertex = 0;
476  if (++iFace > nface) iFace = 1;
477  return false; // Last Edge
478  }else{
479  ++iQVertex;
480  return true; // not Last Edge
481  }
482 }
483 
484 class polyhedron_trd2 : public polyhedron {
485  public:
486  polyhedron_trd2(double Dx1, double Dx2,
487  double Dy1, double Dy2, double Dz);
488  virtual ~polyhedron_trd2(){}
489  public:
490  polyhedron_trd2(const polyhedron_trd2& a_from):polyhedron(a_from){}
492  polyhedron::operator=(a_from);
493  return *this;
494  }
495 
496  //virtual polyhedron& operator = (const polyhedron& from) {
497  // return polyhedron::operator = (from);
498  //}
499 };
500 
501 class polyhedron_arb8 : public polyhedron {
502  public:
503  polyhedron_arb8(double Dz,const double* xy);
504  virtual ~polyhedron_arb8(){}
505  public:
506  polyhedron_arb8(const polyhedron_arb8& a_from):polyhedron(a_from){}
508  polyhedron::operator=(a_from);
509  return *this;
510  }
511  //virtual polyhedron& operator = (const polyhedron& from) {
512  // return polyhedron::operator = (from);
513  //}
514 };
515 
516 class polyhedron_xtru : public polyhedron {
517  public:
518  polyhedron_xtru(int a_npts,int a_nz,
519  double* a_xs,double* a_ys,double* a_zs,
520  bool a_acw = true,
521  bool a_zfb = true);
522  virtual ~polyhedron_xtru(){}
523  public:
524  polyhedron_xtru(const polyhedron_xtru& a_from):polyhedron(a_from){}
526  polyhedron::operator=(a_from);
527  return *this;
528  }
529 };
530 
531 class polyhedron_hype : public polyhedron {
532  public:
533  polyhedron_hype(double a_st_in,double a_st_out,
534  double a_rmin,double a_rmax,double a_dz,
535  int a_nz = 10,int a_nphi = 24);
536  virtual ~polyhedron_hype(){}
537  public:
538  polyhedron_hype(const polyhedron_hype& a_from):polyhedron(a_from){}
540  polyhedron::operator=(a_from);
541  return *this;
542  }
543 };
544 
546  public:
547  polyhedron_trd1(double Dx1, double Dx2,
548  double Dy, double Dz);
549  virtual ~polyhedron_trd1(){}
550  public:
554  return *this;
555  }
556  //virtual polyhedron& operator = (const polyhedron& from) {
557  // return polyhedron::operator = (from);
558  //}
559 };
560 
562  public:
563  polyhedron_box(double Dx, double Dy, double Dz);
564  virtual ~polyhedron_box(){}
565  public:
569  return *this;
570  }
571  //virtual polyhedron& operator = (const polyhedron& from) {
572  // return polyhedron::operator = (from);
573  //}
574 };
575 
576 class polyhedron_trap : public polyhedron {
577 public:
578  polyhedron_trap(double Dz, double Theta, double Phi,
579  double Dy1,
580  double Dx1, double Dx2, double Alp1,
581  double Dy2,
582  double Dx3, double Dx4, double Alp2);
583  virtual ~polyhedron_trap(){}
584  public:
585  polyhedron_trap(const polyhedron_trap& a_from):polyhedron(a_from){}
587  polyhedron::operator=(a_from);
588  return *this;
589  }
590  //virtual polyhedron& operator = (const polyhedron& from) {
591  // return polyhedron::operator = (from);
592  //}
593 };
594 
596 public:
597  polyhedron_para(double Dx, double Dy, double Dz,
598  double Alpha, double Theta, double Phi);
599  virtual ~polyhedron_para(){}
600  public:
604  return *this;
605  }
606  //virtual polyhedron& operator = (const polyhedron& from) {
607  // return polyhedron::operator = (from);
608  //}
609 };
610 
611 class polyhedron_cons : public polyhedron {
612 public:
613  polyhedron_cons(double Rmn1, double Rmx1,
614  double Rmn2, double Rmx2, double Dz,
615  double Phi1, double Dphi,
616  int nstep = 0); //G.Barrand
617  virtual ~polyhedron_cons(){}
618  public:
619  polyhedron_cons(const polyhedron_cons& a_from):polyhedron(a_from){}
621  polyhedron::operator=(a_from);
622  return *this;
623  }
624  //virtual polyhedron& operator = (const polyhedron& from) {
625  // return polyhedron::operator = (from);
626  //}
627 };
628 
630 public:
631  polyhedron_cone(double Rmn1, double Rmx1,
632  double Rmn2, double Rmx2, double Dz,
633  int nstep = 0); //G.Barrand
634  virtual ~polyhedron_cone(){}
635  public:
639  return *this;
640  }
641  //virtual polyhedron& operator = (const polyhedron& from) {
642  // return polyhedron::operator = (from);
643  //}
644 };
645 
647 public:
648  polyhedron_tubs(double Rmin, double Rmax, double Dz,
649  double Phi1, double Dphi,
650  int nstep = 0); //G.Barrand
651  virtual ~polyhedron_tubs(){}
652  public:
656  return *this;
657  }
658  //virtual polyhedron& operator = (const polyhedron& from) {
659  // return polyhedron::operator = (from);
660  //}
661 };
662 
664 public:
665  polyhedron_tube(double Rmin, double Rmax, double Dz,int nstep = 0); //G.Barrand
666  virtual ~polyhedron_tube(){}
667  public:
671  return *this;
672  }
673  //virtual polyhedron& operator = (const polyhedron& from) {
674  // return polyhedron::operator = (from);
675  //}
676 };
677 
678 class polyhedron_pgon : public polyhedron {
679 public:
680  polyhedron_pgon(double phi, double dphi, int npdv, int nz,
681  const double *z,
682  const double *rmin,
683  const double *rmax);
684  virtual ~polyhedron_pgon(){}
685  public:
686  polyhedron_pgon(const polyhedron_pgon& a_from):polyhedron(a_from){}
688  polyhedron::operator=(a_from);
689  return *this;
690  }
691  //virtual polyhedron& operator = (const polyhedron& from) {
692  // return polyhedron::operator = (from);
693  //}
694 };
695 
697 public:
698  polyhedron_pcon(double phi, double dphi, int nz,
699  const double *z,
700  const double *rmin,
701  const double *rmax);
702  virtual ~polyhedron_pcon(){}
703  public:
707  return *this;
708  }
709  //virtual polyhedron& operator = (const polyhedron& from) {
710  // return polyhedron::operator = (from);
711  //}
712 };
713 
715 public:
716  polyhedron_sphere(double rmin, double rmax,
717  double phi, double dphi,
718  double the, double dthe,
719  int nphi = 0,
720  int nthe = 0); //G.Barrand
721  virtual ~polyhedron_sphere(){}
722  public:
725  polyhedron::operator=(a_from);
726  return *this;
727  }
728  //virtual polyhedron& operator = (const polyhedron& from) {
729  // return polyhedron::operator = (from);
730  //}
731 };
732 
733 class polyhedron_torus : public polyhedron {
734 public:
735  polyhedron_torus(double rmin, double rmax, double rtor,
736  double phi, double dphi,
737  int nphi = 0,
738  int nthe = 0); //G.Barrand
739  virtual ~polyhedron_torus(){}
740  public:
743  polyhedron::operator=(a_from);
744  return *this;
745  }
746  //virtual polyhedron& operator = (const polyhedron& from) {
747  // return polyhedron::operator = (from);
748  //}
749 };
750 
751 //G.Barrand : begin
753 #ifdef TOOLS_MEM
755 #endif
756 public:
757  enum Operation { //Must be the same than BooleanProcessor OP_XXX.
758  UNION = 0
761  };
762 private:
763  typedef std::pair<Operation,polyhedron> op_t;
764 public:
766 #ifdef TOOLS_MEM
767  mem::increment(s_class().c_str());
768 #endif
769  }
771 #ifdef TOOLS_MEM
772  mem::decrement(s_class().c_str());
773 #endif
774  }
775 private:
777 #ifdef TOOLS_MEM
778  mem::increment(s_class().c_str());
779 #endif
780  }
781  polyhedronProcessor& operator=(const polyhedronProcessor&){return *this;}
782 public:
783  void push_back(Operation a_op,const polyhedron& a_polyhedron) {
784  m_ops.push_back(op_t(a_op,a_polyhedron));
785  }
787  void clear() { m_ops.clear();}
788  bool is_same_op() const {
789  if(!m_ops.size()) return true;
790  Operation op = m_ops[0].first;
791  std::vector<op_t>::const_iterator it;
792  for(it=m_ops.begin();it!=m_ops.end();++it) {
793  if((*it).first!=op) return false;
794  }
795  return true;
796  }
797 
798 //private:
799  bool execute1(polyhedron&,const std::vector<unsigned int>&);
800 private:
801  std::vector<op_t> m_ops;
802 };
803 //G.Barrand : end
804 
805 //inline const std::string& stype(const polyhedron&) {
806 // static const std::string s_v("tools::hep::polyhedron");
807 // return s_v;
808 //}
809 
810 }}
811 
812 #include "polyhedron.icc"
813 
814 namespace tools {
815 namespace hep {
816 
817 template <class MATRIX>
818 inline void tsf_polyhedron(polyhedron& a_ph,const MATRIX& a_matrix) {
819  typedef typename MATRIX::elem_t T;
820  int nvert = a_ph.GetNoVertices();
821  hep::HVPoint3D* pV = a_ph.GetPV();
822  if (nvert > 0) {
823  T x,y,z;
824  for (int i=1; i<=nvert; i++) {
825  hep::HVPoint3D& p = pV[i];
826  x = T(p.x());
827  y = T(p.y());
828  z = T(p.z());
829  a_matrix.mul_3(x,y,z);
830  p.set_value(x,y,z);
831  }
832  }
833 }
834 
835 }}
836 
837 #endif
838 
839 //--------------------------------------------------------------------//
840 // JFB: //
841 // polyhedron was HepPolyhedron, retrofitted to Open Inventor //
842 // infrastructure: //
843 //--------------------------------------------------------------------//
844 
845 
846 // ********************************************************************
847 // * DISCLAIMER *
848 // * *
849 // * The following disclaimer summarizes all the specific disclaimers *
850 // * of contributors to this software. The specific disclaimers,which *
851 // * govern, are listed with their locations in: *
852 // * http://cern.ch/geant4/license *
853 // * *
854 // * Neither the authors of this software system, nor their employing *
855 // * institutes,nor the agencies providing financial support for this *
856 // * work make any representation or warranty, express or implied, *
857 // * regarding this software system or assume any liability for its *
858 // * use. *
859 // * *
860 // * This code implementation is the intellectual property of the *
861 // * GEANT4 collaboration. *
862 // * By copying, distributing or modifying the Program (or any work *
863 // * based on the Program) you indicate your acceptance of this *
864 // * statement, and all its terms. *
865 // ********************************************************************
866 //
867 //
868 //
869 //
870 // Class Description:
871 // polyhedron is an intermediate class between description of a shape
872 // and visualization systems. It is intended to provide some service like:
873 // - polygonization of shapes with triangulization (quadrilaterization)
874 // of complex polygons;
875 // - calculation of normals for faces and vertices;
876 // - finding result of boolean operation on polyhedra;
877 //
878 // Public constructors:
879 //
880 // polyhedron_box (dx,dy,dz)
881 // - create polyhedron for Box;
882 // polyhedron_trd1 (dx1,dx2,dy,dz)
883 // - create polyhedron for G3 Trd1;
884 // polyhedron_trd2 (dx1,dx2,dy1,dy2,dz)
885 // - create polyhedron for G3 Trd2;
886 // polyhedron_trap (dz,theta,phi, h1,bl1,tl1,alp1, h2,bl2,tl2,alp2)
887 // - create polyhedron for G3 Trap;
888 // polyhedron_para (dx,dy,dz,alpha,theta,phi)
889 // - create polyhedron for G3 Para;
890 // polyhedron_tube (rmin,rmax,dz,nstep=0)
891 // - create polyhedron for G3 Tube;
892 // polyhedron_tubs (rmin,rmax,dz,phi1,dphi,nstep=0)
893 // - create polyhedron for G3 Tubs;
894 // polyhedron_cone (rmin1,rmax1,rmin2,rmax2,dz,nstep=0)
895 // - create polyhedron for G3 Cone;
896 // polyhedron_cons (rmin1,rmax1,rmin2,rmax2,dz,phi1,dphi,nstep=0)
897 // - create polyhedron for G3 Cons;
898 // polyhedron_pgon (phi,dphi,npdv,nz, z(*),rmin(*),rmax(*))
899 // - create polyhedron for G3 Pgon;
900 // polyhedron_pcon (phi,dphi,nz, z(*),rmin(*),rmax(*))
901 // - create polyhedron for G3 Pcon;
902 // polyhedron_sphere (rmin,rmax,phi,dphi,the,dthe,nstep=0)
903 // - create polyhedron for Sphere;
904 // polyhedron_torus (rmin,rmax,rtor,phi,dphi,nstep=0)
905 // - create polyhedron for Torus;
906 // Public functions:
907 //
908 // GetNoVertices () - returns number of vertices;
909 // GetNoFacets () - returns number of faces;
910 // GetNextVertexIndex (index,edgeFlag) - get vertex indeces of the
911 // quadrilaterals in order;
912 // returns false when finished each face;
913 // GetVertex (index) - returns vertex by index;
914 // GetNextVertex (vertex,edgeFlag) - get vertices with edge visibility
915 // of the quadrilaterals in order;
916 // returns false when finished each face;
917 // GetNextVertex (vertex,edgeFlag,normal) - get vertices with edge
918 // visibility and normal of the quadrilaterals
919 // in order; returns false when finished each face;
920 // GetNextEdgeIndeces (i1,i2,edgeFlag) - get indeces of the next edge;
921 // returns false for the last edge;
922 // GetNextEdgeIndeces (i1,i2,edgeFlag,iface1,iface2) - get indeces of
923 // the next edge with indeces of the faces
924 // to which the edge belongs;
925 // returns false for the last edge;
926 // GetNextEdge (p1,p2,edgeFlag) - get next edge;
927 // returns false for the last edge;
928 // GetNextEdge (p1,p2,edgeFlag,iface1,iface2) - get next edge with indeces
929 // of the faces to which the edge belongs;
930 // returns false for the last edge;
931 // GetFacet (index,n,nodes,edgeFlags=0,normals=0) - get face by index;
932 // GetNextFacet (n,nodes,edgeFlags=0,normals=0) - get next face with normals
933 // at the nodes; returns false for the last face;
934 // GetNormal (index) - get normal of face given by index;
935 // GetUnitNormal (index) - get unit normal of face given by index;
936 // GetNextNormal (normal) - get normals of each face in order;
937 // returns false when finished all faces;
938 // GetNextUnitNormal (normal) - get normals of unit length of each face
939 // in order; returns false when finished all faces;
940 // GetSurfaceArea() - get surface area of the polyhedron;
941 // GetVolume() - get volume of the polyhedron;
942 // GetNumberOfRotationSteps() - get number of steps for whole circle;
943 // SetNumberOfRotationSteps (n) - set number of steps for whole circle;
944 // ResetNumberOfRotationSteps() - reset number of steps for whole circle
945 // to default value;
946 // History:
947 //
948 // 20.06.96 Evgeni Chernyaev <Evgueni.Tcherniaev@cern.ch> - initial version
949 //
950 // 23.07.96 John Allison
951 // - added GetNoVertices, GetNoFacets, GetNextVertex, GetNextNormal
952 //
953 // 30.09.96 E.Chernyaev
954 // - added GetNextVertexIndex, GetVertex by Yasuhide Sawada
955 // - added GetNextUnitNormal, GetNextEdgeIndeces, GetNextEdge
956 // - improvements: angles now expected in radians
957 // int -> G4int, double -> G4double
958 // - G4ThreeVector replaced by either G4Point3D or G4Normal3D
959 //
960 // 15.12.96 E.Chernyaev
961 // - private functions G4PolyhedronAlloc, G4PolyhedronPrism renamed
962 // to AllocateMemory and CreatePrism
963 // - added private functions GetNumberOfRotationSteps, RotateEdge,
964 // RotateAroundZ, SetReferences
965 // - rewritten G4PolyhedronCons;
966 // - added G4PolyhedronPara, ...Trap, ...Pgon, ...Pcon, ...Sphere, ...Torus,
967 // so full List of implemented shapes now looks like:
968 // BOX, TRD1, TRD2, TRAP, TUBE, TUBS, CONE, CONS, PARA, PGON, PCON,
969 // SPHERE, TORUS
970 //
971 // 01.06.97 E.Chernyaev
972 // - RotateAroundZ modified and SetSideFacets added to allow Rmin=Rmax
973 // in bodies of revolution
974 //
975 // 24.06.97 J.Allison
976 // - added static private member fNumberOfRotationSteps and static public
977 // functions void SetNumberOfRotationSteps (G4int n) and
978 // void ResetNumberOfRotationSteps (). Modified
979 // GetNumberOfRotationSteps() appropriately. Made all three functions
980 // inline (at end of this .hh file).
981 // Usage:
982 // G4Polyhedron::SetNumberOfRotationSteps
983 // (fpView -> GetViewParameters ().GetNoOfSides ());
984 // pPolyhedron = solid.CreatePolyhedron ();
985 // G4Polyhedron::ResetNumberOfRotationSteps ();
986 //
987 // 19.03.00 E.Chernyaev
988 // - added boolean operations (add, subtract, intersect) on polyhedra;
989 //
990 // 25.05.01 E.Chernyaev
991 // - added GetSurfaceArea() and GetVolume();
992 //
993 
tools::hep::polyhedron::polyhedron
polyhedron(int Nvert=0, int Nface=0)
Definition: polyhedron:205
tools::hep::polyhedron_box::~polyhedron_box
virtual ~polyhedron_box()
Definition: polyhedron:564
tools::hep::polyhedron::ResetNumberOfRotationSteps
void ResetNumberOfRotationSteps()
tools::hep::polyhedronProcessor::clear
void clear()
Definition: polyhedron:787
tools::hep::SbFacet::operator==
friend int operator==(const SbFacet &v1, const SbFacet &v2)
tools::hep::polyhedron_tube::polyhedron_tube
polyhedron_tube(const polyhedron_tube &a_from)
Definition: polyhedron:668
tools::hep::polyhedron_pgon::operator=
polyhedron_pgon & operator=(const polyhedron_pgon &a_from)
Definition: polyhedron:687
tools::hep::polyhedron_cons::polyhedron_cons
polyhedron_cons(double Rmn1, double Rmx1, double Rmn2, double Rmx2, double Dz, double Phi1, double Dphi, int nstep=0)
tools::hep::polyhedron::GetFacet
void GetFacet(int iFace, int &n, HVPoint3D *nodes, int *edgeFlags=0, HVNormal3D *normals=0) const
tools::hep::polyhedronProcessor::execute
bool execute(polyhedron &)
tools::hep::SbFacet
Definition: polyhedron:45
tools::hep::polyhedron::set_polyhedron_sphere
bool set_polyhedron_sphere(double rmin, double rmax, double phi, double dphi, double the, double dthe, int nphi=0, int nthe=0)
tools::hep::polyhedron_sphere::polyhedron_sphere
polyhedron_sphere(double rmin, double rmax, double phi, double dphi, double the, double dthe, int nphi=0, int nthe=0)
tools::hep::SbFacet::operator<<
friend std::ostream & operator<<(std::ostream &, const SbFacet &facet)
tools::hep::polyhedron_pgon::~polyhedron_pgon
virtual ~polyhedron_pgon()
Definition: polyhedron:684
tools::hep::polyhedron::set_polyhedron_tubs
bool set_polyhedron_tubs(double Rmin, double Rmax, double Dz, double Phi1, double Dphi, int nstep)
Definition: polyhedron:387
tools::hep::SbFacet::GetEdge
void GetEdge(int i, int &v, int &f) const
Definition: polyhedron:119
tools::hep::polyhedron::SetReferences
void SetReferences()
tools::hep::polyhedron::GetSurfaceArea
double GetSurfaceArea() const
tools::hep::polyhedron_trd1::~polyhedron_trd1
virtual ~polyhedron_trd1()
Definition: polyhedron:549
tools::hep::polyhedron::GetNextVertexIndex
bool GetNextVertexIndex(int &index, int &edgeFlag) const
Definition: polyhedron:446
tools::hep::polyhedron::intersect
polyhedron intersect(const polyhedron &p) const
tools::hep::polyhedron_xtru::~polyhedron_xtru
virtual ~polyhedron_xtru()
Definition: polyhedron:522
tools::hep::polyhedron_box::polyhedron_box
polyhedron_box(double Dx, double Dy, double Dz)
tools::hep::polyhedron::operator=
polyhedron & operator=(const polyhedron &from)
tools::hep::polyhedron_pcon::operator=
polyhedron_pcon & operator=(const polyhedron_pcon &a_from)
Definition: polyhedron:705
tools::hep::polyhedron_trap::operator=
polyhedron_trap & operator=(const polyhedron_trap &a_from)
Definition: polyhedron:586
tools::hep::polyhedron::nface
int nface
Definition: polyhedron:157
tools::hep::polyhedron_sphere::polyhedron_sphere
polyhedron_sphere(const polyhedron_sphere &a_from)
Definition: polyhedron:723
tools::hep::polyhedronProcessor::execute1
bool execute1(polyhedron &, const std::vector< unsigned int > &)
tools::hep::polyhedron::InvertFacets
void InvertFacets()
tools::hep::polyhedron::GetNextEdge
bool GetNextEdge(HVPoint3D &p1, HVPoint3D &p2, int &edgeFlag) const
tools::vec3d
Definition: vec3d:13
tools::hep::polyhedron::GetPF
SbFacet * GetPF() const
Definition: polyhedron:270
tools::vec3::set_value
void set_value(const T &a0, const T &a1, const T &a2)
Definition: vec3:92
tools::hep::polyhedron::~polyhedron
virtual ~polyhedron()
Definition: polyhedron:217
tools::hep::polyhedron_torus::polyhedron_torus
polyhedron_torus(double rmin, double rmax, double rtor, double phi, double dphi, int nphi=0, int nthe=0)
tools::hep::polyhedron_cons::operator=
polyhedron_cons & operator=(const polyhedron_cons &a_from)
Definition: polyhedron:620
tools::hep::polyhedron::Set
void Set(int Nvert, HVPoint3D *aV, int Nface, SbFacet *aF)
Definition: polyhedron:241
tools::hep::polyhedron::GetNextVertex
bool GetNextVertex(HVPoint3D &vertex, int &edgeFlag) const
tools::hep::polyhedron::set_polyhedron_trd2
bool set_polyhedron_trd2(double Dx1, double Dx2, double Dy1, double Dy2, double Dz)
tools::hep::polyhedron::dump
void dump(std::ostream &) const
tools::hep::polyhedron::Translate
polyhedron & Translate(double, double, double)
tools::hep::polyhedron::set_polyhedron_pgon
bool set_polyhedron_pgon(double phi, double dphi, int npdv, int nz, const double *z, const double *rmin, const double *rmax)
tools::hep::SbFacet::set
void set(int v1, int f1, int v2, int f2, int v3, int f3, int v4, int f4)
Definition: polyhedron:124
tools::hep::polyhedron_torus::operator=
polyhedron_torus & operator=(const polyhedron_torus &a_from)
Definition: polyhedron:742
tools::hep::polyhedron_hype::polyhedron_hype
polyhedron_hype(double a_st_in, double a_st_out, double a_rmin, double a_rmax, double a_dz, int a_nz=10, int a_nphi=24)
tools::hep::polyhedronProcessor::Operation
Operation
Definition: polyhedron:757
tools::hep::polyhedronProcessor::polyhedronProcessor
polyhedronProcessor()
Definition: polyhedron:765
tools::hep::polyhedron::FindNodeNormal
HVNormal3D FindNodeNormal(int iFace, int iNode) const
tools::hep::HVNormal3D
HVPoint3D HVNormal3D
Definition: polyhedron:38
tools::hep::SbFacet::operator=
SbFacet & operator=(const SbFacet &aFrom)
Definition: polyhedron:92
tools::hep::polyhedron_para::operator=
polyhedron_para & operator=(const polyhedron_para &a_from)
Definition: polyhedron:602
tools::hep::polyhedron_pcon::polyhedron_pcon
polyhedron_pcon(const polyhedron_pcon &a_from)
Definition: polyhedron:704
tools::hep::polyhedron_trd2::polyhedron_trd2
polyhedron_trd2(double Dx1, double Dx2, double Dy1, double Dy2, double Dz)
tools::hep::polyhedron::GetVolume
double GetVolume() const
tools::hep::polyhedron::GetNextEdgeIndeces
bool GetNextEdgeIndeces(int &i1, int &i2, int &edgeFlag) const
tools::hep::polyhedron
Definition: polyhedron:142
tools::hep::polyhedron_sphere::operator=
polyhedron_sphere & operator=(const polyhedron_sphere &a_from)
Definition: polyhedron:724
tools::hep::polyhedron_tubs::operator=
polyhedron_tubs & operator=(const polyhedron_tubs &a_from)
Definition: polyhedron:654
tools::hep::polyhedron::set_polyhedron_box
bool set_polyhedron_box(double Dx, double Dy, double Dz)
Definition: polyhedron:357
tools::hep::polyhedron::GetVertexFast
const HVPoint3D & GetVertexFast(int index) const
tools::hep::polyhedron::GetNoFacets
int GetNoFacets() const
Definition: polyhedron:253
tools::vec3::y
const T & y() const
Definition: vec3:85
tools::hep::polyhedron::set_polyhedron_pcon
bool set_polyhedron_pcon(double phi, double dphi, int nz, const double *z, const double *rmin, const double *rmax)
Definition: polyhedron:380
tools::hep::polyhedronProcessor::UNION
@ UNION
Definition: polyhedron:758
tools::hep::SbFacet::Set
void Set(int v[8])
Definition: polyhedron:131
tools::hep::polyhedron_trd2::operator=
polyhedron_trd2 & operator=(const polyhedron_trd2 &a_from)
Definition: polyhedron:491
tools::hep::polyhedron_trap::~polyhedron_trap
virtual ~polyhedron_trap()
Definition: polyhedron:583
tools::hep::polyhedron::Transform
polyhedron & Transform(const rotd &rot, double, double, double)
tools::hep::polyhedronProcessor::~polyhedronProcessor
virtual ~polyhedronProcessor()
Definition: polyhedron:770
tools::hep::polyhedron_cone::~polyhedron_cone
virtual ~polyhedron_cone()
Definition: polyhedron:634
TOOLS_SCLASS
#define TOOLS_SCLASS(a_name)
Definition: S_STRING:41
tools::hep::polyhedron_torus::~polyhedron_torus
virtual ~polyhedron_torus()
Definition: polyhedron:739
tools::hep::polyhedron::set_polyhedron_torus
bool set_polyhedron_torus(double rmin, double rmax, double rtor, double phi, double dphi, int nphi, int nthe)
tools::hep::SbFacet::SbFacet
SbFacet(int v1=0, int f1=0, int v2=0, int f2=0, int v3=0, int f3=0, int v4=0, int f4=0)
Definition: polyhedron:63
tools::hep::polyhedron::polyhedron
polyhedron(const polyhedron &from)
tools::hep::polyhedron::set_polyhedron_hype
bool set_polyhedron_hype(double a_st_in, double a_st_out, double a_rmin, double a_rmax, double a_dz, int a_nz=10, int a_nphi=24)
tools::hep::polyhedron::Transform
polyhedron & Transform(const rotd &rot, const vec3d &trans)
tools::hep::Sb_iabs
int Sb_iabs(int a)
Definition: polyhedron:441
tools::hep::polyhedron_torus
Definition: polyhedron:733
tools::vec3::x
const T & x() const
Definition: vec3:84
tools::hep::polyhedron_box
Definition: polyhedron:561
tools::hep::polyhedron::SetSideFacets
void SetSideFacets(int ii[4], int vv[4], int *kk, double *r, double dphi, int ns, int &kface)
tools::hep::polyhedron::RotateEdge
void RotateEdge(int k1, int k2, double r1, double r2, int v1, int v2, int vEdge, bool ifWholeCircle, int ns, int &kface)
tools::hep::polyhedron_tubs::~polyhedron_tubs
virtual ~polyhedron_tubs()
Definition: polyhedron:651
polyhedron.icc
tools::hep::polyhedron_xtru
Definition: polyhedron:516
tools::hep::HVPoint3D
vec3d HVPoint3D
Definition: polyhedron:37
tools::hep::polyhedron_arb8::polyhedron_arb8
polyhedron_arb8(const polyhedron_arb8 &a_from)
Definition: polyhedron:506
tools::hep::polyhedron::subtract
polyhedron subtract(const polyhedron &p) const
tools::hep::polyhedron_pcon::~polyhedron_pcon
virtual ~polyhedron_pcon()
Definition: polyhedron:702
tools::hep::polyhedron::GetNextEdge
bool GetNextEdge(HVPoint3D &p1, HVPoint3D &p2, int &edgeFlag, int &iface1, int &iface2) const
tools::hep::polyhedron_pcon
Definition: polyhedron:696
tools::hep::polyhedron_pgon
Definition: polyhedron:678
tools::hep::SbFacet::~SbFacet
virtual ~SbFacet()
Definition: polyhedron:72
tools::hep::tsf_polyhedron
void tsf_polyhedron(polyhedron &a_ph, const MATRIX &a_matrix)
Definition: polyhedron:818
tools::hep::polyhedron_sphere::~polyhedron_sphere
virtual ~polyhedron_sphere()
Definition: polyhedron:721
tools::hep::polyhedron::do_not_set_NUMBER_OF_STEPS
static void do_not_set_NUMBER_OF_STEPS(int)
Definition: polyhedron:203
tools::hep::polyhedron::_M_PI
static double _M_PI()
Definition: polyhedron:164
tools::hep::polyhedron::set_polyhedron_cons
bool set_polyhedron_cons(double Rmn1, double Rmx1, double Rmn2, double Rmx2, double Dz, double Phi1, double Dphi, int nstep=0)
tools::hep::polyhedron_xtru::polyhedron_xtru
polyhedron_xtru(const polyhedron_xtru &a_from)
Definition: polyhedron:524
tools::hep::polyhedron::set_polyhedron_arb8
bool set_polyhedron_arb8(double Dz, const double *xy)
tools::hep::SbFacet::operator!=
friend int operator!=(const SbFacet &v1, const SbFacet &v2)
tools::hep::polyhedron_arb8::operator=
polyhedron_arb8 & operator=(const polyhedron_arb8 &a_from)
Definition: polyhedron:507
tools::hep::polyhedron_arb8::polyhedron_arb8
polyhedron_arb8(double Dz, const double *xy)
tools::hep::polyhedron_trd1::polyhedron_trd1
polyhedron_trd1(double Dx1, double Dx2, double Dy, double Dz)
tools::hep::polyhedron::operator!=
friend int operator!=(const polyhedron &v1, const polyhedron &v2)
tools::hep::polyhedron::set_polyhedron_trap
bool set_polyhedron_trap(double Dz, double Theta, double Phi, double Dy1, double Dx1, double Dx2, double Alp1, double Dy2, double Dx3, double Dx4, double Alp2)
tools::hep::polyhedron_box::operator=
polyhedron_box & operator=(const polyhedron_box &a_from)
Definition: polyhedron:567
tools::hep::polyhedron_para::polyhedron_para
polyhedron_para(const polyhedron_para &a_from)
Definition: polyhedron:601
tools::hep::polyhedron_box::polyhedron_box
polyhedron_box(const polyhedron_box &a_from)
Definition: polyhedron:566
tools::hep::polyhedron::set_polyhedron_cone
bool set_polyhedron_cone(double Rmn1, double Rmx1, double Rmn2, double Rmx2, double Dz, int nstep)
Definition: polyhedron:394
tools::hep::polyhedronProcessor
Definition: polyhedron:752
tools::hep::polyhedron::pF
SbFacet * pF
Definition: polyhedron:159
tools::hep::polyhedron_tube::polyhedron_tube
polyhedron_tube(double Rmin, double Rmax, double Dz, int nstep=0)
tools::hep::polyhedron_xtru::operator=
polyhedron_xtru & operator=(const polyhedron_xtru &a_from)
Definition: polyhedron:525
tools
inlined C code : ///////////////////////////////////
Definition: aida_ntuple:26
tools::hep::polyhedron::FindNeighbour
int FindNeighbour(int iFace, int iNode, int iOrder) const
tools::hep::polyhedron::set_polyhedron_para
bool set_polyhedron_para(double Dx, double Dy, double Dz, double Alpha, double Theta, double Phi)
Definition: polyhedron:370
tools::hep::polyhedron::set_polyhedron_tube
bool set_polyhedron_tube(double Rmin, double Rmax, double Dz, int nstep=0)
Definition: polyhedron:347
tools::hep::polyhedron_trd2
Definition: polyhedron:484
tools::hep::polyhedron_para
Definition: polyhedron:595
tools::hep::polyhedron_trap
Definition: polyhedron:576
tools::hep::polyhedron_arb8::~polyhedron_arb8
virtual ~polyhedron_arb8()
Definition: polyhedron:504
tools::hep::polyhedron::GetNoVertices
int GetNoVertices() const
Definition: polyhedron:250
tools::hep::polyhedron_hype::polyhedron_hype
polyhedron_hype(const polyhedron_hype &a_from)
Definition: polyhedron:538
tools::rotd
Definition: rotd:17
tools::hep::polyhedron_tubs::polyhedron_tubs
polyhedron_tubs(double Rmin, double Rmax, double Dz, double Phi1, double Dphi, int nstep=0)
tools::hep::polyhedron::GetNextFacet
bool GetNextFacet(int &n, HVPoint3D *nodes, int *edgeFlags=0, HVNormal3D *normals=0) const
tools::hep::polyhedron_sphere
Definition: polyhedron:714
tools::hep::SbFacet::SbFacet
SbFacet(const SbFacet &aFrom)
Definition: polyhedron:79
tools::hep::polyhedron::GetNextVertex
bool GetNextVertex(HVPoint3D &vertex, int &edgeFlag, HVNormal3D &normal) const
tools::hep::polyhedron::set_polyhedron_eltu
bool set_polyhedron_eltu(double a_dx, double a_dy, double a_dz, int a_nz=10, int a_nphi=24)
tools::hep::polyhedron_pcon::polyhedron_pcon
polyhedron_pcon(double phi, double dphi, int nz, const double *z, const double *rmin, const double *rmax)
tools::hep::polyhedron_torus::polyhedron_torus
polyhedron_torus(const polyhedron_torus &a_from)
Definition: polyhedron:741
tools::hep::polyhedron::set_polyhedron_xtru
bool set_polyhedron_xtru(int a_npts, int a_nz, double *a_xs, double *a_ys, double *a_zs, bool a_acw=true, bool a_zfb=true)
tools::hep::polyhedron::GetNormal
HVNormal3D GetNormal(int iFace) const
tools::hep::polyhedron::isConsistent
bool isConsistent(const char *=0) const
tools::hep::polyhedron_trd2::polyhedron_trd2
polyhedron_trd2(const polyhedron_trd2 &a_from)
Definition: polyhedron:490
tools::hep::polyhedron::GetFacet
void GetFacet(int iFace, int &n, int *iNodes, int *edgeFlags=0, int *iFaces=0) const
tools::hep::polyhedron::GetNextEdgeIndeces
bool GetNextEdgeIndeces(int &i1, int &i2, int &edgeFlag, int &iface1, int &iface2) const
tools::hep::polyhedron_tubs::polyhedron_tubs
polyhedron_tubs(const polyhedron_tubs &a_from)
Definition: polyhedron:653
tools::hep::polyhedron_cons::polyhedron_cons
polyhedron_cons(const polyhedron_cons &a_from)
Definition: polyhedron:619
tools::hep::polyhedron_tube::~polyhedron_tube
virtual ~polyhedron_tube()
Definition: polyhedron:666
tools::hep::polyhedron::GetNextUnitNormal
bool GetNextUnitNormal(HVNormal3D &normal) const
tools::hep::polyhedron::nvert
int nvert
Definition: polyhedron:157
tools::hep::polyhedron_trap::polyhedron_trap
polyhedron_trap(double Dz, double Theta, double Phi, double Dy1, double Dx1, double Dx2, double Alp1, double Dy2, double Dx3, double Dx4, double Alp2)
tools::hep::polyhedron::AllocateMemory
void AllocateMemory(int Nvert, int Nface)
tools::hep::polyhedron_cone::polyhedron_cone
polyhedron_cone(const polyhedron_cone &a_from)
Definition: polyhedron:636
tools::hep::polyhedron_trd2::~polyhedron_trd2
virtual ~polyhedron_trd2()
Definition: polyhedron:488
tools::hep::polyhedron_xtru::polyhedron_xtru
polyhedron_xtru(int a_npts, int a_nz, double *a_xs, double *a_ys, double *a_zs, bool a_acw=true, bool a_zfb=true)
tools::hep::polyhedronProcessor::push_back
void push_back(Operation a_op, const polyhedron &a_polyhedron)
Definition: polyhedron:783
tools::vec3::z
const T & z() const
Definition: vec3:86
tools::hep::polyhedron::SetNumberOfRotationSteps
void SetNumberOfRotationSteps(int n)
tools::hep::polyhedronProcessor::SUBTRACTION
@ SUBTRACTION
Definition: polyhedron:760
tools::hep::polyhedron_trd1::polyhedron_trd1
polyhedron_trd1(const polyhedron_trd1 &a_from)
Definition: polyhedron:551
tools::hep::polyhedron::set_polyhedron_trd1
bool set_polyhedron_trd1(double Dx1, double Dx2, double Dy, double Dz)
Definition: polyhedron:360
tools::hep::polyhedron::operator<<
friend std::ostream & operator<<(std::ostream &, const polyhedron &ph)
tools::hep::polyhedron::GetVertex
HVPoint3D GetVertex(int index) const
tools::hep::polyhedron_para::polyhedron_para
polyhedron_para(double Dx, double Dy, double Dz, double Alpha, double Theta, double Phi)
tools::hep::polyhedron::isEqual
bool isEqual(const polyhedron &p) const
tools::hep::polyhedron::GetNextNormal
bool GetNextNormal(HVNormal3D &normal) const
tools::hep::polyhedron_tube::operator=
polyhedron_tube & operator=(const polyhedron_tube &a_from)
Definition: polyhedron:669
tools::hep::polyhedron::GetUnitNormal
HVNormal3D GetUnitNormal(int iFace) const
tools::hep::polyhedron::GetPV
HVPoint3D * GetPV() const
Definition: polyhedron:269
tools::hep::polyhedron_hype::operator=
polyhedron_hype & operator=(const polyhedron_hype &a_from)
Definition: polyhedron:539
tools::hep::polyhedron::CreatePrism
void CreatePrism()
tools::hep::polyhedron_cone::operator=
polyhedron_cone & operator=(const polyhedron_cone &a_from)
Definition: polyhedron:637
tools::hep::polyhedron_cons::~polyhedron_cons
virtual ~polyhedron_cons()
Definition: polyhedron:617
tools::hep::polyhedron_trd1
Definition: polyhedron:545
tools::hep::polyhedron::NUMBER_OF_STEPS
static int NUMBER_OF_STEPS()
Definition: polyhedron:201
tools::hep::polyhedron::Empty
void Empty()
Definition: polyhedron:246
tools::hep::polyhedron_hype::~polyhedron_hype
virtual ~polyhedron_hype()
Definition: polyhedron:536
tools::hep::polyhedron_arb8
Definition: polyhedron:501
tools::hep::polyhedron_hype
Definition: polyhedron:531
tools::hep::polyhedronProcessor::INTERSECTION
@ INTERSECTION
Definition: polyhedron:759
tools::hep::polyhedron_pgon::polyhedron_pgon
polyhedron_pgon(const polyhedron_pgon &a_from)
Definition: polyhedron:686
tools::hep::polyhedron_trap::polyhedron_trap
polyhedron_trap(const polyhedron_trap &a_from)
Definition: polyhedron:585
tools::hep::polyhedron::operator==
friend int operator==(const polyhedron &v1, const polyhedron &v2)
tools::hep::polyhedron_cone::polyhedron_cone
polyhedron_cone(double Rmn1, double Rmx1, double Rmn2, double Rmx2, double Dz, int nstep=0)
tools::hep::polyhedron_tube
Definition: polyhedron:663
tools::hep::polyhedron::vxy
double vxy(const double *xy, int i, int j)
Definition: polyhedron:352
tools::hep::polyhedron::RotateAroundZ
void RotateAroundZ(int nstep, double phi, double dphi, int np1, int np2, const double *z, double *r, int nodeVis, int edgeVis)
tools::hep::polyhedron_pgon::polyhedron_pgon
polyhedron_pgon(double phi, double dphi, int npdv, int nz, const double *z, const double *rmin, const double *rmax)
tools::hep::polyhedronProcessor::is_same_op
bool is_same_op() const
Definition: polyhedron:788
tools::hep::polyhedron_tubs
Definition: polyhedron:646
tools::hep::SbFacet::isEqual
bool isEqual(const SbFacet &aFrom) const
Definition: polyhedron:104
tools::hep::polyhedron_cone
Definition: polyhedron:629
tools::hep::polyhedron::GetNumberOfRotationSteps
int GetNumberOfRotationSteps()
tools::hep::polyhedron::pV
HVPoint3D * pV
Definition: polyhedron:158
tools::hep::polyhedron_para::~polyhedron_para
virtual ~polyhedron_para()
Definition: polyhedron:599
tools::hep::polyhedron_cons
Definition: polyhedron:611
tools::hep::polyhedron::add
polyhedron add(const polyhedron &p) const
tools::hep::HVVector3D
HVPoint3D HVVector3D
Definition: polyhedron:39
tools::hep::polyhedron_trd1::operator=
polyhedron_trd1 & operator=(const polyhedron_trd1 &a_from)
Definition: polyhedron:552