VGM Version 5.3
Loading...
Searching...
No Matches
Arb8.h
Go to the documentation of this file.
1// $Id$
2
3// -----------------------------------------------------------------------
4// The Geant4GM package of the Virtual Geometry Model
5// Copyright (C) 2007, Ivana Hrivnacova
6// All rights reserved.
7//
8// For the licensing terms see vgm/LICENSE.
9// Contact: ivana@ipno.in2p3.fr
10// -----------------------------------------------------------------------
11
13//
21
22#ifndef GEANT4_GM_ARB8_SOLID_H
23#define GEANT4_GM_ARB8_SOLID_H
24
26
28
29#include "G4ThreeVector.hh"
30#include "globals.hh"
31
32#include <string>
33#include <vector>
34
35class G4TessellatedSolid;
36class G4VFacet;
37
38namespace Geant4GM {
39
40class Arb8 : public BaseVGM::VArb8
41{
42 public:
43 Arb8(
44 const std::string& name, double hz, std::vector<VGM::TwoVector> vertices);
45 virtual ~Arb8();
46
47 // static methods
48 static bool IsTwisted(std::vector<VGM::TwoVector> vertices);
49
50 // methods
51 virtual std::string Name() const;
52 virtual int NofVertices() const;
53 virtual VGM::TwoVector Vertex(int index) const;
54 virtual double TwistAngle(int index) const;
55 virtual double ZHalfLength() const;
56
57 protected:
58 Arb8();
59 Arb8(const Arb8& rhs);
60
61 private:
62 // methods
63 void ReorderVertices(std::vector<G4ThreeVector>& vertices);
64
65 G4VFacet* MakeDownFacet(std::vector<G4ThreeVector> fromVertices, int ind1,
66 int ind2, int ind3) const;
67 G4VFacet* MakeUpFacet(std::vector<G4ThreeVector> fromVertices, int ind1,
68 int ind2, int ind3) const;
69 G4VFacet* MakeSideFacet(G4ThreeVector downVertex0, G4ThreeVector downVertex1,
70 G4ThreeVector upVertex1, G4ThreeVector upVertex0) const;
71
72 // static data members
73 static const int fgkNofVertices;
74 static const double fgkTolerance;
75
76 // data members
77 double fHz;
78 std::vector<VGM::TwoVector> fVertices;
79 G4TessellatedSolid* fTessellatedSolid;
80};
81
82} // namespace Geant4GM
83
84#endif // GEANT4_GM_EXTRUDED_SOLID_H
The ABC for Arb8 solids.
Definition VArb8.h:30
VGM implementation for Geant4 Arb8 solid, the shape is implemented using G4TessellatedSolid.
Definition Arb8.h:41
static bool IsTwisted(std::vector< VGM::TwoVector > vertices)
Definition Arb8.cxx:34
virtual ~Arb8()
Definition Arb8.cxx:178
virtual int NofVertices() const
Return the number of vertices.
Definition Arb8.cxx:292
virtual double ZHalfLength() const
Return the half-length along the z axis in mm.
Definition Arb8.cxx:321
virtual VGM::TwoVector Vertex(int index) const
Return the index-th vertex.
Definition Arb8.cxx:295
virtual double TwistAngle(int index) const
Return the index-th twist angle.
Definition Arb8.cxx:307
virtual std::string Name() const
Return the name of this solid.
Definition Arb8.cxx:286
VGM implementation for Geant4.
Definition Element.h:29
std::pair< double, double > TwoVector
Definition TwoVector.h:28