Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/specific/include/G4VTwistedFaceted.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.                      *
 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 // G4VTwistedFaceted
 27 //
 28 // Class description:
 29 //
 30 //  G4VTwistedFaceted is an abstract base class for twisted boxoids:
 31 //  G4TwistedTrd, G4TwistedTrap and G4TwistedBox
 32 
 33 // Author: 27-Oct-2004 - O.Link (CERN)
 34 // --------------------------------------------------------------------
 35 #ifndef G4VTWISTEDFACETED_HH
 36 #define G4VTWISTEDFACETED_HH 1
 37 
 38 #include "G4VSolid.hh"
 39 #include "G4TwoVector.hh"
 40 #include "G4TwistTrapAlphaSide.hh"
 41 #include "G4TwistTrapParallelSide.hh"
 42 #include "G4TwistBoxSide.hh"
 43 #include "G4TwistTrapFlatSide.hh"
 44 
 45 class G4SolidExtentList;
 46 class G4ClippablePolygon;
 47 
 48 class G4VTwistedFaceted: public G4VSolid
 49 {
 50   public:
 51 
 52     G4VTwistedFaceted(const G4String& pname,    // Name of instance
 53                             G4double PhiTwist,  // twist angle
 54                             G4double pDz,       // half z lenght
 55                             G4double pTheta,  // direction between end planes
 56                             G4double pPhi,    // defined by polar & azim. angles
 57                             G4double pDy1,    // half y length at -pDz
 58                             G4double pDx1,    // half x length at -pDz,-pDy
 59                             G4double pDx2,    // half x length at -pDz,+pDy
 60                             G4double pDy2,    // half y length at +pDz
 61                             G4double pDx3,    // half x length at +pDz,-pDy
 62                             G4double pDx4,    // half x length at +pDz,+pDy
 63                             G4double pAlph    // tilt angle at +pDz
 64                      );
 65 
 66     ~G4VTwistedFaceted() override;
 67 
 68     void ComputeDimensions(      G4VPVParameterisation*,
 69                            const G4int,
 70                            const G4VPhysicalVolume*  ) override;
 71 
 72     void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const override;
 73 
 74     G4bool CalculateExtent(const EAxis              pAxis,
 75                            const G4VoxelLimits&     pVoxelLimit,
 76                            const G4AffineTransform& pTransform,
 77                                  G4double&          pMin,
 78                                  G4double&          pMax ) const override;
 79 
 80     G4double DistanceToIn (const G4ThreeVector& p,
 81                            const G4ThreeVector& v ) const override;
 82 
 83     G4double DistanceToIn (const G4ThreeVector& p ) const override;
 84 
 85     G4double DistanceToOut(const G4ThreeVector& p,
 86                            const G4ThreeVector& v,
 87                            const G4bool         calcnorm  = false,
 88                                  G4bool* validnorm = nullptr,
 89                                  G4ThreeVector*  n = nullptr ) const override;
 90 
 91     G4double DistanceToOut(const G4ThreeVector& p) const override;
 92 
 93     EInside Inside (const G4ThreeVector& p) const override;
 94 
 95     G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const override;
 96 
 97     G4ThreeVector GetPointOnSurface() const override;
 98     G4ThreeVector GetPointInSolid(G4double z) const;
 99 
100     G4double GetCubicVolume() override;
101     G4double GetSurfaceArea() override;
102 
103     void DescribeYourselfTo (G4VGraphicsScene& scene) const override;
104     G4Polyhedron* CreatePolyhedron   () const override;
105     G4Polyhedron* GetPolyhedron      () const override;
106 
107     std::ostream &StreamInfo(std::ostream& os) const override;
108 
109     // accessors
110 
111     inline G4double GetTwistAngle() const { return fPhiTwist; }
112 
113     inline G4double GetDx1   () const { return fDx1   ; }
114     inline G4double GetDx2   () const { return fDx2   ; }
115     inline G4double GetDx3   () const { return fDx3   ; }
116     inline G4double GetDx4   () const { return fDx4   ; }
117     inline G4double GetDy1   () const { return fDy1   ; }
118     inline G4double GetDy2   () const { return fDy2   ; }
119     inline G4double GetDz    () const { return fDz    ; }
120     inline G4double GetPhi   () const { return fPhi   ; }
121     inline G4double GetTheta () const { return fTheta ; }
122     inline G4double GetAlpha () const { return fAlph  ; }
123 
124     inline G4double Xcoef(G4double u, G4double phi, G4double ftg) const;
125       // For calculating the w(u) function
126 
127     inline G4double GetValueA(G4double phi) const;
128     inline G4double GetValueB(G4double phi) const;
129     inline G4double GetValueD(G4double phi) const;
130 
131     G4VisExtent     GetExtent    () const override;
132     G4GeometryType  GetEntityType() const override;
133 
134     G4VTwistedFaceted(__void__&);
135       // Fake default constructor for usage restricted to direct object
136       // persistency for clients requiring preallocation of memory for
137       // persistifiable objects.
138 
139     G4VTwistedFaceted(const G4VTwistedFaceted& rhs);
140     G4VTwistedFaceted& operator=(const G4VTwistedFaceted& rhs);
141       // Copy constructor and assignment operator.
142 
143   protected:
144 
145     mutable G4bool fRebuildPolyhedron = false;
146     mutable G4Polyhedron* fpPolyhedron = nullptr;  // polyhedron for vis
147 
148     G4double fCubicVolume = 0.0; // volume of the solid
149     G4double fSurfaceArea = 0.0; // area of the solid
150 
151   private:
152 
153     double GetLateralFaceArea(const G4TwoVector& p1,
154                               const G4TwoVector& p2,
155                               const G4TwoVector& p3,
156                               const G4TwoVector& p4) const;
157     void CreateSurfaces();
158 
159   private:
160 
161     G4double fTheta;
162     G4double fPhi ;
163 
164     G4double fDy1;
165     G4double fDx1;
166     G4double fDx2;
167 
168     G4double fDy2;
169     G4double fDx3;
170     G4double fDx4;
171 
172     G4double fDz;        // Half-length along the z axis
173 
174     G4double fDx ;       // maximum side in x
175     G4double fDy ;       // maximum side in y
176 
177     G4double fAlph ;
178     G4double fTAlph ;    // std::tan(fAlph)
179 
180     G4double fdeltaX ;
181     G4double fdeltaY ;
182 
183     G4double fPhiTwist;  // twist angle ( dphi in surface equation)
184 
185     G4VTwistSurface* fLowerEndcap ;  // surface of -ve z
186     G4VTwistSurface* fUpperEndcap ;  // surface of +ve z
187 
188     G4VTwistSurface* fSide0 ;    // Twisted Side at phi = 0 deg
189     G4VTwistSurface* fSide90 ;   // Twisted Side at phi = 90 deg
190     G4VTwistSurface* fSide180 ;  // Twisted Side at phi = 180 deg
191     G4VTwistSurface* fSide270 ;  // Twisted Side at phi = 270 deg
192 
193 };
194 
195 //=====================================================================
196 
197 inline
198 G4double G4VTwistedFaceted::GetValueA(G4double phi) const
199 {
200   return ( fDx4 + fDx2  + ( fDx4 - fDx2 ) * ( 2 * phi ) / fPhiTwist  ) ;
201 }
202 
203 inline
204 G4double G4VTwistedFaceted::GetValueD(G4double phi) const
205 {
206   return ( fDx3 + fDx1  + ( fDx3 - fDx1 ) * ( 2 * phi ) / fPhiTwist  ) ;
207 }
208 
209 inline
210 G4double G4VTwistedFaceted::GetValueB(G4double phi) const
211 {
212   return ( fDy2 + fDy1  + ( fDy2 - fDy1 ) * ( 2 * phi ) / fPhiTwist ) ;
213 }
214 
215 inline
216 G4double G4VTwistedFaceted::Xcoef(G4double u, G4double phi, G4double ftg) const
217 {
218   return GetValueA(phi)/2. + (GetValueD(phi)-GetValueA(phi))/4.
219     - u*( ( GetValueD(phi)-GetValueA(phi) ) / ( 2 * GetValueB(phi) ) - ftg );
220 }
221 
222 #endif
223