VGM Version 5.3
Loading...
Searching...
No Matches
Ctubs.cxx
Go to the documentation of this file.
1// $Id$
2
3// -----------------------------------------------------------------------
4// The Geant4GM package of the Virtual Geometry Model
5// Copyright (C) 2014, 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 Ctubs
14// ---------------
15// VGM implementation for cut tubs solid in Geant4.
16//
17// Author: Ivana Hrivnacova; IPN Orsay
18
19#include "ClhepVGM/Units.h"
20
23
24#include "G4CutTubs.hh"
25#include "G4ReflectedSolid.hh"
26#include "geomdefs.hh"
27
28#include <cmath>
29
30using CLHEP::Hep3Vector;
31using CLHEP::HepRotation;
32
33const double Geant4GM::Ctubs::fgkTolerance = 1e-9;
34//_____________________________________________________________________________
35Geant4GM::Ctubs::Ctubs(const std::string& name, double rin, double rout,
36 double hz, double sphi, double dphi, double nxlow, double nylow, double nzlow,
37 double nxhigh, double nyhigh, double nzhigh)
38 : VGM::ISolid(),
39 VGM::ICtubs(),
40 BaseVGM::VCtubs(),
41 fIsReflected(false),
42 fCutTubs(0)
43{
56
58 fCutTubs = new G4CutTubs(name, rin / ClhepVGM::Units::Length(),
61 G4ThreeVector(nxlow / ClhepVGM::Units::Length(),
63 G4ThreeVector(nxhigh / ClhepVGM::Units::Length(),
64 nyhigh / ClhepVGM::Units::Length(), nzhigh / ClhepVGM::Units::Length()));
65
66 Geant4GM::SolidMap::Instance()->AddSolid(this, fCutTubs);
67}
68
69//_____________________________________________________________________________
70Geant4GM::Ctubs::Ctubs(G4CutTubs* ctubs, G4ReflectedSolid* reflCtubs)
71 : VGM::ISolid(),
72 VGM::ICtubs(),
73 BaseVGM::VCtubs(),
74 fIsReflected(false),
75 fCutTubs(ctubs)
76{
78
79 if (reflCtubs) {
80 fIsReflected = true;
81 Geant4GM::SolidMap::Instance()->AddSolid(this, reflCtubs);
82 }
83 else
85}
86
87//_____________________________________________________________________________
88Geant4GM::Ctubs::Ctubs() : VGM::ISolid(), VGM::ICtubs(), BaseVGM::VCtubs()
89{
91}
92
93//_____________________________________________________________________________
95 : VGM::ISolid(rhs), VGM::ICtubs(rhs), BaseVGM::VCtubs(rhs)
96{
98}
99
100//_____________________________________________________________________________
102{
103 //
104}
105
106//_____________________________________________________________________________
107std::string Geant4GM::Ctubs::Name() const { return fCutTubs->GetName(); }
108
109//_____________________________________________________________________________
111{
112 return fCutTubs->GetInnerRadius() * ClhepVGM::Units::Length();
113}
114
115//_____________________________________________________________________________
117{
118 return fCutTubs->GetOuterRadius() * ClhepVGM::Units::Length();
119}
120
121//_____________________________________________________________________________
123{
124 return fCutTubs->GetZHalfLength() * ClhepVGM::Units::Length();
125}
126
127//_____________________________________________________________________________
129{
130 return fCutTubs->GetStartPhiAngle() * ClhepVGM::Units::Angle();
131}
132
133//_____________________________________________________________________________
135{
136 return fCutTubs->GetDeltaPhiAngle() * ClhepVGM::Units::Angle();
137}
138//_____________________________________________________________________________
140{
141 if (!fIsReflected)
142 return fCutTubs->GetLowNorm().x() * ClhepVGM::Units::Length();
143 else
144 return fCutTubs->GetHighNorm().x() * ClhepVGM::Units::Length();
145}
146
147//_____________________________________________________________________________
149{
150 if (!fIsReflected)
151 return fCutTubs->GetLowNorm().y() * ClhepVGM::Units::Length();
152 else
153 return fCutTubs->GetHighNorm().y() * ClhepVGM::Units::Length();
154}
155
156//_____________________________________________________________________________
158{
159 if (!fIsReflected)
160 return fCutTubs->GetLowNorm().z() * ClhepVGM::Units::Length();
161 else
162 return -fCutTubs->GetHighNorm().z() * ClhepVGM::Units::Length();
163}
164
165//_____________________________________________________________________________
167{
168 if (!fIsReflected)
169 return fCutTubs->GetHighNorm().x() * ClhepVGM::Units::Length();
170 else
171 return fCutTubs->GetLowNorm().x() * ClhepVGM::Units::Length();
172}
173
174//_____________________________________________________________________________
176{
177 if (!fIsReflected)
178 return fCutTubs->GetHighNorm().y() * ClhepVGM::Units::Length();
179 else
180 return fCutTubs->GetLowNorm().y() * ClhepVGM::Units::Length();
181}
182
183//_____________________________________________________________________________
185{
186 if (!fIsReflected)
187 return fCutTubs->GetHighNorm().z() * ClhepVGM::Units::Length();
188 else
189 return -fCutTubs->GetLowNorm().z() * ClhepVGM::Units::Length();
190}
static double Length()
Return CLHEP default length unit in VGM units.
Definition Units.cxx:83
static double Angle()
Return CLHEP default angle unit in VGM units.
Definition Units.cxx:92
VGM implementation for cut tubs solid in Geant4.
Definition Ctubs.h:38
virtual double NyHigh() const
Y-component of the normal unit vector to the cut plane in +z.
Definition Ctubs.cxx:175
virtual double NxHigh() const
X-component of the normal unit vector to the cut plane in +z.
Definition Ctubs.cxx:166
virtual double OuterRadius() const
Return the outside radius in mm.
Definition Ctubs.cxx:116
virtual double ZHalfLength() const
Return the half-length along the z axis in m.
Definition Ctubs.cxx:122
virtual double NzLow() const
Z-component of the normal unit vector to the cut plane in -z.
Definition Ctubs.cxx:157
virtual double NzHigh() const
Z-component of the normal unit vector to the cut plane in +z.
Definition Ctubs.cxx:184
virtual ~Ctubs()
Definition Ctubs.cxx:101
virtual double DeltaPhi() const
Return the opening angle of the segment in deg.
Definition Ctubs.cxx:134
virtual double NyLow() const
Y-component of the normal unit vector to the cut plane in -z.
Definition Ctubs.cxx:148
virtual double NxLow() const
X-component of the normal unit vector to the cut plane in -z.
Definition Ctubs.cxx:139
virtual std::string Name() const
Return the name of this solid.
Definition Ctubs.cxx:107
virtual double StartPhi() const
Return the starting angle of the segment in deg.
Definition Ctubs.cxx:128
virtual double InnerRadius() const
Return the inside radius in mm.
Definition Ctubs.cxx:110
static SolidMap * Instance()
Definition SolidMap.cxx:28
void AddSolid(VGM::ISolid *, G4VSolid *)
Definition SolidMap.cxx:59
BaseVGM utilities.
Definition utilities.h:23
VGM interfaces.
Definition VMedium.h:28