// @(#)root/geom:$Id$
// Author: Andrei Gheata   24/10/01

/*************************************************************************
 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

#ifndef ROOT_TGeoArb8
#define ROOT_TGeoArb8

#ifndef ROOT_TGeoBBox
#include "TGeoBBox.h"
#endif


////////////////////////////////////////////////////////////////////////////
//                                                                        //
// TGeoArb8 - a arbitrary trapezoid with less than 8 vertices standing on //
//   two paralel planes perpendicular to Z axis. Parameters :             //
//            - dz - half length in Z;                                    //
//            - xy[8][2] - vector of (x,y) coordinates of vertices        //
//               - first four points (xy[i][j], i<4, j<2) are the (x,y)   //
//                 coordinates of the vertices sitting on the -dz plane;  //
//               - last four points (xy[i][j], i>=4, j<2) are the (x,y)   //
//                 coordinates of the vertices sitting on the +dz plane;  //
//   The order of defining the vertices of an arb8 is the following :     //
//      - point 0 is connected with points 1,3,4                          //
//      - point 1 is connected with points 0,2,5                          //
//      - point 2 is connected with points 1,3,6                          //
//      - point 3 is connected with points 0,2,7                          //
//      - point 4 is connected with points 0,5,7                          //
//      - point 5 is connected with points 1,4,6                          //
//      - point 6 is connected with points 2,5,7                          //
//      - point 7 is connected with points 3,4,6                          //
//   Points can be identical in order to create shapes with less than     //
//   8 vertices.                                                          //
//                                                                        //
////////////////////////////////////////////////////////////////////////////

class TGeoArb8 : public TGeoBBox
{
protected:
   enum EGeoArb8Type {
//      kArb8Trd1 = BIT(25), // trd1 type
//      kArb8Trd2 = BIT(26), // trd2 type
      kArb8Trap = BIT(27), // planar surface trapezoid
      kArb8Tra  = BIT(28)  // general twisted trapezoid
   };
   // data members
   Double_t              fDz;          // half length in Z
   Double_t             *fTwist;       //! [4] tangents of twist angles
   Double_t              fXY[8][2];    // list of vertices

   TGeoArb8(const TGeoArb8&);
   TGeoArb8& operator=(const TGeoArb8&);

public:
   // constructors
   TGeoArb8();
   TGeoArb8(Double_t dz, Double_t *vertices=0);
   TGeoArb8(const char *name, Double_t dz, Double_t *vertices=0);
   // destructor
   virtual ~TGeoArb8();
   // methods
   virtual Double_t      Capacity() const;
   virtual void          ComputeBBox();
   virtual void          ComputeNormal(const Double_t *point, const Double_t *dir, Double_t *norm);
   virtual void          ComputeNormal_v(const Double_t *points, const Double_t *dirs, Double_t *norms, Int_t vecsize);
   void                  ComputeTwist();
   virtual Bool_t        Contains(const Double_t *point) const;
   virtual void          Contains_v(const Double_t *points, Bool_t *inside, Int_t vecsize) const;
   Double_t              DistToPlane(const Double_t *point, const Double_t *dir, Int_t ipl, Bool_t in) const;
   virtual Double_t      DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1,
                                   Double_t step=TGeoShape::Big(), Double_t *safe=0) const;
   virtual void          DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const;
   virtual Double_t      DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1,
                                   Double_t step=TGeoShape::Big(), Double_t *safe=0) const;
   virtual void          DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const;
   virtual TGeoVolume   *Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
                                Double_t start, Double_t step);
   virtual Double_t      GetAxisRange(Int_t iaxis, Double_t &xlo, Double_t &xhi) const;
   virtual void          GetBoundingCylinder(Double_t *param) const;
   virtual Int_t         GetByteCount() const {return 100;}
   Double_t              GetClosestEdge(const Double_t *point, Double_t *vert, Int_t &isegment) const;
   virtual Bool_t        GetPointsOnFacet(Int_t /*index*/, Int_t /*npoints*/, Double_t * /*array*/) const;
   Double_t              GetDz() const {return fDz;}
   virtual Int_t         GetFittingBox(const TGeoBBox *parambox, TGeoMatrix *mat, Double_t &dx, Double_t &dy, Double_t &dz) const;
   virtual TGeoShape    *GetMakeRuntimeShape(TGeoShape * /*mother*/, TGeoMatrix * /*mat*/) const {return 0;}
   static void           GetPlaneNormal(Double_t *p1, Double_t *p2, Double_t *p3, Double_t *norm);
   Double_t             *GetVertices() {return &fXY[0][0];}
   Double_t              GetTwist(Int_t iseg) const;
   virtual Bool_t        IsCylType() const {return kFALSE;}
   static Bool_t         IsSamePoint(const Double_t *p1, const Double_t *p2) {return (TMath::Abs(p1[0]-p2[0])<1.E-16 && TMath::Abs(p1[1]-p2[1])<1.E-16)?kTRUE:kFALSE;}
   static Bool_t         InsidePolygon(Double_t x, Double_t y, Double_t *pts);
   virtual void          InspectShape() const;
   Bool_t                IsTwisted() const {return (fTwist==0)?kFALSE:kTRUE;}
   Double_t              SafetyToFace(const Double_t *point, Int_t iseg, Bool_t in) const;
   virtual Double_t      Safety(const Double_t *point, Bool_t in=kTRUE) const;
   virtual void          Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const;
   virtual void          SavePrimitive(std::ostream &out, Option_t *option = "");
   void                  SetPlaneVertices(Double_t zpl, Double_t *vertices) const;
   virtual void          SetVertex(Int_t vnum, Double_t x, Double_t y);
   virtual void          SetDimensions(Double_t *param);
   void                  SetDz(Double_t dz) {fDz = dz;}
   virtual void          SetPoints(Double_t *points) const;
   virtual void          SetPoints(Float_t *points) const;
   virtual void          Sizeof3D() const;

   ClassDef(TGeoArb8, 1)         // arbitrary trapezoid with 8 vertices
};

////////////////////////////////////////////////////////////////////////////
//                                                                        //
// TGeoTrap                                                               //
//                                                                        //
// Trap is a general trapezoid, i.e. one for which the faces perpendicular//
// to z are trapezia and their centres are not the same x, y. It has 11   //
// parameters: the half length in z, the polar angles from the centre of  //
// the face at low z to that at high z, H1 the half length in y at low z, //
// LB1 the half length in x at low z and y low edge, LB2 the half length  //
// in x at low z and y high edge, TH1 the angle w.r.t. the y axis from the//
// centre of low y edge to the centre of the high y edge, and H2, LB2,    //
// LH2, TH2, the corresponding quantities at high z.                      //
//                                                                        //
////////////////////////////////////////////////////////////////////////////

class TGeoTrap : public TGeoArb8
{
protected:
   // data members
   Double_t              fTheta; // theta angle
   Double_t              fPhi;   // phi angle
   Double_t              fH1;    // half length in y at low z
   Double_t              fBl1;   // half length in x at low z and y low edge
   Double_t              fTl1;   // half length in x at low z and y high edge
   Double_t              fAlpha1;// angle between centers of x edges an y axis at low z
   Double_t              fH2;    // half length in y at high z
   Double_t              fBl2;   // half length in x at high z and y low edge
   Double_t              fTl2;   // half length in x at high z and y high edge
   Double_t              fAlpha2;// angle between centers of x edges an y axis at low z

public:
   // constructors
   TGeoTrap();
   TGeoTrap(Double_t dz, Double_t theta, Double_t phi);
   TGeoTrap(Double_t dz, Double_t theta, Double_t phi, Double_t h1,
            Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
            Double_t tl2, Double_t alpha2);
   TGeoTrap(const char *name, Double_t dz, Double_t theta, Double_t phi, Double_t h1,
            Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
            Double_t tl2, Double_t alpha2);
   // destructor
   virtual ~TGeoTrap();
   virtual Double_t      DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1,
                                   Double_t step=TGeoShape::Big(), Double_t *safe=0) const;
   virtual void          DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const;
   virtual Double_t      DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1,
                                   Double_t step=TGeoShape::Big(), Double_t *safe=0) const;
   virtual void          DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const;
   virtual TGeoVolume   *Divide(TGeoVolume *voldiv, const char *divname, Int_t iaxis, Int_t ndiv,
                                Double_t start, Double_t step);
   Double_t              GetTheta() const {return fTheta;}
   Double_t              GetPhi() const   {return fPhi;}
   Double_t              GetH1() const    {return fH1;}
   Double_t              GetBl1() const   {return fBl1;}
   Double_t              GetTl1() const   {return fTl1;}
   Double_t              GetAlpha1() const   {return fAlpha1;}
   Double_t              GetH2() const    {return fH2;}
   Double_t              GetBl2() const   {return fBl2;}
   Double_t              GetTl2() const   {return fTl2;}
   Double_t              GetAlpha2() const   {return fAlpha2;}
   virtual TGeoShape    *GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const;
   virtual void          SetDimensions(Double_t *param);
   virtual Double_t      Safety(const Double_t *point, Bool_t in=kTRUE) const;
   virtual void          Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const;
   virtual void          SavePrimitive(std::ostream &out, Option_t *option = "");

   ClassDef(TGeoTrap, 1)         // G3 TRAP shape
};

////////////////////////////////////////////////////////////////////////////
//                                                                        //
// TGeoGtra                                                               //
//                                                                        //
// Gtra is a twisted general trapezoid, i.e. one for which the faces perpendicular//
// to z are trapezia and their centres are not the same x, y. It has 12   //
// parameters: the half length in z, the polar angles from the centre of  //
// the face at low z to that at high z, the twist angle, H1 the half length in y at low z, //
// LB1 the half length in x at low z and y low edge, LB2 the half length  //
// in x at low z and y high edge, TH1 the angle w.r.t. the y axis from the//
// centre of low y edge to the centre of the high y edge, and H2, LB2,    //
// LH2, TH2, the corresponding quantities at high z.                      //
//                                                                        //
////////////////////////////////////////////////////////////////////////////

class TGeoGtra : public TGeoTrap
{
protected:
   // data members
   Double_t          fTwistAngle; // twist angle in degrees
public:
   // constructors
   TGeoGtra();
   TGeoGtra(Double_t dz, Double_t theta, Double_t phi, Double_t twist, Double_t h1,
            Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
            Double_t tl2, Double_t alpha2);
   TGeoGtra(const char *name, Double_t dz, Double_t theta, Double_t phi, Double_t twist, Double_t h1,
            Double_t bl1, Double_t tl1, Double_t alpha1, Double_t h2, Double_t bl2,
            Double_t tl2, Double_t alpha2);
   // destructor
   virtual ~TGeoGtra();
   virtual Double_t      DistFromInside(const Double_t *point, const Double_t *dir, Int_t iact=1,
                                   Double_t step=TGeoShape::Big(), Double_t *safe=0) const;
   virtual void          DistFromInside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const;
   virtual Double_t      DistFromOutside(const Double_t *point, const Double_t *dir, Int_t iact=1,
                                   Double_t step=TGeoShape::Big(), Double_t *safe=0) const;
   virtual void          DistFromOutside_v(const Double_t *points, const Double_t *dirs, Double_t *dists, Int_t vecsize, Double_t *step) const;
   virtual TGeoShape    *GetMakeRuntimeShape(TGeoShape *mother, TGeoMatrix *mat) const;
   Double_t              GetTwistAngle() const {return fTwistAngle;}
   virtual Double_t      Safety(const Double_t *point, Bool_t in=kTRUE) const;
   virtual void          Safety_v(const Double_t *points, const Bool_t *inside, Double_t *safe, Int_t vecsize) const;
   virtual void          SetDimensions(Double_t *param);
   virtual void          SavePrimitive(std::ostream &out, Option_t *option = "");

   ClassDef(TGeoGtra, 1)         // G3 GTRA shape
};

#endif
 TGeoArb8.h:1
 TGeoArb8.h:2
 TGeoArb8.h:3
 TGeoArb8.h:4
 TGeoArb8.h:5
 TGeoArb8.h:6
 TGeoArb8.h:7
 TGeoArb8.h:8
 TGeoArb8.h:9
 TGeoArb8.h:10
 TGeoArb8.h:11
 TGeoArb8.h:12
 TGeoArb8.h:13
 TGeoArb8.h:14
 TGeoArb8.h:15
 TGeoArb8.h:16
 TGeoArb8.h:17
 TGeoArb8.h:18
 TGeoArb8.h:19
 TGeoArb8.h:20
 TGeoArb8.h:21
 TGeoArb8.h:22
 TGeoArb8.h:23
 TGeoArb8.h:24
 TGeoArb8.h:25
 TGeoArb8.h:26
 TGeoArb8.h:27
 TGeoArb8.h:28
 TGeoArb8.h:29
 TGeoArb8.h:30
 TGeoArb8.h:31
 TGeoArb8.h:32
 TGeoArb8.h:33
 TGeoArb8.h:34
 TGeoArb8.h:35
 TGeoArb8.h:36
 TGeoArb8.h:37
 TGeoArb8.h:38
 TGeoArb8.h:39
 TGeoArb8.h:40
 TGeoArb8.h:41
 TGeoArb8.h:42
 TGeoArb8.h:43
 TGeoArb8.h:44
 TGeoArb8.h:45
 TGeoArb8.h:46
 TGeoArb8.h:47
 TGeoArb8.h:48
 TGeoArb8.h:49
 TGeoArb8.h:50
 TGeoArb8.h:51
 TGeoArb8.h:52
 TGeoArb8.h:53
 TGeoArb8.h:54
 TGeoArb8.h:55
 TGeoArb8.h:56
 TGeoArb8.h:57
 TGeoArb8.h:58
 TGeoArb8.h:59
 TGeoArb8.h:60
 TGeoArb8.h:61
 TGeoArb8.h:62
 TGeoArb8.h:63
 TGeoArb8.h:64
 TGeoArb8.h:65
 TGeoArb8.h:66
 TGeoArb8.h:67
 TGeoArb8.h:68
 TGeoArb8.h:69
 TGeoArb8.h:70
 TGeoArb8.h:71
 TGeoArb8.h:72
 TGeoArb8.h:73
 TGeoArb8.h:74
 TGeoArb8.h:75
 TGeoArb8.h:76
 TGeoArb8.h:77
 TGeoArb8.h:78
 TGeoArb8.h:79
 TGeoArb8.h:80
 TGeoArb8.h:81
 TGeoArb8.h:82
 TGeoArb8.h:83
 TGeoArb8.h:84
 TGeoArb8.h:85
 TGeoArb8.h:86
 TGeoArb8.h:87
 TGeoArb8.h:88
 TGeoArb8.h:89
 TGeoArb8.h:90
 TGeoArb8.h:91
 TGeoArb8.h:92
 TGeoArb8.h:93
 TGeoArb8.h:94
 TGeoArb8.h:95
 TGeoArb8.h:96
 TGeoArb8.h:97
 TGeoArb8.h:98
 TGeoArb8.h:99
 TGeoArb8.h:100
 TGeoArb8.h:101
 TGeoArb8.h:102
 TGeoArb8.h:103
 TGeoArb8.h:104
 TGeoArb8.h:105
 TGeoArb8.h:106
 TGeoArb8.h:107
 TGeoArb8.h:108
 TGeoArb8.h:109
 TGeoArb8.h:110
 TGeoArb8.h:111
 TGeoArb8.h:112
 TGeoArb8.h:113
 TGeoArb8.h:114
 TGeoArb8.h:115
 TGeoArb8.h:116
 TGeoArb8.h:117
 TGeoArb8.h:118
 TGeoArb8.h:119
 TGeoArb8.h:120
 TGeoArb8.h:121
 TGeoArb8.h:122
 TGeoArb8.h:123
 TGeoArb8.h:124
 TGeoArb8.h:125
 TGeoArb8.h:126
 TGeoArb8.h:127
 TGeoArb8.h:128
 TGeoArb8.h:129
 TGeoArb8.h:130
 TGeoArb8.h:131
 TGeoArb8.h:132
 TGeoArb8.h:133
 TGeoArb8.h:134
 TGeoArb8.h:135
 TGeoArb8.h:136
 TGeoArb8.h:137
 TGeoArb8.h:138
 TGeoArb8.h:139
 TGeoArb8.h:140
 TGeoArb8.h:141
 TGeoArb8.h:142
 TGeoArb8.h:143
 TGeoArb8.h:144
 TGeoArb8.h:145
 TGeoArb8.h:146
 TGeoArb8.h:147
 TGeoArb8.h:148
 TGeoArb8.h:149
 TGeoArb8.h:150
 TGeoArb8.h:151
 TGeoArb8.h:152
 TGeoArb8.h:153
 TGeoArb8.h:154
 TGeoArb8.h:155
 TGeoArb8.h:156
 TGeoArb8.h:157
 TGeoArb8.h:158
 TGeoArb8.h:159
 TGeoArb8.h:160
 TGeoArb8.h:161
 TGeoArb8.h:162
 TGeoArb8.h:163
 TGeoArb8.h:164
 TGeoArb8.h:165
 TGeoArb8.h:166
 TGeoArb8.h:167
 TGeoArb8.h:168
 TGeoArb8.h:169
 TGeoArb8.h:170
 TGeoArb8.h:171
 TGeoArb8.h:172
 TGeoArb8.h:173
 TGeoArb8.h:174
 TGeoArb8.h:175
 TGeoArb8.h:176
 TGeoArb8.h:177
 TGeoArb8.h:178
 TGeoArb8.h:179
 TGeoArb8.h:180
 TGeoArb8.h:181
 TGeoArb8.h:182
 TGeoArb8.h:183
 TGeoArb8.h:184
 TGeoArb8.h:185
 TGeoArb8.h:186
 TGeoArb8.h:187
 TGeoArb8.h:188
 TGeoArb8.h:189
 TGeoArb8.h:190
 TGeoArb8.h:191
 TGeoArb8.h:192
 TGeoArb8.h:193
 TGeoArb8.h:194
 TGeoArb8.h:195
 TGeoArb8.h:196
 TGeoArb8.h:197
 TGeoArb8.h:198
 TGeoArb8.h:199
 TGeoArb8.h:200
 TGeoArb8.h:201
 TGeoArb8.h:202
 TGeoArb8.h:203
 TGeoArb8.h:204
 TGeoArb8.h:205
 TGeoArb8.h:206
 TGeoArb8.h:207
 TGeoArb8.h:208
 TGeoArb8.h:209
 TGeoArb8.h:210
 TGeoArb8.h:211
 TGeoArb8.h:212
 TGeoArb8.h:213
 TGeoArb8.h:214
 TGeoArb8.h:215
 TGeoArb8.h:216
 TGeoArb8.h:217
 TGeoArb8.h:218
 TGeoArb8.h:219
 TGeoArb8.h:220
 TGeoArb8.h:221
 TGeoArb8.h:222
 TGeoArb8.h:223
 TGeoArb8.h:224
 TGeoArb8.h:225
 TGeoArb8.h:226
 TGeoArb8.h:227
 TGeoArb8.h:228
 TGeoArb8.h:229
 TGeoArb8.h:230
 TGeoArb8.h:231
 TGeoArb8.h:232