VGM Version 5.3
Loading...
Searching...
No Matches
TessellatedSolid.cxx
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
12//
13// Class TessellatedSolid
14// --------------------
15// VGM implementation for Geant4 TessellatedSolid solid.
16//
17// Author: Ivana Hrivnacova; IPN Orsay
18
21
22#include "ClhepVGM/Units.h"
23
25
26#include "G4QuadrangularFacet.hh"
27#include "G4ReflectedSolid.hh"
28#include "G4TessellatedSolid.hh"
29#include "G4TriangularFacet.hh"
30
31#include <iostream>
32#include <math.h>
33#include <set>
34
35//_____________________________________________________________________________
37 const std::string& name, std::vector<std::vector<VGM::ThreeVector> > facets)
38 : VGM::ISolid(),
39 VGM::ITessellatedSolid(),
40 BaseVGM::VTessellatedSolid(),
41 fName(name),
42 fTessellatedSolid(0)
43{
47
48 fTessellatedSolid = new G4TessellatedSolid(name);
49
50 // Add triangular facets
51 //
52 for (G4int i = 0; i < G4int(facets.size()); i++) {
53
54 std::vector<VGM::ThreeVector> facet = facets[i];
55
56 // check number of vertices
57 if (facet.size() != 3 && facet.size() != 4) {
58 std::cerr << "+++ Error +++" << std::endl;
59 std::cerr << " Number of vertices in a facet = " << facet.size()
60 << " has to be == 3 or 4" << std::endl;
61 exit(1);
62 }
63
64 if (facet.size() == 3) {
65 VGM::ThreeVector vertex0 = facet[0];
66 VGM::ThreeVector vertex1 = facet[1];
67 VGM::ThreeVector vertex2 = facet[2];
68
69 fTessellatedSolid->AddFacet(new G4TriangularFacet(
70 G4ThreeVector(vertex0[VGM::kDx] / ClhepVGM::Units::Length(),
73 G4ThreeVector(vertex1[VGM::kDx] / ClhepVGM::Units::Length(),
76 G4ThreeVector(vertex2[VGM::kDx] / ClhepVGM::Units::Length(),
79 ABSOLUTE));
80 }
81 else {
82 // Quadrangular facet
83 VGM::ThreeVector vertex0 = facet[0];
84 VGM::ThreeVector vertex1 = facet[1];
85 VGM::ThreeVector vertex2 = facet[2];
86 VGM::ThreeVector vertex3 = facet[3];
87
88 fTessellatedSolid->AddFacet(new G4QuadrangularFacet(
89 G4ThreeVector(vertex0[VGM::kDx] / ClhepVGM::Units::Length(),
92 G4ThreeVector(vertex1[VGM::kDx] / ClhepVGM::Units::Length(),
95 G4ThreeVector(vertex2[VGM::kDx] / ClhepVGM::Units::Length(),
98 G4ThreeVector(vertex3[VGM::kDx] / ClhepVGM::Units::Length(),
100 vertex3[VGM::kDz] / ClhepVGM::Units::Length()),
101 ABSOLUTE));
102 }
103 }
104 fTessellatedSolid->SetSolidClosed(true);
105
106 // G4cout << *fTessellatedSolid << G4endl;
107
108 Geant4GM::SolidMap::Instance()->AddSolid(this, fTessellatedSolid);
109}
110
111//_____________________________________________________________________________
113 G4TessellatedSolid* tessellated, G4ReflectedSolid* reflTessellated)
114 : VGM::ISolid(),
115 VGM::ITessellatedSolid(),
116 BaseVGM::VTessellatedSolid(),
117 fName(tessellated->GetName()),
118 fIsReflected(false),
119 fTessellatedSolid(tessellated)
120{
122
123 if (reflTessellated) {
124 fIsReflected = true;
125 Geant4GM::SolidMap::Instance()->AddSolid(this, reflTessellated);
126 }
127 else {
128 Geant4GM::SolidMap::Instance()->AddSolid(this, tessellated);
129 }
130}
131
132//_____________________________________________________________________________
134 : VGM::ISolid(),
135 VGM::ITessellatedSolid(),
136 BaseVGM::VTessellatedSolid(),
137 fName(),
138 fIsReflected(false),
139 fTessellatedSolid(0)
140{
142}
143
144//_____________________________________________________________________________
146 : VGM::ISolid(rhs),
147 VGM::ITessellatedSolid(rhs),
148 BaseVGM::VTessellatedSolid(rhs),
149 fName(),
150 fIsReflected(false),
151 fTessellatedSolid(0)
152{
154}
155
156//_____________________________________________________________________________
161
162//_____________________________________________________________________________
163void Geant4GM::TessellatedSolid::CheckFacetIndex(int ifacet) const
164{
165 if (ifacet < 0 || ifacet > NofFacets()) {
166 std::cerr << "+++ Error +++" << std::endl;
167 std::cerr << " Wrong facet index: " << ifacet << std::endl;
168 exit(1);
169 }
170
171 G4VFacet* facet = fTessellatedSolid->GetFacet(ifacet);
172 if (!facet) {
173 std::cerr << "+++ Error +++" << std::endl;
174 std::cerr << " Facet with index: " << ifacet << " not found."
175 << std::endl;
176 exit(1);
177 }
178}
179
180//_____________________________________________________________________________
181void Geant4GM::TessellatedSolid::CheckVertexIndex(int ifacet, int index) const
182{
183 CheckFacetIndex(ifacet);
184
185 if (index < 0 || index > NofVertices(ifacet)) {
186 std::cerr << "+++ Error +++" << std::endl;
187 std::cerr << " Wrong vertex index: " << index << " in " << ifacet
188 << " th facet." << std::endl;
189 exit(1);
190 }
191}
192
193//_____________________________________________________________________________
195{
196 return fTessellatedSolid->GetName();
197}
198
199//_____________________________________________________________________________
201{
202 return fTessellatedSolid->GetNumberOfFacets();
203}
204
205//_____________________________________________________________________________
207{
208 CheckFacetIndex(ifacet);
209
210 G4VFacet* facet = fTessellatedSolid->GetFacet(ifacet);
211
212 return facet->GetNumberOfVertices();
213}
214
215//_____________________________________________________________________________
217{
218 CheckVertexIndex(ifacet, index);
219
220 G4VFacet* facet = fTessellatedSolid->GetFacet(ifacet);
221
222 VGM::ThreeVector vertex;
223 vertex.push_back(facet->GetVertex(index).x() * ClhepVGM::Units::Length());
224 vertex.push_back(facet->GetVertex(index).y() * ClhepVGM::Units::Length());
225 vertex.push_back(facet->GetVertex(index).z() * ClhepVGM::Units::Length());
226
227 if (fIsReflected) vertex[VGM::kDz] = -vertex[VGM::kDz];
228
229 return vertex;
230}
static double Length()
Return CLHEP default length unit in VGM units.
Definition Units.cxx:83
static SolidMap * Instance()
Definition SolidMap.cxx:28
void AddSolid(VGM::ISolid *, G4VSolid *)
Definition SolidMap.cxx:59
VGM implementation for Geant4 tessellated solid.
virtual int NofVertices(int ifacet) const
Return the number of vertices in the the ifacet-th facet.
virtual VGM::ThreeVector Vertex(int ifacet, int index) const
Return the index-th vertex in the ifacet-th facet.
virtual int NofFacets() const
Return the number of facets.
virtual std::string Name() const
Return the name of this solid.
BaseVGM utilities.
Definition utilities.h:23
VGM interfaces.
Definition VMedium.h:28
std::vector< double > ThreeVector
Definition ThreeVector.h:27
@ kDx
Definition Transform.h:44
@ kDz
Definition Transform.h:46
@ kDy
Definition Transform.h:45