G4TriangularFacet.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration and of QinetiQ Ltd,   *
00020 // * subject to DEFCON 705 IPR conditions.                            *
00021 // * By using,  copying,  modifying or  distributing the software (or *
00022 // * any work based  on the software)  you  agree  to acknowledge its *
00023 // * use  in  resulting  scientific  publications,  and indicate your *
00024 // * acceptance of all terms of the Geant4 Software license.          *
00025 // ********************************************************************
00026 //
00027 // $Id: G4TriangularFacet.hh 67011 2013-01-29 16:17:41Z gcosmo $
00028 //
00029 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00030 //
00031 // Class G4TriangularFacet
00032 //
00033 // Class description:
00034 //
00035 //   The G4TriangularFacet class is used for the contruction of
00036 //   G4TessellatedSolid.
00037 //   It is defined by three fVertices, which shall be supplied in anti-clockwise
00038 //   order looking from the outsider of the solid where it belongs.
00039 //   Its constructor:
00040 //   
00041 //      G4TriangularFacet (const G4ThreeVector Pt0, const G4ThreeVector vt1,
00042 //                         const G4ThreeVector vt2, G4FacetVertexType);
00043 //
00044 //   takes 4 parameters to define the three fVertices:
00045 //      1) G4FacetvertexType = "ABSOLUTE": in this case Pt0, vt1 and vt2 are 
00046 //         the 3 fVertices in anti-clockwise order looking from the outsider.
00047 //      2) G4FacetvertexType = "RELATIVE": in this case the first vertex is Pt0,
00048 //         the second vertex is Pt0+vt1 and the third vertex is Pt0+vt2, all  
00049 //         in anti-clockwise order when looking from the outsider.
00050 
00051 // CHANGE HISTORY
00052 // --------------
00053 //
00054 // 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
00055 // 12 October 2012, M Gayer, CERN, - Reviewed optimized implementation.
00056 //
00058 #ifndef G4TriangularFacet_hh
00059 #define G4TriangularFacet_hh 1
00060 
00061 #include "G4VFacet.hh"
00062 #include "G4Types.hh"
00063 #include "G4ThreeVector.hh"
00064 
00065 class G4TriangularFacet : public G4VFacet
00066 {
00067   public:  // with desctiption
00068 
00069     G4TriangularFacet ();
00070    ~G4TriangularFacet ();
00071 
00072     G4TriangularFacet (const G4ThreeVector &vt0, const G4ThreeVector &vt1,
00073                        const G4ThreeVector &vt2, G4FacetVertexType);
00074     G4TriangularFacet (const G4TriangularFacet &right);
00075 
00076     G4TriangularFacet &operator=(const G4TriangularFacet &right);    
00077 
00078     G4VFacet *GetClone ();
00079     G4TriangularFacet *GetFlippedFacet ();
00080 
00081     G4ThreeVector Distance (const G4ThreeVector &p);
00082     G4double Distance (const G4ThreeVector &p, G4double minDist);
00083     G4double Distance (const G4ThreeVector &p, G4double minDist,
00084                        const G4bool outgoing);
00085     G4double Extent   (const G4ThreeVector axis);
00086     G4bool Intersect  (const G4ThreeVector &p, const G4ThreeVector &v,
00087                        const G4bool outgoing, G4double &distance,
00088                              G4double &distFromSurface, G4ThreeVector &normal);
00089     G4double GetArea ();
00090     G4ThreeVector GetPointOnFace () const;
00091 
00092     G4ThreeVector GetSurfaceNormal () const;
00093     void SetSurfaceNormal (G4ThreeVector normal);
00094 
00095     G4GeometryType GetEntityType () const;
00096 
00097     inline G4bool IsDefined () const;
00098     inline G4int GetNumberOfVertices () const;
00099     inline G4ThreeVector GetVertex (G4int i) const;
00100     inline void SetVertex (G4int i, const G4ThreeVector &val);
00101 
00102     inline G4ThreeVector GetCircumcentre () const;
00103     inline G4double GetRadius () const;
00104 
00105     inline G4int AllocatedMemory();
00106 
00107     inline G4int GetVertexIndex (G4int i) const;
00108     inline void SetVertexIndex (G4int i, G4int j);
00109     inline void SetVertices(std::vector<G4ThreeVector> *v);
00110 
00111   private:
00112 
00113     void CopyFrom(const G4TriangularFacet &rhs);
00114 
00115     G4ThreeVector fSurfaceNormal;
00116     G4double fArea;
00117     G4ThreeVector fCircumcentre;
00118     G4double fRadius;
00119     G4int fIndices[3];
00120 
00121     std::vector<G4ThreeVector> *fVertices;
00122 
00123     G4double fA, fB, fC;
00124     G4double fDet;
00125     G4double fSqrDist;
00126     G4ThreeVector fE1, fE2;
00127     G4bool fIsDefined;
00128 };
00129 
00131 // Inlined Methods
00133 
00134 inline G4bool G4TriangularFacet::IsDefined () const
00135 {
00136   return fIsDefined;
00137 }
00138 
00139 inline G4int G4TriangularFacet::GetNumberOfVertices () const
00140 {
00141   return 3;
00142 }
00143 
00144 inline G4ThreeVector G4TriangularFacet::GetVertex (G4int i) const
00145 {      
00146   G4int indice = fIndices[i];
00147   return indice < 0 ? (*fVertices)[i] : (*fVertices)[indice];
00148 }
00149 
00150 inline void G4TriangularFacet::SetVertex (G4int i, const G4ThreeVector &val)
00151 {
00152   (*fVertices)[i] = val;
00153 }
00154 
00155 inline G4ThreeVector G4TriangularFacet::GetCircumcentre () const
00156 {
00157   return fCircumcentre;
00158 }
00159 
00160 inline G4double G4TriangularFacet::GetRadius () const
00161 {
00162   return fRadius;
00163 }
00164 
00165 inline G4int G4TriangularFacet::AllocatedMemory()
00166 {
00167   G4int size = sizeof(*this);
00168   size += GetNumberOfVertices() * sizeof(G4ThreeVector);
00169   return size;
00170 }
00171 
00172 inline G4int G4TriangularFacet::GetVertexIndex (G4int i) const
00173 {
00174   if (i < 3) return fIndices[i];
00175   else       return 999999999;
00176 }
00177 
00178 inline void G4TriangularFacet::SetVertexIndex (G4int i, G4int j)
00179 {
00180   fIndices[i] = j;
00181 }
00182 
00183 inline void G4TriangularFacet::SetVertices(std::vector<G4ThreeVector> *v)
00184 {
00185   if (fIndices[0] < 0 && fVertices)
00186   {
00187     delete fVertices;
00188     fVertices = 0;
00189   }
00190   fVertices = v;
00191 }
00192 
00193 #endif

Generated on Mon May 27 17:50:03 2013 for Geant4 by  doxygen 1.4.7