32 #include <boost/python.hpp> 33 #include "G4TwistedTubs.hh" 45 = &G4TwistedTubs::GetEndInnerRadius;
47 = &G4TwistedTubs::GetEndInnerRadius;
51 = &G4TwistedTubs::GetEndOuterRadius;
53 = &G4TwistedTubs::GetEndOuterRadius;
58 G4double twistedangle,
64 return new G4TwistedTubs(name, twistedangle, endinnerrad,
65 endouterrad, halfzlen, dphi);
69 G4double twistedangle,
76 return new G4TwistedTubs(name, twistedangle, endinnerrad,
77 endouterrad, halfzlen, nseg, totphi);
81 G4double twistedangle,
84 G4double negativeEndz,
85 G4double positiveEndz,
88 return new G4TwistedTubs(name, twistedangle, innerrad, outerrad,
89 negativeEndz, positiveEndz, dphi);
93 G4double twistedangle,
96 G4double negativeEndz,
97 G4double positiveEndz,
101 return new G4TwistedTubs(name, twistedangle, innerrad, outerrad,
102 negativeEndz, positiveEndz, nseg, totphi);
115 class_<G4TwistedTubs, G4TwistedTubs*, bases<G4VSolid> >
116 (
"G4TwistedTubs",
"twisted tube solid class", no_init)
118 .def(init<
const G4String&, G4double, G4double, G4double,
119 G4double, G4double>())
120 .def(init<
const G4String&, G4double, G4double, G4double,
121 G4double, G4int, G4double>())
122 .def(init<
const G4String&, G4double, G4double, G4double,
123 G4double, G4double, G4double>())
124 .def(init<
const G4String&, G4double, G4double, G4double,
125 G4double, G4double, G4int, G4double>())
127 .def(
"GetDPhi", &G4TwistedTubs::GetDPhi)
128 .def(
"GetPhiTwist", &G4TwistedTubs::GetPhiTwist)
129 .def(
"GetInnerRadius", &G4TwistedTubs::GetInnerRadius)
130 .def(
"GetOuterRadius", &G4TwistedTubs::GetOuterRadius)
131 .def(
"GetInnerStereo", &G4TwistedTubs::GetInnerStereo)
132 .def(
"GetOuterStereo", &G4TwistedTubs::GetOuterStereo)
133 .def(
"GetZHalfLength", &G4TwistedTubs::GetZHalfLength)
134 .def(
"GetKappa", &G4TwistedTubs::GetKappa)
135 .def(
"GetTanInnerStereo", &G4TwistedTubs::GetTanInnerStereo)
136 .def(
"GetTanInnerStereo2", &G4TwistedTubs::GetTanInnerStereo2)
137 .def(
"GetTanOuterStereo", &G4TwistedTubs::GetTanOuterStereo)
138 .def(
"GetTanOuterStereo2", &G4TwistedTubs::GetTanOuterStereo2)
139 .def(
"GetEndZ", &G4TwistedTubs::GetEndZ)
140 .def(
"GetEndPhi", &G4TwistedTubs::GetEndPhi)
146 .def(self_ns::str(
self))
151 return_value_policy<manage_new_object>());
153 return_value_policy<manage_new_object>());
155 return_value_policy<manage_new_object>());
157 return_value_policy<manage_new_object>());
G4double(G4TwistedTubs::* f1_GetEndInnerRadius)(G4int) const
G4double(G4TwistedTubs::* f2_GetEndInnerRadius)(G4int) const
void export_G4TwistedTubs()
G4TwistedTubs * f4_CreateTwistedTubs(const G4String &name, G4double twistedangle, G4double innerrad, G4double outerrad, G4double negativeEndz, G4double positiveEndz, G4int nseg, G4double totphi)
G4double(G4TwistedTubs::* f2_GetEndOuterRadius)() const
G4TwistedTubs * f1_CreateTwistedTubs(const G4String &name, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4double dphi)
G4TwistedTubs * f3_CreateTwistedTubs(const G4String &name, G4double twistedangle, G4double innerrad, G4double outerrad, G4double negativeEndz, G4double positiveEndz, G4double dphi)
G4TwistedTubs * f2_CreateTwistedTubs(const G4String &name, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4int nseg, G4double totphi)
G4double(G4TwistedTubs::* f1_GetEndOuterRadius)(G4int) const