Belle II Software  release-06-01-15
BelleCrystal.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include "G4Types.hh"
12 
13 #include "G4CSGSolid.hh"
14 #include "G4Polyhedron.hh"
15 #include <vector>
16 
17 namespace Belle2 {
22  namespace ECL {
23  struct Plane_t {
24  G4ThreeVector n;// Normal unit vector (a,b,c)
25  double d; // offset (d)
26  // => n.x*x + n.y*y + n.z*z + d = 0
27  };
28 
29  struct Point_t {
30  double x, y;
31  };
32 
34  class BelleCrystal : public G4CSGSolid {
35  public:
37  explicit BelleCrystal(const G4String& pName);
39  BelleCrystal(const G4String& pName, int, const G4ThreeVector*);
40 
42  virtual ~BelleCrystal() ;
43 
44  Plane_t GetSidePlane(G4int n) const {return fPlanes[n];}
45 
46  // Methods for solid
47 
48  G4double GetCubicVolume();
49  G4double GetSurfaceArea();
52  void ComputeDimensions(G4VPVParameterisation* p,
53  const G4int n,
54  const G4VPhysicalVolume* pRep);
55 
57  G4bool CalculateExtent(const EAxis pAxis,
58  const G4VoxelLimits& pVoxelLimit,
59  const G4AffineTransform& pTransform,
60  G4double& pMin, G4double& pMax) const;
61 
62  EInside Inside(const G4ThreeVector& p) const;
63 
64  G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
65 
66  G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
67 
68  G4double DistanceToIn(const G4ThreeVector& p) const;
69 
70  G4double DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v,
71  const G4bool calcNorm = false,
72  G4bool* validNorm = 0, G4ThreeVector* n = 0) const;
73 
74  G4double DistanceToOut(const G4ThreeVector& p) const;
75 
76  G4GeometryType GetEntityType() const;
77 
78  G4ThreeVector GetPointOnSurface() const;
79 
80  G4VSolid* Clone() const;
81 
82  std::ostream& StreamInfo(std::ostream& os) const;
83 
84  // Visualisation functions
85 
86  void DescribeYourselfTo(G4VGraphicsScene& scene) const;
87  G4Polyhedron* CreatePolyhedron() const;
88 
90  void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
91 
92  public: // without description
93 
99  explicit BelleCrystal(__void__&);
100 
102  BelleCrystal(const BelleCrystal& rhs);
104  BelleCrystal& operator=(const BelleCrystal& rhs);
105 
106  G4ThreeVector vertex(unsigned int i) const;
107 
108  protected: // with description
109 
110  G4bool MakePlane(const G4ThreeVector& p1,
111  const G4ThreeVector& p2,
112  const G4ThreeVector& p3,
113  const G4ThreeVector& p4,
114  Plane_t& plane) const;
115 
116  private:
117  G4ThreeVector GetPointOnTriangle(int) const;
118  double area(int, double&) const;
119  double getvolarea() const;
120  const unsigned int* ivertx(unsigned int i) const;
121  private:
122  unsigned int nsides;
123  double fDz;
124  std::vector<Plane_t> fPlanes;
125  std::vector<Point_t> fx;
126 
127  mutable std::vector<double> fareas;
128  };
129 
131  class PolyhedronBelleCrystal: public G4Polyhedron {
132  public:
133  PolyhedronBelleCrystal(int, const G4ThreeVector*);
134  virtual ~PolyhedronBelleCrystal();
135  };
136 
137  }
139 }
a Belle crystal in Geant4
Definition: BelleCrystal.h:34
const unsigned int * ivertx(unsigned int i) const
get the ith vertex
G4ThreeVector GetPointOnTriangle(int) const
Returns a random point on the surface of one of the faces.
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
Two vectors define an axis-parallel bounding box for the shape.
double area(int, double &) const
triangle area
G4double GetCubicVolume()
get the Cubic volume
BelleCrystal & operator=(const BelleCrystal &rhs)
assignment operator
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
calculate the extent of the volume
virtual ~BelleCrystal()
Destructor.
G4double GetSurfaceArea()
get the surface area
BelleCrystal(const G4String &pName)
Constructor for "nominal" BelleCrystal.
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
compute the dimensions
unsigned int nsides
the number of sides
Definition: BelleCrystal.h:122
Belle crystal in polyhedron.
Definition: BelleCrystal.h:131
PolyhedronBelleCrystal(int, const G4ThreeVector *)
constructor
virtual ~PolyhedronBelleCrystal()
destructor
Abstract base class for different kinds of events.