VGM Version 5.3
Loading...
Searching...
No Matches
Arb8.cxx
Go to the documentation of this file.
1// $Id$
2
3// -----------------------------------------------------------------------
4// The RootGM 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
12//
13// Class RootGM::Arb8
14// --------------------
15// VGM implementation for Root Arb8 solid.
16//
17// Author: Ivana Hrivnacova; IPN Orsay
18
19#include "RootGM/solids/Arb8.h"
20#include "RootGM/common/Units.h"
22
23#include "TGeoArb8.h"
24#include "TMath.h"
25
26#include <cstdlib>
27#include <iostream>
28#include <math.h>
29
31
32//_____________________________________________________________________________
34 const std::string& name, double hz, std::vector<VGM::TwoVector> vertices)
35 : VGM::ISolid(), VGM::IArb8(), BaseVGM::VArb8(), fArb8(0)
36{
57
58 fArb8 = new TGeoArb8(hz / RootGM::Units::Length());
59 fArb8->SetName(name.c_str());
60
61 // set vertices
62 for (unsigned int i = 0; i < vertices.size(); ++i) {
63 fArb8->SetVertex(i, vertices[i].first / RootGM::Units::Length(),
64 vertices[i].second / RootGM::Units::Length());
65 }
66
68}
69
70//_____________________________________________________________________________
71RootGM::Arb8::Arb8(TGeoArb8* arb8)
72 : VGM::ISolid(), VGM::IArb8(), BaseVGM::VArb8(), fArb8(arb8)
73{
75
77}
78
79//_____________________________________________________________________________
80RootGM::Arb8::Arb8() : VGM::ISolid(), VGM::IArb8(), BaseVGM::VArb8(), fArb8(0)
81{
83}
84
85//_____________________________________________________________________________
87 : VGM::ISolid(rhs), VGM::IArb8(rhs), BaseVGM::VArb8(rhs), fArb8(0)
88{
90}
91
92//_____________________________________________________________________________
94{
95 //
96}
97
98//_____________________________________________________________________________
99std::string RootGM::Arb8::Name() const { return fArb8->GetName(); }
100
101//_____________________________________________________________________________
102int RootGM::Arb8::NofVertices() const { return fgkNofVertices; }
103
104//_____________________________________________________________________________
106{
107 if (index < 0 || index >= fgkNofVertices) {
108 std::cerr << "+++ Error +++" << std::endl;
109 std::cerr << " Wrong vertex index: " << index << std::endl;
110 exit(1);
111 }
112
113 Double_t* xy = fArb8->GetVertices();
114 /*
115 for ( int i=0; i< 2*fgkNofVertices; i++ )
116 std::cout << xy[i] << ", ";
117 std::cout << std::endl;
118 */
119 // CHECK !!!
120 return VGM::TwoVector(xy[2 * index] * RootGM::Units::Length(),
121 xy[2 * index + 1] * RootGM::Units::Length());
122}
123
124//_____________________________________________________________________________
125double RootGM::Arb8::TwistAngle(int index) const
126{
127 if (index < 0 || index >= 4) {
128 std::cerr << "+++ Error +++" << std::endl;
129 std::cerr << " Wrong twist angle index: " << index << std::endl;
130 exit(1);
131 }
132
133 Double_t twistAngle = TMath::ATan(fArb8->GetTwist(index)) * TMath::RadToDeg();
134
135 return twistAngle * RootGM::Units::Angle();
136}
137
138//_____________________________________________________________________________
140{
141 return fArb8->GetDZ() * RootGM::Units::Length();
142}
VGM implementation for Root Arb8 solid.
Definition Arb8.h:33
virtual double ZHalfLength() const
Return the half-length along the z axis in mm.
Definition Arb8.cxx:139
virtual ~Arb8()
Definition Arb8.cxx:93
virtual VGM::TwoVector Vertex(int index) const
Return the index-th vertex.
Definition Arb8.cxx:105
virtual std::string Name() const
Return the name of this solid.
Definition Arb8.cxx:99
virtual int NofVertices() const
Return the number of vertices.
Definition Arb8.cxx:102
virtual double TwistAngle(int index) const
Return the index-th twist angle.
Definition Arb8.cxx:125
static const int fgkNofVertices
Definition Arb8.h:51
TGeoArb8 * fArb8
Definition Arb8.h:53
void AddSolid(VGM::ISolid *, TGeoShape *)
Definition SolidMap.cxx:59
static SolidMap * Instance()
Definition SolidMap.cxx:28
static double Length()
Return Root length unit in VGM units.
Definition Units.h:84
static double Angle()
Return Root angle unit in VGM units.
Definition Units.h:85
BaseVGM utilities.
Definition utilities.h:23
VGM interfaces.
Definition VMedium.h:28
std::pair< double, double > TwoVector
Definition TwoVector.h:28