Geant4  10.02.p01
G4DrawVoxels.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4DrawVoxels.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 //
30 // class G4DrawVoxels
31 //
32 // Implementation
33 //
34 // Define G4DrawVoxelsDebug for debugging information on G4cout
35 //
36 // History:
37 // 03/08/1999 The G4VisAttributes have been made member data for
38 // lifetime reasons / visualisation L.G
39 // 29/07/1999 first comitted version L.G.
40 // --------------------------------------------------------------------
41 
42 #include "G4DrawVoxels.hh"
43 #include "G4AffineTransform.hh"
44 #include "G4SmartVoxelHeader.hh"
45 #include "G4LogicalVolume.hh"
46 #include "G4VSolid.hh"
47 #include "G4VVisManager.hh"
48 #include "G4Colour.hh"
51 
52 #define voxel_width 0
53 
54 // Private Constructor
55 //
57 {
62 }
63 
64 // Destructor
65 //
67 {
68 }
69 
70 // Methods that allow changing colors of the drawing
71 //
73  G4VisAttributes& VA_voxelY,
74  G4VisAttributes& VA_voxelZ)
75 {
76  fVoxelsVisAttributes[0]=VA_voxelX;
77  fVoxelsVisAttributes[1]=VA_voxelY;
78  fVoxelsVisAttributes[2]=VA_voxelZ;
79 }
80 
82 {
83  fBoundingBoxVisAttributes=VA_boundingbox;
84 }
85 
86 // ***************************************************************
87 
88 void
90  const G4SmartVoxelHeader* header,
91  G4VoxelLimits& limit,
92  G4PlacedPolyhedronList* ppl) const
93 {
94  // Let's draw the selected voxelisation now !
95 
96  G4VSolid* solid=lv->GetSolid();
97 
99  G4double xmax=0,xmin=0,ymax=0,ymin=0,zmax=0,zmin=0;
100 
101  if (lv->GetNoDaughters()<=0)
102  {
103  return;
104  }
105 
106  // Let's get the data for the voxelisation
107 
108  solid->CalculateExtent(kXAxis,limit,G4AffineTransform(),xmin,xmax);
109  // G4AffineTransform() is identity
110  solid->CalculateExtent(kYAxis,limit,G4AffineTransform(),ymin,ymax);
111  // extents according to the axis of the local frame
112  solid->CalculateExtent(kZAxis,limit,G4AffineTransform(),zmin,zmax);
113  dx=xmax-xmin;
114  dy=ymax-ymin;
115  dz=zmax-zmin;
116 
117  // Preparing the colored bounding polyhedronBox for the pVolume
118  //
119  G4PolyhedronBox bounding_polyhedronBox(dx*0.5,dy*0.5,dz*0.5);
120  bounding_polyhedronBox.SetVisAttributes(&fBoundingBoxVisAttributes);
121  G4ThreeVector t_centerofBoundingBox((xmin+xmax)*0.5,
122  (ymin+ymax)*0.5,
123  (zmin+zmax)*0.5);
124 
125  ppl->push_back(G4PlacedPolyhedron(bounding_polyhedronBox,
126  G4Translate3D(t_centerofBoundingBox)));
127 
128  G4ThreeVector t_FirstCenterofVoxelPlane;
129  const G4VisAttributes* voxelsVisAttributes=0;
130 
131  G4ThreeVector unit_translation_vector;
132  G4ThreeVector current_translation_vector;
133 
134  switch(header->GetAxis())
135  {
136  case kXAxis:
137  dx=voxel_width;
138  unit_translation_vector=G4ThreeVector(1,0,0);
139  t_FirstCenterofVoxelPlane=G4ThreeVector(xmin,(ymin+ymax)*0.5,
140  (zmin+zmax)*0.5);
141  voxelsVisAttributes=&fVoxelsVisAttributes[0];
142  break;
143  case kYAxis:
144  dy=voxel_width;
145  t_FirstCenterofVoxelPlane=G4ThreeVector((xmin+xmax)*0.5,ymin,
146  (zmin+zmax)*0.5);
147  unit_translation_vector=G4ThreeVector(0,1,0);
148  voxelsVisAttributes=&fVoxelsVisAttributes[1];
149  break;
150  case kZAxis:
151  dz=voxel_width;
152  t_FirstCenterofVoxelPlane=G4ThreeVector((xmin+xmax)*0.5,
153  (ymin+ymax)*0.5,zmin);
154  unit_translation_vector=G4ThreeVector(0,0,1);
155  voxelsVisAttributes=&fVoxelsVisAttributes[2];
156  break;
157  default:
158  break;
159  };
160 
161  G4PolyhedronBox voxel_plane(dx*0.5,dy*0.5,dz*0.5);
162  voxel_plane.SetVisAttributes(voxelsVisAttributes);
163 
164  G4SmartVoxelProxy* slice=header->GetSlice(0);
165  G4int slice_no=0,no_slices=header->GetNoSlices();
166  G4double beginning=header->GetMinExtent(),
167  step=(header->GetMaxExtent()-beginning)/no_slices;
168 
169  while (slice_no<no_slices)
170  {
171  if (slice->IsHeader())
172  {
173  G4VoxelLimits newlimit(limit);
174  newlimit.AddLimit(header->GetAxis(),beginning+step*slice_no,
175  beginning+step*(slice->GetHeader()->GetMaxEquivalentSliceNo()+1));
176  ComputeVoxelPolyhedra(lv,slice->GetHeader(),newlimit,ppl);
177  }
178  current_translation_vector=unit_translation_vector;
179  current_translation_vector*=step*slice_no;
180 
181  ppl->push_back(G4PlacedPolyhedron(voxel_plane,
182  G4Translate3D(current_translation_vector
183  +t_FirstCenterofVoxelPlane)));
184  slice_no=(slice->IsHeader()
185  ? slice->GetHeader()->GetMaxEquivalentSliceNo()+1
186  : slice->GetNode()->GetMaxEquivalentSliceNo()+1);
187  if (slice_no<no_slices) { slice=header->GetSlice(slice_no); }
188  }
189 }
190 
191 // ########################################################################
192 
195 {
197  G4VoxelLimits limits; // Working object for recursive call.
198  ComputeVoxelPolyhedra(lv,lv->GetVoxelHeader(),limits,pplist);
199  return pplist; //it s up to the calling program to destroy it then!
200 }
201 
203 {
205 
206  if (lv->GetNoDaughters()<=0)
207  {
208  return;
209  }
210 
211  // Computing the transformation according to the world volume
212  // (the drawing is directly in the world volume while the axis
213  // are relative to the mother volume of lv's daughter.)
214 
215  G4TouchableHistoryHandle aTouchable =
217  GetNavigatorForTracking()->CreateTouchableHistoryHandle();
218  G4AffineTransform globTransform =
219  aTouchable->GetHistory()->GetTopTransform().Inverse();
220  G4Transform3D transf3D(globTransform.NetRotation(),
221  globTransform.NetTranslation());
222 
224  if(pVVisManager)
225  {
226  // Drawing the bounding and voxel polyhedra for the pVolume
227  //
228  for (size_t i=0;i<pplist->size();i++)
229  {
230  pVVisManager->Draw((*pplist)[i].GetPolyhedron(),
231  (*pplist)[i].GetTransform()*transf3D);
232  }
233  }
234  else
235  {
236  G4Exception("G4DrawVoxels::DrawVoxels()",
237  "GeomNav1002", JustWarning,
238  "Pointer to visualization manager is null!");
239  }
240  delete pplist;
241 }
virtual void Draw(const G4Circle &, const G4Transform3D &objectTransformation=G4Transform3D())=0
G4SmartVoxelHeader * GetVoxelHeader() const
G4SmartVoxelHeader * GetHeader() const
void SetColour(const G4Colour &)
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const =0
G4int GetMaxEquivalentSliceNo() const
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
static G4VVisManager * GetConcreteInstance()
G4VisAttributes fVoxelsVisAttributes[3]
Definition: G4DrawVoxels.hh:79
void SetVoxelsVisAttributes(G4VisAttributes &, G4VisAttributes &, G4VisAttributes &)
Definition: G4DrawVoxels.cc:72
G4PlacedPolyhedronList * CreatePlacedPolyhedra(const G4LogicalVolume *) const
G4ThreeVector NetTranslation() const
int G4int
Definition: G4Types.hh:78
G4int GetMaxEquivalentSliceNo() const
void SetBoundingBoxVisAttributes(G4VisAttributes &)
Definition: G4DrawVoxels.cc:81
EAxis GetAxis() const
void AddLimit(const EAxis pAxis, const G4double pMin, const G4double pMax)
G4double GetMinExtent() const
void ComputeVoxelPolyhedra(const G4LogicalVolume *, const G4SmartVoxelHeader *, G4VoxelLimits &, G4PlacedPolyhedronList *) const
Definition: G4DrawVoxels.cc:89
HepGeom::Transform3D G4Transform3D
G4int GetNoSlices() const
G4double GetMaxExtent() const
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:80
static G4TransportationManager * GetTransportationManager()
G4RotationMatrix NetRotation() const
G4SmartVoxelNode * GetNode() const
HepGeom::Translate3D G4Translate3D
G4bool IsHeader() const
G4SmartVoxelProxy * GetSlice(G4int n) const
double G4double
Definition: G4Types.hh:76
#define voxel_width
Definition: G4DrawVoxels.cc:52
G4VisAttributes fBoundingBoxVisAttributes
Definition: G4DrawVoxels.hh:80
G4VSolid * GetSolid() const
void DrawVoxels(const G4LogicalVolume *lv) const
std::vector< G4PlacedPolyhedron > G4PlacedPolyhedronList