Belle II Software development
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
17namespace Belle2 {
22 namespace ECL {
24 struct Plane_t {
25 G4ThreeVector n;
26 double d;
27 // => n.x*x + n.y*y + n.z*z + d = 0
28 };
29
31 struct Point_t {
32 double x;
33 double y;
34 };
35
37 class BelleCrystal : public G4CSGSolid {
38 public:
40 explicit BelleCrystal(const G4String& pName);
42 BelleCrystal(const G4String& pName, int, const G4ThreeVector*);
43
45 virtual ~BelleCrystal() ;
46
48 Plane_t GetSidePlane(G4int n) const {return fPlanes[n];}
49
50 // Methods for solid
51
52 G4double GetCubicVolume();
53 G4double GetSurfaceArea();
56 void ComputeDimensions(G4VPVParameterisation* p,
57 const G4int n,
58 const G4VPhysicalVolume* pRep);
59
61 G4bool CalculateExtent(const EAxis pAxis,
62 const G4VoxelLimits& pVoxelLimit,
63 const G4AffineTransform& pTransform,
64 G4double& pMin, G4double& pMax) const;
65
67 EInside Inside(const G4ThreeVector& p) const;
68
70 G4ThreeVector SurfaceNormal(const G4ThreeVector& p) const;
71
73 G4double DistanceToIn(const G4ThreeVector& p, const G4ThreeVector& v) const;
74
76 G4double DistanceToIn(const G4ThreeVector& p) const;
77
79 G4double DistanceToOut(const G4ThreeVector& p, const G4ThreeVector& v,
80 const G4bool calcNorm = false,
81 G4bool* validNorm = 0, G4ThreeVector* n = 0) const;
82
84 G4double DistanceToOut(const G4ThreeVector& p) const;
85
87 G4GeometryType GetEntityType() const;
88
90 G4ThreeVector GetPointOnSurface() const;
91
93 G4VSolid* Clone() const;
94
96 std::ostream& StreamInfo(std::ostream& os) const;
97
99 void DescribeYourselfTo(G4VGraphicsScene& scene) const;
101 G4Polyhedron* CreatePolyhedron() const;
102
104 void BoundingLimits(G4ThreeVector& pMin, G4ThreeVector& pMax) const;
105
106 public: // without description
107
113 explicit BelleCrystal(__void__&);
114
116 BelleCrystal(const BelleCrystal& rhs);
120 G4ThreeVector vertex(unsigned int i) const;
121
122 protected: // with description
123
132 G4bool MakePlane(const G4ThreeVector& p1,
133 const G4ThreeVector& p2,
134 const G4ThreeVector& p3,
135 const G4ThreeVector& p4,
136 Plane_t& plane) const;
137
138 private:
139 G4ThreeVector GetPointOnTriangle(int) const;
140 double area(int, double&) const;
141 double getvolarea() const;
142 const unsigned int* ivertx(unsigned int i) const;
143 private:
144 unsigned int nsides;
145 double fDz;
146 std::vector<Plane_t> fPlanes;
147 std::vector<Point_t> fx;
149 mutable std::vector<double> fareas;
150 };
151
153 class PolyhedronBelleCrystal: public G4Polyhedron {
154 public:
155 PolyhedronBelleCrystal(int, const G4ThreeVector*);
156 virtual ~PolyhedronBelleCrystal();
157 };
158
159 }
161}
a Belle crystal in Geant4
Definition: BelleCrystal.h:37
std::vector< Plane_t > fPlanes
vector of planes
Definition: BelleCrystal.h:146
G4GeometryType GetEntityType() const
Get entity type.
void DescribeYourselfTo(G4VGraphicsScene &scene) const
Visualisation functions.
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.
G4Polyhedron * CreatePolyhedron() const
create polyhedron
G4bool MakePlane(const G4ThreeVector &p1, const G4ThreeVector &p2, const G4ThreeVector &p3, const G4ThreeVector &p4, Plane_t &plane) const
Calculate the coef's of the plane p1->p2->p3->p4->p1 where the ThreeVectors 1-4 are in anti-clockwise...
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
Plane_t GetSidePlane(G4int n) const
Get side plane.
Definition: BelleCrystal.h:48
G4ThreeVector GetPointOnSurface() const
get point on surface
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
calculate the extent of the volume
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
Calculate side nearest to p, and return normal.
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Calculate exact shortest distance to any boundary from outside.
virtual ~BelleCrystal()
Destructor.
G4VSolid * Clone() const
Make a clone of the object.
std::ostream & StreamInfo(std::ostream &os) const
Stream object contents to an output stream.
std::vector< double > fareas
vector of area values
Definition: BelleCrystal.h:149
double getvolarea() const
get volume area
EInside Inside(const G4ThreeVector &p) const
Return whether point inside/outside/on surface, using tolerance.
std::vector< Point_t > fx
vector of points
Definition: BelleCrystal.h:147
G4double GetSurfaceArea()
get the surface area
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
compute the dimensions
unsigned int nsides
the number of sides
Definition: BelleCrystal.h:144
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const
Calculate exact shortest distance to any boundary from inside.
Belle crystal in polyhedron.
Definition: BelleCrystal.h:153
virtual ~PolyhedronBelleCrystal()
destructor
Abstract base class for different kinds of events.
struct for plane
Definition: BelleCrystal.h:24
G4ThreeVector n
Normal unit vector (a,b,c)
Definition: BelleCrystal.h:25
double d
offset (d)
Definition: BelleCrystal.h:26
struct for Point
Definition: BelleCrystal.h:31
double y
y coordinate
Definition: BelleCrystal.h:33
double x
x coordinate
Definition: BelleCrystal.h:32