Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TwistTrapFlatSide.hh
Go to the documentation of this file.
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// G4FlatTrapSurface
27//
28// Class description:
29//
30// Class describing a flat boundary surface for a trapezoid.
31
32// Author: 27-Oct-2004 - O.Link ([email protected])
33// --------------------------------------------------------------------
34#ifndef G4TWISTTRAPFLATSIDE_HH
35#define G4TWISTTRAPFLATSIDE_HH
36
37#include "G4VTwistSurface.hh"
38
40{
41 public: // with description
42
43 G4TwistTrapFlatSide( const G4String& name,
44 G4double PhiTwist,
45 G4double pDx1,
46 G4double pDx2,
47 G4double pDy,
48 G4double pDz,
49 G4double pAlpha,
50 G4double pPhi,
51 G4double pTheta,
52 G4int handedness );
53 virtual ~G4TwistTrapFlatSide();
54
55 virtual G4ThreeVector GetNormal(const G4ThreeVector& /* xx */ ,
56 G4bool isGlobal = false);
57 virtual G4int DistanceToSurface(const G4ThreeVector& gp,
58 const G4ThreeVector& gv,
59 G4ThreeVector gxx[],
60 G4double distance[],
61 G4int areacode[],
62 G4bool isvalid[],
63 EValidate validate = kValidateWithTol);
64
65 virtual G4int DistanceToSurface(const G4ThreeVector& gp,
66 G4ThreeVector gxx[],
67 G4double distance[],
68 G4int areacode[]);
69
70
72 G4bool isGlobal = false);
75 virtual G4double GetSurfaceArea();
76 virtual void GetFacets( G4int m, G4int n, G4double xyz[][3],
77 G4int faces[][4], G4int iside );
78
79 public: // without description
80
81 G4TwistTrapFlatSide(__void__&);
82 // Fake default constructor for usage restricted to direct object
83 // persistency for clients requiring preallocation of memory for
84 // persistifiable objects.
85
86 protected: // with description
87
88 virtual G4int GetAreaCode(const G4ThreeVector& xx,
89 G4bool withTol = true);
90
91 private:
92
93 virtual void SetCorners();
94 virtual void SetBoundaries();
95
96 inline double xAxisMax(G4double u, G4double fTanAlpha) const;
97
98 private:
99
100 G4double fDx1;
101 G4double fDx2;
102 G4double fDy;
103 G4double fDz;
104 G4double fPhiTwist;
105 G4double fAlpha;
106 G4double fTAlph;
107 G4double fPhi;
108 G4double fTheta;
109 G4double fdeltaX;
110 G4double fdeltaY;
111};
112
113//========================================================
114// inline functions
115//========================================================
116
117inline
118G4double G4TwistTrapFlatSide::xAxisMax(G4double u, G4double fTanAlpha) const
119{
120 return ( ( fDx2 + fDx1 )/2. + u*(fDx2 - fDx1)/(2.*fDy) + u *fTanAlpha ) ;
121}
122
123inline G4ThreeVector
125{
126 G4ThreeVector SurfPoint ( x,y,0);
127
128 if (isGlobal) { return (fRot*SurfPoint + fTrans); }
129 return SurfPoint;
130}
131
132inline
134{
135 return -xAxisMax(y, -fTAlph) ;
136}
137
138inline
140{
141 return xAxisMax(y, fTAlph) ;
142}
143
144inline
146{
147 return 2*(fDx1 + fDx2)*fDy ;
148}
149
150#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual G4int GetAreaCode(const G4ThreeVector &xx, G4bool withTol=true)
virtual void GetFacets(G4int m, G4int n, G4double xyz[][3], G4int faces[][4], G4int iside)
virtual G4ThreeVector GetNormal(const G4ThreeVector &, G4bool isGlobal=false)
virtual G4double GetSurfaceArea()
virtual G4double GetBoundaryMin(G4double u)
virtual G4double GetBoundaryMax(G4double u)
virtual G4ThreeVector SurfacePoint(G4double x, G4double y, G4bool isGlobal=false)
virtual G4int DistanceToSurface(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector gxx[], G4double distance[], G4int areacode[], G4bool isvalid[], EValidate validate=kValidateWithTol)
G4RotationMatrix fRot
G4ThreeVector fTrans