Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/specific/include/G4QuadrangularFacet.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 and of QinetiQ Ltd,   *
 20 // * subject to DEFCON 705 IPR conditions.                            *
 21 // * By using,  copying,  modifying or  distributing the software (or *
 22 // * any work based  on the software)  you  agree  to acknowledge its *
 23 // * use  in  resulting  scientific  publications,  and indicate your *
 24 // * acceptance of all terms of the Geant4 Software license.          *
 25 // ********************************************************************
 26 //
 27 // G4QuadrangularFacet
 28 //
 29 // Class description:
 30 //
 31 //   The G4QuadrangularFacet class is used for the contruction of
 32 //   G4TessellatedSolid.
 33 //   It is defined by four fVertices, which shall be in the same plane and be
 34 //   supplied in anti-clockwise order looking from the outsider of the solid
 35 //   where it belongs. Its constructor
 36 //   
 37 //     G4QuadrangularFacet (const G4ThreeVector& Pt0, const G4ThreeVector& vt1,
 38 //                          const G4ThreeVector& vt2, const G4ThreeVector& vt3,
 39 //                          G4FacetVertexType);
 40 //
 41 //   takes 5 parameters to define the four fVertices:
 42 //     1) G4FacetvertexType = "ABSOLUTE": in this case Pt0, vt1, vt2 and vt3
 43 //        are the four fVertices required in anti-clockwise order when looking
 44 //        from the outsider.
 45 //     2) G4FacetvertexType = "RELATIVE": in this case the first vertex is Pt0,
 46 //        the second vertex is Pt0+vt, the third vertex is Pt0+vt2 and 
 47 //        the fourth vertex is Pt0+vt3, in anti-clockwise order when looking 
 48 //        from the outsider.
 49 
 50 // 31 October 2004, P R Truscott, QinetiQ Ltd, UK - Created.
 51 // 12 October 2012, M Gayer, CERN, - Reviewed optimized implementation.
 52 // --------------------------------------------------------------------
 53 #ifndef G4QUADRANGULARFACET_HH
 54 #define G4QUADRANGULARFACET_HH
 55 
 56 #include "G4VFacet.hh"
 57 #include "G4Types.hh"
 58 #include "G4ThreeVector.hh"
 59 #include "G4TriangularFacet.hh"
 60 
 61 class G4QuadrangularFacet : public G4VFacet
 62 {
 63   public:
 64 
 65     G4QuadrangularFacet (const G4ThreeVector& Pt0, const G4ThreeVector& vt1,
 66                          const G4ThreeVector& vt2, const G4ThreeVector& vt3,
 67                                G4FacetVertexType);
 68     G4QuadrangularFacet (const G4QuadrangularFacet& right);
 69    ~G4QuadrangularFacet () override;
 70 
 71     G4QuadrangularFacet& operator=(const G4QuadrangularFacet& right);    
 72 
 73     G4VFacet* GetClone () override;
 74 
 75     G4ThreeVector Distance (const G4ThreeVector& p);
 76     G4double Distance (const G4ThreeVector& p, G4double minDist) override;
 77     G4double Distance (const G4ThreeVector& p, G4double minDist,
 78                        const G4bool outgoing) override;
 79     G4double Extent   (const G4ThreeVector axis) override;
 80     G4bool Intersect  (const G4ThreeVector& p, const G4ThreeVector& v,
 81                        const G4bool outgoing, G4double& distance,
 82                              G4double& distFromSurface,
 83                              G4ThreeVector& normal) override;
 84     G4ThreeVector GetSurfaceNormal () const override;
 85 
 86     G4double GetArea () const override;
 87     G4ThreeVector GetPointOnFace () const override;
 88 
 89     G4GeometryType GetEntityType () const override;
 90 
 91     inline G4bool IsDefined () const override;
 92     inline G4int GetNumberOfVertices () const override;
 93     inline G4ThreeVector GetVertex (G4int i) const override;
 94     inline void SetVertex (G4int i, const G4ThreeVector& val) override;
 95     inline void SetVertices(std::vector<G4ThreeVector>* v) override;
 96 
 97     inline G4double GetRadius () const override;
 98     inline G4ThreeVector GetCircumcentre () const override;
 99 
100   private:
101 
102     inline G4int GetVertexIndex (G4int i) const override;
103     inline void SetVertexIndex (G4int i, G4int val) override;
104 
105     inline G4int AllocatedMemory() override;
106 
107   private:
108 
109     G4double fRadius = 0.0;
110 
111     G4ThreeVector fCircumcentre;
112 
113     G4TriangularFacet fFacet1, fFacet2;
114 };
115 
116 // --------------------------------------------------------------------
117 // Inlined Methods
118 // --------------------------------------------------------------------
119 
120 inline G4int G4QuadrangularFacet::GetNumberOfVertices () const
121 {
122   return 4;
123 }
124 
125 inline G4ThreeVector G4QuadrangularFacet::GetVertex (G4int i) const
126 {
127   return i == 3 ? fFacet2.GetVertex(2) : fFacet1.GetVertex(i);
128 }
129 
130 
131 inline G4double G4QuadrangularFacet::GetRadius () const
132 {
133   return fRadius;
134 }
135 
136 inline G4ThreeVector G4QuadrangularFacet::GetCircumcentre () const
137 {
138   return fCircumcentre;
139 }
140 
141 inline void G4QuadrangularFacet::SetVertex (G4int i, const G4ThreeVector &val)
142 {
143   switch (i)
144   {
145     case 0:
146       fFacet1.SetVertex(0, val);
147       fFacet2.SetVertex(0, val);
148       break;
149     case 1:
150       fFacet1.SetVertex(1, val);
151       break;
152     case 2:
153       fFacet1.SetVertex(2, val);
154       fFacet2.SetVertex(1, val);
155       break;
156     case 3:
157       fFacet2.SetVertex(2, val);
158       break;
159     }
160 }
161 
162 inline void G4QuadrangularFacet::SetVertices(std::vector<G4ThreeVector>* v)
163 {
164   fFacet1.SetVertices(v);
165   fFacet2.SetVertices(v);
166 }
167 
168 inline G4bool G4QuadrangularFacet::IsDefined () const
169 {
170   return fFacet1.IsDefined();
171 }
172 
173 inline G4int G4QuadrangularFacet::GetVertexIndex (G4int i) const
174 {
175   return i == 3 ? fFacet2.GetVertexIndex(2) : fFacet1.GetVertexIndex(i);
176 }
177 
178 
179 inline void G4QuadrangularFacet::SetVertexIndex (G4int i, G4int val)
180 {
181   switch (i)
182   {
183     case 0:
184       fFacet1.SetVertexIndex(0, val);
185       fFacet2.SetVertexIndex(0, val);
186       break;
187     case 1:
188       fFacet1.SetVertexIndex(1, val);
189       break;
190     case 2:
191       fFacet1.SetVertexIndex(2, val);
192       fFacet2.SetVertexIndex(1, val);
193       break;
194     case 3:
195       fFacet2.SetVertexIndex(2, val);
196       break;
197   }
198 }
199 
200 inline G4int G4QuadrangularFacet::AllocatedMemory()
201 {
202   return sizeof(*this) + fFacet1.AllocatedMemory() + fFacet2.AllocatedMemory();
203 }
204 
205 #endif
206