Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/geometry/solids/specific/include/G4TwistBoxSide.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 // G4TwistBoxSide
 27 //
 28 // Class description:
 29 //
 30 // G4TwistBoxSide describes a twisted boundary surface for a trapezoid.
 31 
 32 // Author: 27-Oct-2004 - O.Link (Oliver.Link@cern.ch)
 33 // --------------------------------------------------------------------
 34 #ifndef G4TWISTBOXSIDE_HH
 35 #define G4TWISTBOXSIDE_HH
 36 
 37 #include "G4VTwistSurface.hh"
 38 
 39 #include <vector>
 40 
 41 class G4TwistBoxSide : public G4VTwistSurface
 42 {
 43   public:
 44    
 45     G4TwistBoxSide(const G4String& name,
 46                          G4double  PhiTwist, // twist angle
 47                          G4double  pDz,      // half z lenght
 48                          G4double  pTheta,   // direction between end planes
 49                          G4double  pPhi,     // by polar and azimutal angles
 50                          G4double  pDy1,     // half y length at -pDz
 51                          G4double  pDx1,     // half x length at -pDz,-pDy
 52                          G4double  pDx2,     // half x length at -pDz,+pDy
 53                          G4double  pDy2,     // half y length at +pDz
 54                          G4double  pDx3,     // half x length at +pDz,-pDy
 55                          G4double  pDx4,     // half x length at +pDz,+pDy
 56                          G4double  pAlph,    // tilt angle at +pDz
 57                          G4double  AngleSide // parity
 58                    );
 59   
 60     ~G4TwistBoxSide() override;
 61    
 62     G4ThreeVector  GetNormal(const G4ThreeVector& xx,
 63                                    G4bool isGlobal = false) override ;   
 64    
 65     G4int DistanceToSurface(const G4ThreeVector& gp,
 66                             const G4ThreeVector& gv,
 67                                   G4ThreeVector  gxx[],
 68                                   G4double  distance[],
 69                                   G4int     areacode[],
 70                                   G4bool    isvalid[],
 71                             EValidate validate = kValidateWithTol) override;
 72                                                   
 73     G4int DistanceToSurface(const G4ThreeVector& gp,
 74                                   G4ThreeVector  gxx[],
 75                                   G4double       distance[],
 76                                   G4int          areacode[]) override;
 77 
 78     G4TwistBoxSide(__void__&);
 79       // Fake default constructor for usage restricted to direct object
 80       // persistency for clients requiring preallocation of memory for
 81       // persistifiable objects.
 82 
 83   private:
 84 
 85     G4int GetAreaCode(const G4ThreeVector& xx, 
 86                             G4bool withTol = true) override;
 87     void SetCorners() override;
 88     void SetBoundaries() override;
 89 
 90     void GetPhiUAtX(const G4ThreeVector& p, G4double& phi, G4double& u);
 91     G4ThreeVector ProjectPoint(const G4ThreeVector& p,
 92                                      G4bool isglobal = false);
 93 
 94     inline G4ThreeVector SurfacePoint(G4double phi, G4double u,
 95                                       G4bool isGlobal = false) override;
 96     inline G4double GetBoundaryMin(G4double phi) override;
 97     inline G4double GetBoundaryMax(G4double phi) override;
 98     inline G4double GetSurfaceArea() override;
 99     void GetFacets( G4int m, G4int n, G4double xyz[][3],
100                     G4int faces[][4], G4int iside ) override;
101 
102     inline G4double GetValueA(G4double phi);
103     inline G4double GetValueB(G4double phi);
104     inline G4ThreeVector NormAng(G4double phi, G4double u);
105     inline G4double Xcoef(G4double u, G4double phi);
106       // To calculate the w(u) function
107 
108   private:
109 
110     G4double fTheta;   
111     G4double fPhi ;
112 
113     G4double fDy1;   
114     G4double fDx1;     
115     G4double fDx2;     
116 
117     G4double fDy2;   
118     G4double fDx3;     
119     G4double fDx4;     
120 
121     G4double fDz;         // Half-length along the z axis
122 
123     G4double fAlph;
124     G4double fTAlph;      // std::tan(fAlph)
125     
126     G4double fPhiTwist;   // twist angle ( dphi in surface equation)
127 
128     G4double fAngleSide;
129 
130     G4double fdeltaX;
131     G4double fdeltaY;
132 
133     G4double fDx4plus2;   // fDx4 + fDx2  == a2/2 + a1/2
134     G4double fDx4minus2;  // fDx4 - fDx2          -
135     G4double fDx3plus1;   // fDx3 + fDx1  == d2/2 + d1/2
136     G4double fDx3minus1;  // fDx3 - fDx1          -
137     G4double fDy2plus1;   // fDy2 + fDy1  == b2/2 + b1/2
138     G4double fDy2minus1;  // fDy2 - fDy1          -
139     G4double fa1md1;      // 2 fDx2 - 2 fDx1  == a1 - d1
140     G4double fa2md2;      // 2 fDx4 - 2 fDx3
141 };   
142 
143 //========================================================
144 // inline functions
145 //========================================================
146 
147 inline
148 G4double G4TwistBoxSide::GetValueA(G4double phi)
149 {
150   return ( fDx4plus2 + fDx4minus2 * ( 2 * phi ) / fPhiTwist  ) ;
151 }
152 
153 
154 inline 
155 G4double G4TwistBoxSide::GetValueB(G4double phi) 
156 {
157   return ( fDy2plus1 + fDy2minus1 * ( 2 * phi ) / fPhiTwist ) ;
158 }
159 
160 inline
161 G4double G4TwistBoxSide::Xcoef(G4double u, G4double phi)
162 {
163   
164   return GetValueA(phi)/2. + u*fTAlph    ;
165 
166 }
167 
168 inline G4ThreeVector
169 G4TwistBoxSide::SurfacePoint( G4double phi, G4double u, G4bool isGlobal ) 
170 {
171   // function to calculate a point on the surface, given by parameters phi,u
172 
173   G4ThreeVector SurfPoint ( Xcoef(u,phi) * std::cos(phi)
174                           - u * std::sin(phi) + fdeltaX*phi/fPhiTwist,
175                             Xcoef(u,phi) * std::sin(phi)
176                           + u * std::cos(phi) + fdeltaY*phi/fPhiTwist,
177                             2*fDz*phi/fPhiTwist  );
178 
179   if (isGlobal) { return (fRot * SurfPoint + fTrans); }
180   return SurfPoint;
181 }
182 
183 inline
184 G4double G4TwistBoxSide::GetBoundaryMin(G4double phi)
185 {
186   return -0.5*GetValueB(phi) ;
187 }
188 
189 inline
190 G4double G4TwistBoxSide::GetBoundaryMax(G4double phi)
191 {
192   return 0.5*GetValueB(phi) ;
193 }
194 
195 inline
196 G4double G4TwistBoxSide::GetSurfaceArea()
197 {
198   return (fDz*(std::sqrt(16*fDy1*fDy1
199          + (fa1md1 + 4*fDy1*fTAlph)*(fa1md1 + 4*fDy1*fTAlph))
200          + std::sqrt(16*fDy1*fDy1 + (fa2md2 + 4*fDy1*fTAlph)
201          * (fa2md2 + 4*fDy1*fTAlph))))/2. ;
202 }
203 
204 inline
205 G4ThreeVector G4TwistBoxSide::NormAng( G4double phi, G4double u ) 
206 {
207   // function to calculate the norm at a given point on the surface
208   // replace a1-d1
209 
210   G4ThreeVector nvec( 4*fDz*(std::cos(phi) + fTAlph*std::sin(phi)) ,
211                       4*fDz*(-(fTAlph*std::cos(phi)) + std::sin(phi)),
212                       (fDx2 + fDx4)*fPhiTwist*fTAlph
213                                    + 2*fDx4minus2*(-1 + fTAlph*phi)
214                                    + 2*fPhiTwist*(1 + fTAlph*fTAlph)*u
215                                 - 2*(fdeltaX - fdeltaY*fTAlph)*std::cos(phi)
216                                 - 2*(fdeltaY + fdeltaX*fTAlph)*std::sin(phi) );
217   return nvec.unit();
218 }
219 
220 #endif
221