Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/visualization/HepRep/src/G4HepRepFileSceneHandler.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //
 28 // Joseph Perl  27th January 2002
 29 // A base class for a scene handler to export geometry and trajectories
 30 // to the HepRep xml file format.
 31 
 32 #include "G4HepRepFileSceneHandler.hh"
 33 #include "G4HepRepFile.hh"
 34 #include "G4HepRepMessenger.hh"
 35 #include "G4UIcommand.hh"
 36 
 37 #include "G4PhysicalConstants.hh"
 38 #include "G4SystemOfUnits.hh"
 39 #include "G4Version.hh"
 40 #include "G4VSolid.hh"
 41 #include "G4PhysicalVolumeModel.hh"
 42 #include "G4VPhysicalVolume.hh"
 43 #include "G4LogicalVolume.hh"
 44 #include "G4RotationMatrix.hh"
 45 #include "G4Material.hh"
 46 #include "G4Polymarker.hh"
 47 #include "G4Polyline.hh"
 48 #include "G4Text.hh"
 49 #include "G4Circle.hh"
 50 #include "G4Square.hh"
 51 #include "G4Polyhedron.hh"
 52 #include "G4VTrajectory.hh"
 53 #include "G4VTrajectoryPoint.hh"
 54 #include "G4TrajectoriesModel.hh"
 55 #include "G4VHit.hh"
 56 #include "G4AttDef.hh"
 57 #include "G4AttValue.hh"
 58 #include "G4AttCheck.hh"
 59 #include "G4VisManager.hh"
 60 #include "G4VisTrajContext.hh"
 61 #include "G4VTrajectoryModel.hh"
 62 
 63 // HepRep
 64 #include "G4HepRepFileXMLWriter.hh"
 65 
 66 G4int G4HepRepFileSceneHandler::fSceneIdCount = 0;
 67 // Counter for HepRep scene handlers.
 68 
 69 G4HepRepFileSceneHandler::G4HepRepFileSceneHandler(G4VGraphicsSystem& system,
 70                                                    const G4String& name)
 71   : G4VSceneHandler(system, fSceneIdCount++, name)
 72 {
 73   hepRepXMLWriter = ((G4HepRepFile*) (&system))->GetHepRepXMLWriter();
 74   fileCounter     = 0;
 75 
 76   inPrimitives2D       = false;
 77   warnedAbout3DText    = false;
 78   warnedAbout2DMarkers = false;
 79   haveVisible          = false;
 80   drawingTraj          = false;
 81   doneInitTraj         = false;
 82   drawingHit           = false;
 83   doneInitHit          = false;
 84   trajContext          = 0;
 85   trajAttValues        = 0;
 86   trajAttDefs          = 0;
 87   hitAttValues         = 0;
 88   hitAttDefs           = 0;
 89 }
 90 
 91 G4HepRepFileSceneHandler::~G4HepRepFileSceneHandler() {}
 92 
 93 void G4HepRepFileSceneHandler::BeginModeling()
 94 {
 95   G4VisManager* visManager        = G4VisManager::GetInstance();
 96   const G4VTrajectoryModel* model = visManager->CurrentTrajDrawModel();
 97   trajContext                     = &model->GetContext();
 98 
 99   G4VSceneHandler::BeginModeling();  // Required: see G4VSceneHandler.hh.
100 }
101 
102 void G4HepRepFileSceneHandler::EndModeling()
103 {
104   G4VSceneHandler::EndModeling();  // Required: see G4VSceneHandler.hh.
105 }
106 
107 void G4HepRepFileSceneHandler::BeginPrimitives2D(
108   const G4Transform3D& objectTransformation)
109 {
110 #ifdef G4HEPREPFILEDEBUG
111   G4cout << "G4HepRepFileSceneHandler::BeginPrimitives2D() " << G4endl;
112 #endif
113   inPrimitives2D = true;
114   G4VSceneHandler::BeginPrimitives2D(objectTransformation);
115 }
116 
117 void G4HepRepFileSceneHandler::EndPrimitives2D()
118 {
119 #ifdef G4HEPREPFILEDEBUG
120   G4cout << "G4HepRepFileSceneHandler::EndPrimitives2D() " << G4endl;
121 #endif
122   G4VSceneHandler::EndPrimitives2D();
123   inPrimitives2D = false;
124 }
125 
126 #ifdef G4HEPREPFILEDEBUG
127 void G4HepRepFileSceneHandler::PrintThings()
128 {
129   G4cout << "  with transformation " << fObjectTransformation;
130   G4PhysicalVolumeModel* pPVModel =
131     dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
132   if(pPVModel)
133   {
134     G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
135     G4LogicalVolume* pCurrentLV   = pPVModel->GetCurrentLV();
136     G4int currentDepth            = pPVModel->GetCurrentDepth();
137     G4cout << "\n  current physical volume: " << pCurrentPV->GetName()
138            << "\n  current logical volume: " << pCurrentLV->GetName()
139            << "\n  current depth of geometry tree: " << currentDepth;
140   }
141   G4cout << G4endl;
142 }
143 #endif
144 
145 void G4HepRepFileSceneHandler::AddSolid(const G4Box& box)
146 {
147 #ifdef G4HEPREPFILEDEBUG
148   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Box& box) called for "
149          << box.GetName() << G4endl;
150   PrintThings();
151 #endif
152 
153   if(drawingTraj)
154     return;
155 
156   if(drawingHit)
157     InitHit();
158 
159   haveVisible = false;
160   AddHepRepInstance("Prism", NULL);
161 
162   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
163 
164   // Get and check applicable vis attributes.
165   fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
166   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
167      messenger->getCullInvisibles())
168     return;
169 
170   hepRepXMLWriter->addPrimitive();
171 
172   G4double dx = box.GetXHalfLength();
173   G4double dy = box.GetYHalfLength();
174   G4double dz = box.GetZHalfLength();
175 
176   G4Point3D vertex1(G4Point3D(dx, dy, -dz));
177   G4Point3D vertex2(G4Point3D(dx, -dy, -dz));
178   G4Point3D vertex3(G4Point3D(-dx, -dy, -dz));
179   G4Point3D vertex4(G4Point3D(-dx, dy, -dz));
180   G4Point3D vertex5(G4Point3D(dx, dy, dz));
181   G4Point3D vertex6(G4Point3D(dx, -dy, dz));
182   G4Point3D vertex7(G4Point3D(-dx, -dy, dz));
183   G4Point3D vertex8(G4Point3D(-dx, dy, dz));
184 
185   vertex1 = (fObjectTransformation) *vertex1;
186   vertex2 = (fObjectTransformation) *vertex2;
187   vertex3 = (fObjectTransformation) *vertex3;
188   vertex4 = (fObjectTransformation) *vertex4;
189   vertex5 = (fObjectTransformation) *vertex5;
190   vertex6 = (fObjectTransformation) *vertex6;
191   vertex7 = (fObjectTransformation) *vertex7;
192   vertex8 = (fObjectTransformation) *vertex8;
193 
194   hepRepXMLWriter->addPoint(vertex1.x(), vertex1.y(), vertex1.z());
195   hepRepXMLWriter->addPoint(vertex2.x(), vertex2.y(), vertex2.z());
196   hepRepXMLWriter->addPoint(vertex3.x(), vertex3.y(), vertex3.z());
197   hepRepXMLWriter->addPoint(vertex4.x(), vertex4.y(), vertex4.z());
198   hepRepXMLWriter->addPoint(vertex5.x(), vertex5.y(), vertex5.z());
199   hepRepXMLWriter->addPoint(vertex6.x(), vertex6.y(), vertex6.z());
200   hepRepXMLWriter->addPoint(vertex7.x(), vertex7.y(), vertex7.z());
201   hepRepXMLWriter->addPoint(vertex8.x(), vertex8.y(), vertex8.z());
202 }
203 
204 void G4HepRepFileSceneHandler::AddSolid(const G4Cons& cons)
205 {
206 #ifdef G4HEPREPFILEDEBUG
207   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Cons& cons) called for "
208          << cons.GetName() << G4endl;
209   PrintThings();
210 #endif
211 
212   // HepRApp does not correctly represent the end faces of cones at
213   // non-standard angles, let the base class convert these solids to polygons.
214   G4RotationMatrix r = fObjectTransformation.getRotation();
215   G4bool linedUpWithAnAxis =
216     (std::fabs(r.phiX()) <= .001 || std::fabs(r.phiY()) <= .001 ||
217      std::fabs(r.phiZ()) <= .001 || std::fabs(r.phiX() - pi) <= .001 ||
218      std::fabs(r.phiY() - pi) <= .001 || std::fabs(r.phiZ() - pi) <= .001);
219   // G4cout << "Angle X:" << r.phiX() << ", Angle Y:" << r.phiY() << ", Angle
220   // Z:" << r.phiZ() << G4endl; G4cout << "linedUpWithAnAxis:" <<
221   // linedUpWithAnAxis << G4endl;
222 
223   // HepRep does not have a primitive for a cut cone,
224   // so if this cone is cut, let the base class convert this
225   // solid to polygons.
226   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
227   if(cons.GetDeltaPhiAngle() < twopi || !linedUpWithAnAxis ||
228      messenger->renderCylAsPolygons())
229   {
230     G4VSceneHandler::AddSolid(cons);  // Invoke default action.
231   }
232   else
233   {
234     if(drawingTraj)
235       return;
236 
237     if(drawingHit)
238       InitHit();
239 
240     haveVisible = false;
241     AddHepRepInstance("Cylinder", NULL);
242 
243     // Get and check applicable vis attributes.
244     fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
245     if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
246        messenger->getCullInvisibles())
247       return;
248 
249     G4Point3D vertex1(G4Point3D(0., 0., -cons.GetZHalfLength()));
250     G4Point3D vertex2(G4Point3D(0., 0., cons.GetZHalfLength()));
251 
252     vertex1 = (fObjectTransformation) *vertex1;
253     vertex2 = (fObjectTransformation) *vertex2;
254 
255     // Outer cylinder.
256     hepRepXMLWriter->addPrimitive();
257     hepRepXMLWriter->addAttValue("Radius1", messenger->getScale() *
258                                               cons.GetOuterRadiusMinusZ());
259     hepRepXMLWriter->addAttValue("Radius2", messenger->getScale() *
260                                               cons.GetOuterRadiusPlusZ());
261     hepRepXMLWriter->addPoint(vertex1.x(), vertex1.y(), vertex1.z());
262     hepRepXMLWriter->addPoint(vertex2.x(), vertex2.y(), vertex2.z());
263 
264     // Inner cylinder.
265     hepRepXMLWriter->addPrimitive();
266     hepRepXMLWriter->addAttValue("Radius1", messenger->getScale() *
267                                               cons.GetInnerRadiusMinusZ());
268     hepRepXMLWriter->addAttValue("Radius2", messenger->getScale() *
269                                               cons.GetInnerRadiusPlusZ());
270     hepRepXMLWriter->addPoint(vertex1.x(), vertex1.y(), vertex1.z());
271     hepRepXMLWriter->addPoint(vertex2.x(), vertex2.y(), vertex2.z());
272   }
273 }
274 
275 void G4HepRepFileSceneHandler::AddSolid(const G4Tubs& tubs)
276 {
277 #ifdef G4HEPREPFILEDEBUG
278   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Tubs& tubs) called for "
279          << tubs.GetName() << G4endl;
280   PrintThings();
281 #endif
282 
283   // HepRApp does not correctly represent the end faces of cylinders at
284   // non-standard angles, let the base class convert these solids to polygons.
285   G4RotationMatrix r = fObjectTransformation.getRotation();
286   G4bool linedUpWithAnAxis =
287     (std::fabs(r.phiX()) <= .001 || std::fabs(r.phiY()) <= .001 ||
288      std::fabs(r.phiZ()) <= .001 || std::fabs(r.phiX() - pi) <= .001 ||
289      std::fabs(r.phiY() - pi) <= .001 || std::fabs(r.phiZ() - pi) <= .001);
290   // G4cout << "Angle X:" << r.phiX() << ", Angle Y:" << r.phiY() << ", Angle
291   // Z:" << r.phiZ() << G4endl; G4cout << "linedUpWithAnAxis:" <<
292   // linedUpWithAnAxis << G4endl;
293 
294   // HepRep does not have a primitive for a cut cylinder,
295   // so if this cylinder is cut, let the base class convert this
296   // solid to polygons.
297   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
298   if(tubs.GetDeltaPhiAngle() < twopi || !linedUpWithAnAxis ||
299      messenger->renderCylAsPolygons())
300   {
301     G4VSceneHandler::AddSolid(tubs);  // Invoke default action.
302   }
303   else
304   {
305     if(drawingTraj)
306       return;
307 
308     if(drawingHit)
309       InitHit();
310 
311     haveVisible = false;
312     AddHepRepInstance("Cylinder", NULL);
313 
314     // Get and check applicable vis attributes.
315     fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
316     if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
317        messenger->getCullInvisibles())
318       return;
319 
320     G4Point3D vertex1(G4Point3D(0., 0., -tubs.GetZHalfLength()));
321     G4Point3D vertex2(G4Point3D(0., 0., tubs.GetZHalfLength()));
322 
323     vertex1 = (fObjectTransformation) *vertex1;
324     vertex2 = (fObjectTransformation) *vertex2;
325 
326     // Outer cylinder.
327     hepRepXMLWriter->addPrimitive();
328     hepRepXMLWriter->addAttValue("Radius1",
329                                  messenger->getScale() * tubs.GetOuterRadius());
330     hepRepXMLWriter->addAttValue("Radius2",
331                                  messenger->getScale() * tubs.GetOuterRadius());
332     hepRepXMLWriter->addPoint(vertex1.x(), vertex1.y(), vertex1.z());
333     hepRepXMLWriter->addPoint(vertex2.x(), vertex2.y(), vertex2.z());
334 
335     // Inner cylinder.
336     if(tubs.GetInnerRadius() != 0.)
337     {
338       hepRepXMLWriter->addPrimitive();
339       hepRepXMLWriter->addAttValue("Radius1", messenger->getScale() *
340                                                 tubs.GetInnerRadius());
341       hepRepXMLWriter->addAttValue("Radius2", messenger->getScale() *
342                                                 tubs.GetInnerRadius());
343       hepRepXMLWriter->addPoint(vertex1.x(), vertex1.y(), vertex1.z());
344       hepRepXMLWriter->addPoint(vertex2.x(), vertex2.y(), vertex2.z());
345     }
346   }
347 }
348 
349 void G4HepRepFileSceneHandler::AddSolid(const G4Trd& trd)
350 {
351 #ifdef G4HEPREPFILEDEBUG
352   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Trd& trd) called for "
353          << trd.GetName() << G4endl;
354   PrintThings();
355 #endif
356 
357   if(drawingTraj)
358     return;
359 
360   if(drawingHit)
361     InitHit();
362 
363   haveVisible = false;
364   AddHepRepInstance("Prism", NULL);
365 
366   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
367 
368   // Get and check applicable vis attributes.
369   fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
370   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
371      messenger->getCullInvisibles())
372     return;
373 
374   hepRepXMLWriter->addPrimitive();
375 
376   G4double dx1 = trd.GetXHalfLength1();
377   G4double dy1 = trd.GetYHalfLength1();
378   G4double dx2 = trd.GetXHalfLength2();
379   G4double dy2 = trd.GetYHalfLength2();
380   G4double dz  = trd.GetZHalfLength();
381 
382   G4Point3D vertex1(G4Point3D(dx1, dy1, -dz));
383   G4Point3D vertex2(G4Point3D(dx1, -dy1, -dz));
384   G4Point3D vertex3(G4Point3D(-dx1, -dy1, -dz));
385   G4Point3D vertex4(G4Point3D(-dx1, dy1, -dz));
386   G4Point3D vertex5(G4Point3D(dx2, dy2, dz));
387   G4Point3D vertex6(G4Point3D(dx2, -dy2, dz));
388   G4Point3D vertex7(G4Point3D(-dx2, -dy2, dz));
389   G4Point3D vertex8(G4Point3D(-dx2, dy2, dz));
390 
391   vertex1 = (fObjectTransformation) *vertex1;
392   vertex2 = (fObjectTransformation) *vertex2;
393   vertex3 = (fObjectTransformation) *vertex3;
394   vertex4 = (fObjectTransformation) *vertex4;
395   vertex5 = (fObjectTransformation) *vertex5;
396   vertex6 = (fObjectTransformation) *vertex6;
397   vertex7 = (fObjectTransformation) *vertex7;
398   vertex8 = (fObjectTransformation) *vertex8;
399 
400   hepRepXMLWriter->addPoint(vertex1.x(), vertex1.y(), vertex1.z());
401   hepRepXMLWriter->addPoint(vertex2.x(), vertex2.y(), vertex2.z());
402   hepRepXMLWriter->addPoint(vertex3.x(), vertex3.y(), vertex3.z());
403   hepRepXMLWriter->addPoint(vertex4.x(), vertex4.y(), vertex4.z());
404   hepRepXMLWriter->addPoint(vertex5.x(), vertex5.y(), vertex5.z());
405   hepRepXMLWriter->addPoint(vertex6.x(), vertex6.y(), vertex6.z());
406   hepRepXMLWriter->addPoint(vertex7.x(), vertex7.y(), vertex7.z());
407   hepRepXMLWriter->addPoint(vertex8.x(), vertex8.y(), vertex8.z());
408 }
409 
410 void G4HepRepFileSceneHandler::AddSolid(const G4Trap& trap)
411 {
412 #ifdef G4HEPREPFILEDEBUG
413   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Trap& trap) called for "
414          << trap.GetName() << G4endl;
415   PrintThings();
416 #endif
417   G4VSceneHandler::AddSolid(trap);  // Invoke default action.
418 }
419 
420 void G4HepRepFileSceneHandler::AddSolid(const G4Sphere& sphere)
421 {
422 #ifdef G4HEPREPFILEDEBUG
423   G4cout
424     << "G4HepRepFileSceneHandler::AddSolid(const G4Sphere& sphere) called for "
425     << sphere.GetName() << G4endl;
426   PrintThings();
427 #endif
428   G4VSceneHandler::AddSolid(sphere);  // Invoke default action.
429 }
430 
431 void G4HepRepFileSceneHandler::AddSolid(const G4Para& para)
432 {
433 #ifdef G4HEPREPFILEDEBUG
434   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Para& para) called for "
435          << para.GetName() << G4endl;
436   PrintThings();
437 #endif
438   G4VSceneHandler::AddSolid(para);  // Invoke default action.
439 }
440 
441 void G4HepRepFileSceneHandler::AddSolid(const G4Torus& torus)
442 {
443 #ifdef G4HEPREPFILEDEBUG
444   G4cout
445     << "G4HepRepFileSceneHandler::AddSolid(const G4Torus& torus) called for "
446     << torus.GetName() << G4endl;
447   PrintThings();
448 #endif
449   G4VSceneHandler::AddSolid(torus);  // Invoke default action.
450 }
451 
452 void G4HepRepFileSceneHandler::AddSolid(const G4Polycone& polycone)
453 {
454 #ifdef G4HEPREPFILEDEBUG
455   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Polycone& polycone) "
456             "called for "
457          << polycone.GetName() << G4endl;
458   PrintThings();
459 #endif
460   G4VSceneHandler::AddSolid(polycone);  // Invoke default action.
461 }
462 
463 void G4HepRepFileSceneHandler::AddSolid(const G4Polyhedra& polyhedra)
464 {
465 #ifdef G4HEPREPFILEDEBUG
466   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Polyhedra& polyhedra) "
467             "called for "
468          << polyhedra.GetName() << G4endl;
469   PrintThings();
470 #endif
471   G4VSceneHandler::AddSolid(polyhedra);  // Invoke default action.
472 }
473 
474 void G4HepRepFileSceneHandler::AddSolid(const G4Orb& orb)
475 {
476 #ifdef G4HEPREPFILEDEBUG
477   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Orb& orb) called for "
478          << orb.GetName() << G4endl;
479   PrintThings();
480 #endif
481   G4VSceneHandler::AddSolid(orb);  // Invoke default action.
482 }
483 
484 void G4HepRepFileSceneHandler::AddSolid(const G4Ellipsoid& ellipsoid)
485 {
486 #ifdef G4HEPREPFILEDEBUG
487   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4Ellipsoid& ellipsoid) "
488             "called for "
489          << ellipsoid.GetName() << G4endl;
490   PrintThings();
491 #endif
492   G4VSceneHandler::AddSolid(ellipsoid);  // Invoke default action.
493 }
494 
495 void G4HepRepFileSceneHandler::AddSolid(const G4TessellatedSolid& tess)
496 {
497 #ifdef G4HEPREPFILEDEBUG
498   G4cout << "G4HepRepFileSceneHandler::AddSolid(const G4TessellatedSolid& ) "
499             "called for "
500          << tess.GetName() << G4endl;
501   PrintThings();
502 #endif
503   G4VSceneHandler::AddSolid(tess);  // Invoke default action.
504 }
505 
506 void G4HepRepFileSceneHandler::AddSolid(const G4VSolid& solid)
507 {
508 #ifdef G4HEPREPFILEDEBUG
509   G4cout
510     << "G4HepRepFileSceneHandler::AddSolid(const G4Solid& solid) called for "
511     << solid.GetName() << G4endl;
512   PrintThings();
513 #endif
514   G4VSceneHandler::AddSolid(solid);  // Invoke default action.
515 }
516 
517 void G4HepRepFileSceneHandler::AddCompound(const G4VTrajectory& traj)
518 {
519 #ifdef G4HEPREPFILEDEBUG
520   G4cout << "G4HepRepFileSceneHandler::AddCompound(const G4VTrajectory&) "
521          << G4endl;
522 #endif
523 
524   G4TrajectoriesModel* pTrModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
525   if(!pTrModel)
526     G4Exception("G4HepRepFileSceneHandler::AddCompound(const G4VTrajectory&)",
527                 "vis-HepRep0001", FatalException, "Not a G4TrajectoriesModel.");
528 
529   // Pointers to hold trajectory attribute values and definitions.
530   std::vector<G4AttValue>* rawTrajAttValues = traj.CreateAttValues();
531   trajAttValues                             = new std::vector<G4AttValue>;
532   trajAttDefs                               = new std::map<G4String, G4AttDef>;
533 
534   // Iterators to use with attribute values and definitions.
535   std::vector<G4AttValue>::iterator iAttVal;
536   std::map<G4String, G4AttDef>::const_iterator iAttDef;
537   G4int i;
538 
539   // Get trajectory attributes and definitions in standard HepRep style
540   // (uniform units, 3Vectors decomposed).
541   if(rawTrajAttValues)
542   {
543     G4bool error = G4AttCheck(rawTrajAttValues, traj.GetAttDefs())
544                      .Standard(trajAttValues, trajAttDefs);
545     if(error)
546     {
547       G4cout
548         << "G4HepRepFileSceneHandler::AddCompound(traj):"
549            "\nERROR found during conversion to standard trajectory attributes."
550         << G4endl;
551     }
552 #ifdef G4HEPREPFILEDEBUG
553     G4cout << "G4HepRepFileSceneHandler::AddCompound(traj): standardised "
554               "attributes:\n"
555            << G4AttCheck(trajAttValues, trajAttDefs) << G4endl;
556 #endif
557     delete rawTrajAttValues;
558   }
559 
560   // Open the HepRep output file if it is not already open.
561   CheckFileOpen();
562 
563   // Add the Event Data Type if it hasn't already been added.
564   if(strcmp("Event Data", hepRepXMLWriter->prevTypeName[0]) != 0)
565   {
566     hepRepXMLWriter->addType("Event Data", 0);
567     hepRepXMLWriter->addInstance();
568   }
569 
570   // Add the Trajectories Type.
571   G4String previousName = hepRepXMLWriter->prevTypeName[1];
572   hepRepXMLWriter->addType("Trajectories", 1);
573 
574   // If this is the first trajectory of this event,
575   // specify attribute values common to all trajectories.
576   if(strcmp("Trajectories", previousName) != 0)
577   {
578     hepRepXMLWriter->addAttValue("Layer", 100);
579 
580     // Take all Trajectory attDefs from first trajectory.
581     // Would rather be able to get these attDefs without needing a reference
582     // from any particular trajectory, but don't know how to do that. Write out
583     // trajectory attribute definitions.
584     if(trajAttValues && trajAttDefs)
585     {
586       for(iAttVal = trajAttValues->begin(); iAttVal != trajAttValues->end();
587           ++iAttVal)
588       {
589         iAttDef = trajAttDefs->find(iAttVal->GetName());
590         if(iAttDef != trajAttDefs->end())
591         {
592           // Protect against incorrect use of Category.  Anything value other
593           // than the standard ones will be considered to be in the physics
594           // category.
595           G4String category = iAttDef->second.GetCategory();
596           if(strcmp(category, "Draw") != 0 &&
597              strcmp(category, "Physics") != 0 &&
598              strcmp(category, "Association") != 0 &&
599              strcmp(category, "PickAction") != 0)
600             category = "Physics";
601           hepRepXMLWriter->addAttDef(iAttVal->GetName(),
602                                      iAttDef->second.GetDesc(), category,
603                                      iAttDef->second.GetExtra());
604         }
605       }
606     }
607 
608     // Take all TrajectoryPoint attDefs from first point of first trajectory.
609     // Would rather be able to get these attDefs without needing a reference
610     // from any particular point, but don't know how to do that.
611     if((trajContext->GetDrawStepPts() || trajContext->GetDrawAuxPts()) &&
612        traj.GetPointEntries() > 0)
613     {
614       G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(0);
615 
616       // Pointers to hold trajectory point attribute values and definitions.
617       std::vector<G4AttValue>* rawPointAttValues =
618         aTrajectoryPoint->CreateAttValues();
619       std::vector<G4AttValue>* pointAttValues = new std::vector<G4AttValue>;
620       std::map<G4String, G4AttDef>* pointAttDefs =
621         new std::map<G4String, G4AttDef>;
622 
623       // Get first trajectory point's attributes and definitions in standard
624       // HepRep style (uniform units, 3Vectors decomposed).
625       if(rawPointAttValues)
626       {
627         G4bool error =
628           G4AttCheck(rawPointAttValues, aTrajectoryPoint->GetAttDefs())
629             .Standard(pointAttValues, pointAttDefs);
630         if(error)
631         {
632           G4cout << "G4HepRepFileSceneHandler::AddCompound(traj):"
633                     "\nERROR found during conversion to standard first point "
634                     "attributes."
635                  << G4endl;
636         }
637 
638         // Write out point attribute definitions.
639         if(pointAttValues && pointAttDefs)
640         {
641           for(iAttVal = pointAttValues->begin();
642               iAttVal != pointAttValues->end(); ++iAttVal)
643           {
644             iAttDef = pointAttDefs->find(iAttVal->GetName());
645             if(iAttDef != pointAttDefs->end())
646             {
647               // Protect against incorrect use of Category.  Anything value
648               // other than the standard ones will be considered to be in the
649               // physics category.
650               G4String category = iAttDef->second.GetCategory();
651               if(strcmp(category, "Draw") != 0 &&
652                  strcmp(category, "Physics") != 0 &&
653                  strcmp(category, "Association") != 0 &&
654                  strcmp(category, "PickAction") != 0)
655                 category = "Physics";
656               // Do not write out the Aux or Pos attribute.  Aux does not
657               // conform to the HepRep rule that each object can have only one
658               // instance of a given AttValue. Both of these attributes are
659               // redundant to actual position information of the point.
660               if(strcmp(iAttVal->GetName(), "Aux-X") != 0 &&
661                  strcmp(iAttVal->GetName(), "Aux-Y") != 0 &&
662                  strcmp(iAttVal->GetName(), "Aux-Z") != 0 &&
663                  strcmp(iAttVal->GetName(), "Pos-X") != 0 &&
664                  strcmp(iAttVal->GetName(), "Pos-Y") != 0 &&
665                  strcmp(iAttVal->GetName(), "Pos-Z") != 0)
666                 hepRepXMLWriter->addAttDef(iAttVal->GetName(),
667                                            iAttDef->second.GetDesc(), category,
668                                            iAttDef->second.GetExtra());
669             }
670           }
671         }
672         delete rawPointAttValues;
673       }
674 
675       // Clean up point attributes.
676       if(pointAttValues)
677         delete pointAttValues;
678       if(pointAttDefs)
679         delete pointAttDefs;
680     }
681   }  // end of special treatment for when this is the first trajectory.
682 
683   // Now that we have written out all of the attributes that are based on the
684   // trajectory's particulars, call base class to deconstruct trajectory into
685   // polyline and/or points (or nothing if trajectory is to be filtered out). If
686   // base class calls for drawing points, no points will actually be drawn there
687   // since we instead need to do point drawing from here (in order to obtain the
688   // points attributes, not available from AddPrimitive(...point).  Instead,
689   // such a call will just serve to set the flag that tells us that point
690   // drawing was requested for this trajectory (depends on several factors
691   // including trajContext and filtering).
692   drawingTraj  = true;
693   doneInitTraj = false;
694   G4VSceneHandler::AddCompound(traj);
695   drawingTraj = false;
696 
697   // Draw step points.
698   if(trajContext->GetDrawStepPts())
699   {
700     if(!doneInitTraj)
701       InitTrajectory();
702     // Create Trajectory Points as a subType of Trajectories.
703     // Note that we should create this heprep type even if there are no actual
704     // points. This allows the user to tell that points don't exist (admittedly
705     // odd) rather than that they were omitted by the drawing mode.
706     previousName = hepRepXMLWriter->prevTypeName[2];
707     hepRepXMLWriter->addType("Trajectory Step Points", 2);
708 
709     float redness;
710     float greenness;
711     float blueness;
712     G4int markSize;
713     G4bool visible;
714     G4bool square;
715     G4Colour colour = trajContext->GetStepPtsColour();
716     redness         = colour.GetRed();
717     greenness       = colour.GetGreen();
718     blueness        = colour.GetBlue();
719     markSize        = (G4int) trajContext->GetStepPtsSize();
720     visible         = (G4int) trajContext->GetStepPtsVisible();
721     square          = (trajContext->GetStepPtsType() == G4Polymarker::squares);
722 
723     // Avoiding drawing anything black on black.
724     if(redness == 0. && greenness == 0. && blueness == 0.)
725     {
726       redness   = 1.;
727       greenness = 1.;
728       blueness  = 1.;
729     }
730 
731     // Specify attributes common to all trajectory points.
732     if(strcmp("Trajectory Step Points", previousName) != 0)
733     {
734       hepRepXMLWriter->addAttValue("DrawAs", "Point");
735       hepRepXMLWriter->addAttValue("MarkColor", redness, greenness, blueness);
736       hepRepXMLWriter->addAttValue("MarkSize", markSize);
737       hepRepXMLWriter->addAttValue("Layer", 110);
738       hepRepXMLWriter->addAttValue("Visibility", visible);
739       if(square)
740         hepRepXMLWriter->addAttValue("MarkName", "square");
741       else
742         hepRepXMLWriter->addAttValue("MarkName", "dot");
743     }
744 
745     // Loop over all points on this trajectory.
746     for(i = 0; i < traj.GetPointEntries(); i++)
747     {
748       G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(i);
749 
750       // Each point is a separate instance of the type Trajectory Points.
751       hepRepXMLWriter->addInstance();
752 
753       // Pointers to hold trajectory point attribute values and definitions.
754       std::vector<G4AttValue>* rawPointAttValues =
755         aTrajectoryPoint->CreateAttValues();
756       std::vector<G4AttValue>* pointAttValues = new std::vector<G4AttValue>;
757       std::map<G4String, G4AttDef>* pointAttDefs =
758         new std::map<G4String, G4AttDef>;
759 
760       // Get trajectory point attributes and definitions in standard HepRep
761       // style (uniform units, 3Vectors decomposed).
762       if(rawPointAttValues)
763       {
764         G4bool error =
765           G4AttCheck(rawPointAttValues, aTrajectoryPoint->GetAttDefs())
766             .Standard(pointAttValues, pointAttDefs);
767         if(error)
768         {
769           G4cout
770             << "G4HepRepFileSceneHandler::AddCompound(traj):"
771                "\nERROR found during conversion to standard point attributes."
772             << G4endl;
773         }
774 
775         // Write out point attribute values.
776         if(pointAttValues)
777         {
778           for(iAttVal = pointAttValues->begin();
779               iAttVal != pointAttValues->end(); ++iAttVal)
780             // Do not write out the Aux or Pos attribute.  Aux does not conform
781             // to the HepRep rule that each object can have only one instance of
782             // a given AttValue. Both of these attributes are redundant to
783             // actual position information of the point.
784             if(strcmp(iAttVal->GetName(), "Aux-X") != 0 &&
785                strcmp(iAttVal->GetName(), "Aux-Y") != 0 &&
786                strcmp(iAttVal->GetName(), "Aux-Z") != 0 &&
787                strcmp(iAttVal->GetName(), "Pos-X") != 0 &&
788                strcmp(iAttVal->GetName(), "Pos-Y") != 0 &&
789                strcmp(iAttVal->GetName(), "Pos-Z") != 0)
790               hepRepXMLWriter->addAttValue(iAttVal->GetName(),
791                                            iAttVal->GetValue());
792         }
793       }
794 
795       // Clean up point attributes.
796       delete pointAttDefs;
797       delete pointAttValues;
798       delete rawPointAttValues;
799 
800       // Each trajectory point is made of a single primitive, a point.
801       hepRepXMLWriter->addPrimitive();
802       G4Point3D vertex = aTrajectoryPoint->GetPosition();
803       hepRepXMLWriter->addPoint(vertex.x(), vertex.y(), vertex.z());
804     }
805   }
806 
807   // Draw Auxiliary Points
808   if(trajContext->GetDrawAuxPts())
809   {
810     if(!doneInitTraj)
811       InitTrajectory();
812     // Create Trajectory Points as a subType of Trajectories.
813     // Note that we should create this heprep type even if there are no actual
814     // points. This allows the user to tell that points don't exist (admittedly
815     // odd) rather than that they were omitted by the drawing mode.
816     previousName = hepRepXMLWriter->prevTypeName[2];
817     hepRepXMLWriter->addType("Trajectory Auxiliary Points", 2);
818 
819     float redness;
820     float greenness;
821     float blueness;
822     G4int markSize;
823     G4bool visible;
824     G4bool square;
825     G4Colour colour = trajContext->GetAuxPtsColour();
826     redness         = colour.GetRed();
827     greenness       = colour.GetGreen();
828     blueness        = colour.GetBlue();
829     markSize        = (G4int) trajContext->GetAuxPtsSize();
830     visible         = (G4int) trajContext->GetAuxPtsVisible();
831     square          = (trajContext->GetAuxPtsType() == G4Polymarker::squares);
832 
833     // Avoiding drawing anything black on black.
834     if(redness == 0. && greenness == 0. && blueness == 0.)
835     {
836       redness   = 1.;
837       greenness = 1.;
838       blueness  = 1.;
839     }
840 
841     // Specify attributes common to all trajectory points.
842     if(strcmp("Trajectory Auxiliary Points", previousName) != 0)
843     {
844       hepRepXMLWriter->addAttValue("DrawAs", "Point");
845       hepRepXMLWriter->addAttValue("MarkColor", redness, greenness, blueness);
846       hepRepXMLWriter->addAttValue("MarkSize", markSize);
847       hepRepXMLWriter->addAttValue("Layer", 110);
848       hepRepXMLWriter->addAttValue("Visibility", visible);
849       if(square)
850         hepRepXMLWriter->addAttValue("MarkName", "Square");
851       else
852         hepRepXMLWriter->addAttValue("MarkName", "Dot");
853     }
854 
855     // Loop over all points on this trajectory.
856     for(i = 0; i < traj.GetPointEntries(); i++)
857     {
858       G4VTrajectoryPoint* aTrajectoryPoint = traj.GetPoint(i);
859 
860       // Each point is a separate instance of the type Trajectory Points.
861       hepRepXMLWriter->addInstance();
862 
863       // Pointers to hold trajectory point attribute values and definitions.
864       std::vector<G4AttValue>* rawPointAttValues =
865         aTrajectoryPoint->CreateAttValues();
866       std::vector<G4AttValue>* pointAttValues = new std::vector<G4AttValue>;
867       std::map<G4String, G4AttDef>* pointAttDefs =
868         new std::map<G4String, G4AttDef>;
869 
870       // Get trajectory point attributes and definitions in standard HepRep
871       // style (uniform units, 3Vectors decomposed).
872       if(rawPointAttValues)
873       {
874         G4bool error =
875           G4AttCheck(rawPointAttValues, aTrajectoryPoint->GetAttDefs())
876             .Standard(pointAttValues, pointAttDefs);
877         if(error)
878         {
879           G4cout
880             << "G4HepRepFileSceneHandler::AddCompound(traj):"
881                "\nERROR found during conversion to standard point attributes."
882             << G4endl;
883         }
884 
885         // Write out point attribute values.
886         if(pointAttValues)
887         {
888           for(iAttVal = pointAttValues->begin();
889               iAttVal != pointAttValues->end(); ++iAttVal)
890             // Do not write out the Aux or Pos attribute.  Aux does not conform
891             // to the HepRep rule that each object can have only one instance of
892             // a given AttValue. Both of these attributes are redundant to
893             // actual position information of the point.
894             if(strcmp(iAttVal->GetName(), "Aux-X") != 0 &&
895                strcmp(iAttVal->GetName(), "Aux-Y") != 0 &&
896                strcmp(iAttVal->GetName(), "Aux-Z") != 0 &&
897                strcmp(iAttVal->GetName(), "Pos-X") != 0 &&
898                strcmp(iAttVal->GetName(), "Pos-Y") != 0 &&
899                strcmp(iAttVal->GetName(), "Pos-Z") != 0)
900               hepRepXMLWriter->addAttValue(iAttVal->GetName(),
901                                            iAttVal->GetValue());
902         }
903       }
904 
905       // Clean up point attributes.
906       delete pointAttDefs;
907       delete pointAttValues;
908       delete rawPointAttValues;
909 
910       // Each trajectory point is made of a single primitive, a point.
911       G4Point3D vertex = aTrajectoryPoint->GetPosition();
912 
913       // Loop over auxiliary points associated with this Trajectory Point.
914       const std::vector<G4ThreeVector>* auxiliaries =
915         aTrajectoryPoint->GetAuxiliaryPoints();
916       if(0 != auxiliaries)
917       {
918         for(size_t iAux = 0; iAux < auxiliaries->size(); ++iAux)
919         {
920           const G4ThreeVector auxPos((*auxiliaries)[iAux]);
921           hepRepXMLWriter->addPrimitive();
922           hepRepXMLWriter->addPoint(auxPos.x(), auxPos.y(), auxPos.z());
923         }
924       }
925     }
926   }
927 }
928 
929 void G4HepRepFileSceneHandler::AddCompound(const G4VHit& hit)
930 {
931 #ifdef G4HEPREPFILEDEBUG
932   G4cout << "G4HepRepFileSceneHandler::AddCompound(G4VHit&) " << G4endl;
933 #endif
934 
935   // Pointers to hold hit attribute values and definitions.
936   std::vector<G4AttValue>* rawHitAttValues = hit.CreateAttValues();
937   hitAttValues                             = new std::vector<G4AttValue>;
938   hitAttDefs                               = new std::map<G4String, G4AttDef>;
939 
940   // Iterators to use with attribute values and definitions.
941   std::vector<G4AttValue>::iterator iAttVal;
942   std::map<G4String, G4AttDef>::const_iterator iAttDef;
943 
944   // Get hit attributes and definitions in standard HepRep style
945   // (uniform units, 3Vectors decomposed).
946   if(rawHitAttValues)
947   {
948     G4bool error = G4AttCheck(rawHitAttValues, hit.GetAttDefs())
949                      .Standard(hitAttValues, hitAttDefs);
950     if(error)
951     {
952       G4cout << "G4HepRepFileSceneHandler::AddCompound(hit):"
953                 "\nERROR found during conversion to standard hit attributes."
954              << G4endl;
955     }
956 #ifdef G4HEPREPFILEDEBUG
957     G4cout << "G4HepRepFileSceneHandler::AddCompound(hit): standardised "
958               "attributes:\n"
959            << G4AttCheck(hitAttValues, hitAttDefs) << G4endl;
960 #endif
961     delete rawHitAttValues;
962   }
963 
964   // Open the HepRep output file if it is not already open.
965   CheckFileOpen();
966 
967   // Add the Event Data Type if it hasn't already been added.
968   if(strcmp("Event Data", hepRepXMLWriter->prevTypeName[0]) != 0)
969   {
970     hepRepXMLWriter->addType("Event Data", 0);
971     hepRepXMLWriter->addInstance();
972   }
973 
974   // Find out the current HitType.
975   G4String hitType = "Hits";
976   if(hitAttValues)
977   {
978     G4bool found = false;
979     for(iAttVal = hitAttValues->begin();
980         iAttVal != hitAttValues->end() && !found; ++iAttVal)
981     {
982       if(strcmp(iAttVal->GetName(), "HitType") == 0)
983       {
984         hitType = iAttVal->GetValue();
985         found   = true;
986       }
987     }
988   }
989 
990   // Add the Hits Type.
991   G4String previousName = hepRepXMLWriter->prevTypeName[1];
992   hepRepXMLWriter->addType(hitType, 1);
993 
994   // If this is the first hit of this event,
995   // specify attribute values common to all hits.
996   if(strcmp(hitType, previousName) != 0)
997   {
998     hepRepXMLWriter->addAttValue("Layer", 130);
999 
1000     // Take all Hit attDefs from first hit.
1001     // Would rather be able to get these attDefs without needing a reference
1002     // from any particular hit, but don't know how to do that. Write out hit
1003     // attribute definitions.
1004     if(hitAttValues && hitAttDefs)
1005     {
1006       for(iAttVal = hitAttValues->begin(); iAttVal != hitAttValues->end();
1007           ++iAttVal)
1008       {
1009         iAttDef = hitAttDefs->find(iAttVal->GetName());
1010         if(iAttDef != hitAttDefs->end())
1011         {
1012           // Protect against incorrect use of Category.  Anything value other
1013           // than the standard ones will be considered to be in the physics
1014           // category.
1015           G4String category = iAttDef->second.GetCategory();
1016           if(strcmp(category, "Draw") != 0 &&
1017              strcmp(category, "Physics") != 0 &&
1018              strcmp(category, "Association") != 0 &&
1019              strcmp(category, "PickAction") != 0)
1020             category = "Physics";
1021           hepRepXMLWriter->addAttDef(iAttVal->GetName(),
1022                                      iAttDef->second.GetDesc(), category,
1023                                      iAttDef->second.GetExtra());
1024         }
1025       }
1026     }
1027   }  // end of special treatment for when this is the first hit.
1028 
1029   // Now that we have written out all of the attributes that are based on the
1030   // hit's particulars, call base class to deconstruct hit into a primitives.
1031   drawingHit  = true;
1032   doneInitHit = false;
1033   G4VSceneHandler::AddCompound(hit);  // Invoke default action.
1034   drawingHit = false;
1035 }
1036 
1037 void G4HepRepFileSceneHandler::InitTrajectory()
1038 {
1039   if(!doneInitTraj)
1040   {
1041     // For every trajectory, add an instance of Type Trajectory.
1042     hepRepXMLWriter->addInstance();
1043 
1044     // Write out the trajectory's attribute values.
1045     if(trajAttValues)
1046     {
1047       std::vector<G4AttValue>::iterator iAttVal;
1048       for(iAttVal = trajAttValues->begin(); iAttVal != trajAttValues->end();
1049           ++iAttVal)
1050         hepRepXMLWriter->addAttValue(iAttVal->GetName(), iAttVal->GetValue());
1051       delete trajAttValues;
1052     }
1053 
1054     // Clean up trajectory attributes.
1055     if(trajAttDefs)
1056       delete trajAttDefs;
1057 
1058     doneInitTraj = true;
1059   }
1060 }
1061 
1062 void G4HepRepFileSceneHandler::InitHit()
1063 {
1064   if(!doneInitHit)
1065   {
1066     // For every hit, add an instance of Type Hit.
1067     hepRepXMLWriter->addInstance();
1068 
1069     // Write out the hit's attribute values.
1070     if(hitAttValues)
1071     {
1072       std::vector<G4AttValue>::iterator iAttVal;
1073       for(iAttVal = hitAttValues->begin(); iAttVal != hitAttValues->end();
1074           ++iAttVal)
1075         hepRepXMLWriter->addAttValue(iAttVal->GetName(), iAttVal->GetValue());
1076       delete hitAttValues;
1077     }
1078 
1079     // Clean up hit attributes.
1080     if(hitAttDefs)
1081       delete hitAttDefs;
1082 
1083     doneInitHit = true;
1084   }
1085 }
1086 
1087 void G4HepRepFileSceneHandler::AddPrimitive(const G4Polyline& polyline)
1088 {
1089 #ifdef G4HEPREPFILEDEBUG
1090   G4cout << "G4HepRepFileSceneHandler::AddPrimitive(const G4Polyline& "
1091             "polyline) called:"
1092             "\n  polyline: "
1093          << polyline << G4endl;
1094   PrintThings();
1095 #endif
1096 
1097   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1098 
1099   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
1100      messenger->getCullInvisibles())
1101     return;
1102 
1103   if(inPrimitives2D)
1104   {
1105     if(!warnedAbout2DMarkers)
1106     {
1107       G4cout << "HepRepFile does not currently support 2D lines." << G4endl;
1108       warnedAbout2DMarkers = true;
1109     }
1110     return;
1111   }
1112 
1113   if(drawingTraj)
1114     InitTrajectory();
1115 
1116   if(drawingHit)
1117     InitHit();
1118 
1119   haveVisible = true;
1120   AddHepRepInstance("Line", polyline);
1121 
1122   hepRepXMLWriter->addPrimitive();
1123 
1124   for(size_t i = 0; i < polyline.size(); i++)
1125   {
1126     G4Point3D vertex = (fObjectTransformation) *polyline[i];
1127     hepRepXMLWriter->addPoint(vertex.x(), vertex.y(), vertex.z());
1128   }
1129 }
1130 
1131 void G4HepRepFileSceneHandler::AddPrimitive(const G4Polymarker& line)
1132 {
1133 #ifdef G4HEPREPFILEDEBUG
1134   G4cout
1135     << "G4HepRepFileSceneHandler::AddPrimitive(const G4Polymarker& line) called"
1136     << G4endl;
1137   PrintThings();
1138 #endif
1139 
1140   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1141 
1142   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
1143      messenger->getCullInvisibles())
1144     return;
1145 
1146   if(inPrimitives2D)
1147   {
1148     if(!warnedAbout2DMarkers)
1149     {
1150       G4cout << "HepRepFile does not currently support 2D lines." << G4endl;
1151       warnedAbout2DMarkers = true;
1152     }
1153     return;
1154   }
1155 
1156   MarkerSizeType sizeType;
1157   G4double size = GetMarkerSize(line, sizeType);
1158   if(sizeType == world)
1159     size = 4.;
1160 
1161   if(drawingTraj)
1162     return;
1163 
1164   if(drawingHit)
1165     InitHit();
1166 
1167   haveVisible = true;
1168   AddHepRepInstance("Point", line);
1169 
1170   hepRepXMLWriter->addAttValue("MarkName", "Dot");
1171   hepRepXMLWriter->addAttValue("MarkSize", (G4int) size);
1172 
1173   hepRepXMLWriter->addPrimitive();
1174 
1175   for(size_t i = 0; i < line.size(); i++)
1176   {
1177     G4Point3D vertex = (fObjectTransformation) *line[i];
1178     hepRepXMLWriter->addPoint(vertex.x(), vertex.y(), vertex.z());
1179   }
1180 }
1181 
1182 void G4HepRepFileSceneHandler::AddPrimitive(const G4Text& text)
1183 {
1184 #ifdef G4HEPREPFILEDEBUG
1185   G4cout << "G4HepRepFileSceneHandler::AddPrimitive(const G4Text& text) called:"
1186             "\n  text: "
1187          << text.GetText() << G4endl;
1188   PrintThings();
1189 #endif
1190 
1191   if(!inPrimitives2D)
1192   {
1193     if(!warnedAbout3DText)
1194     {
1195       G4cout << "HepRepFile does not currently support 3D text." << G4endl;
1196       G4cout
1197         << "HepRep browsers can directly display text attributes on request."
1198         << G4endl;
1199       G4cout << "See Application Developers Guide for how to attach attributes "
1200                 "to viewable objects."
1201              << G4endl;
1202       warnedAbout3DText = true;
1203     }
1204     return;
1205   }
1206 
1207   MarkerSizeType sizeType;
1208   G4double size = GetMarkerSize(text, sizeType);
1209   if(sizeType == world)
1210     size = 12.;
1211 
1212   haveVisible = true;
1213   AddHepRepInstance("Text", text);
1214 
1215   hepRepXMLWriter->addAttValue("VAlign", "Top");
1216   hepRepXMLWriter->addAttValue("HAlign", "Left");
1217   hepRepXMLWriter->addAttValue("FontName", "Arial");
1218   hepRepXMLWriter->addAttValue("FontStyle", "Plain");
1219   hepRepXMLWriter->addAttValue("FontSize", (G4int) size);
1220   hepRepXMLWriter->addAttValue("FontHasBanner", "TRUE");
1221   hepRepXMLWriter->addAttValue("FontBannerColor", "0,0,0");
1222 
1223   const G4Colour& colour = GetTextColour(text);
1224   float redness          = colour.GetRed();
1225   float greenness        = colour.GetGreen();
1226   float blueness         = colour.GetBlue();
1227 
1228   // Avoiding drawing anything black on black.
1229   if(redness == 0. && greenness == 0. && blueness == 0.)
1230   {
1231     redness   = 1.;
1232     greenness = 1.;
1233     blueness  = 1.;
1234   }
1235   hepRepXMLWriter->addAttValue("FontColor", redness, greenness, blueness);
1236 
1237   hepRepXMLWriter->addPrimitive();
1238 
1239   hepRepXMLWriter->addAttValue("Text", text.GetText());
1240   hepRepXMLWriter->addAttValue("VPos", .99 - text.GetYOffset());
1241   hepRepXMLWriter->addAttValue("HPos", text.GetXOffset());
1242 }
1243 
1244 void G4HepRepFileSceneHandler::AddPrimitive(const G4Circle& circle)
1245 {
1246 #ifdef G4HEPREPFILEDEBUG
1247   G4cout
1248     << "G4HepRepFileSceneHandler::AddPrimitive(const G4Circle& circle) called:"
1249        "\n  radius: "
1250     << circle.GetWorldRadius() << G4endl;
1251   PrintThings();
1252 #endif
1253 
1254   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1255 
1256   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
1257      messenger->getCullInvisibles())
1258     return;
1259 
1260   if(inPrimitives2D)
1261   {
1262     if(!warnedAbout2DMarkers)
1263     {
1264       G4cout << "HepRepFile does not currently support 2D circles." << G4endl;
1265       warnedAbout2DMarkers = true;
1266     }
1267     return;
1268   }
1269 
1270   MarkerSizeType sizeType;
1271   G4double size = GetMarkerSize(circle, sizeType);
1272   if(sizeType == world)
1273     size = 4.;
1274 
1275   if(drawingTraj)
1276     return;
1277 
1278   if(drawingHit)
1279     InitHit();
1280 
1281   haveVisible = true;
1282   AddHepRepInstance("Point", circle);
1283 
1284   hepRepXMLWriter->addAttValue("MarkName", "Dot");
1285   hepRepXMLWriter->addAttValue("MarkSize", (G4int) size);
1286 
1287   hepRepXMLWriter->addPrimitive();
1288 
1289   G4Point3D center = (fObjectTransformation) *circle.GetPosition();
1290   hepRepXMLWriter->addPoint(center.x(), center.y(), center.z());
1291 }
1292 
1293 void G4HepRepFileSceneHandler::AddPrimitive(const G4Square& square)
1294 {
1295 #ifdef G4HEPREPFILEDEBUG
1296   G4cout
1297     << "G4HepRepFileSceneHandler::AddPrimitive(const G4Square& square) called:"
1298        "\n  side: "
1299     << square.GetWorldRadius() << G4endl;
1300   PrintThings();
1301 #endif
1302 
1303   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1304 
1305   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
1306      messenger->getCullInvisibles())
1307     return;
1308 
1309   if(inPrimitives2D)
1310   {
1311     if(!warnedAbout2DMarkers)
1312     {
1313       G4cout << "HepRepFile does not currently support 2D squares." << G4endl;
1314       warnedAbout2DMarkers = true;
1315     }
1316     return;
1317   }
1318 
1319   MarkerSizeType sizeType;
1320   G4double size = GetMarkerSize(square, sizeType);
1321   if(sizeType == world)
1322     size = 4.;
1323 
1324   if(drawingTraj)
1325     return;
1326 
1327   if(drawingHit)
1328     InitHit();
1329 
1330   haveVisible = true;
1331   AddHepRepInstance("Point", square);
1332 
1333   hepRepXMLWriter->addAttValue("MarkName", "Square");
1334   hepRepXMLWriter->addAttValue("MarkSize", (G4int) size);
1335 
1336   hepRepXMLWriter->addPrimitive();
1337 
1338   G4Point3D center = (fObjectTransformation) *square.GetPosition();
1339   hepRepXMLWriter->addPoint(center.x(), center.y(), center.z());
1340 }
1341 
1342 void G4HepRepFileSceneHandler::AddPrimitive(const G4Polyhedron& polyhedron)
1343 {
1344 #ifdef G4HEPREPFILEDEBUG
1345   G4cout << "G4HepRepFileSceneHandler::AddPrimitive(const G4Polyhedron& "
1346             "polyhedron) called."
1347          << G4endl;
1348   PrintThings();
1349 #endif
1350 
1351   G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1352 
1353   if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
1354      messenger->getCullInvisibles())
1355     return;
1356 
1357   if(polyhedron.GetNoFacets() == 0)
1358     return;
1359 
1360   if(drawingTraj)
1361     return;
1362 
1363   if(drawingHit)
1364     InitHit();
1365 
1366   haveVisible = true;
1367   AddHepRepInstance("Polygon", polyhedron);
1368 
1369   G4Normal3D surfaceNormal;
1370   G4Point3D vertex;
1371 
1372   G4bool notLastFace;
1373   do
1374   {
1375     hepRepXMLWriter->addPrimitive();
1376     notLastFace = polyhedron.GetNextNormal(surfaceNormal);
1377 
1378     G4int edgeFlag = 1;
1379     G4bool notLastEdge;
1380     do
1381     {
1382       notLastEdge = polyhedron.GetNextVertex(vertex, edgeFlag);
1383       vertex      = (fObjectTransformation) *vertex;
1384       hepRepXMLWriter->addPoint(vertex.x(), vertex.y(), vertex.z());
1385     } while(notLastEdge);
1386   } while(notLastFace);
1387 }
1388 
1389 G4HepRepFileXMLWriter* G4HepRepFileSceneHandler::GetHepRepXMLWriter()
1390 {
1391   return hepRepXMLWriter;
1392 }
1393 
1394 void G4HepRepFileSceneHandler::AddHepRepInstance(const char* primName,
1395                                                  const G4Visible visible)
1396 {
1397 #ifdef G4HEPREPFILEDEBUG
1398   G4cout << "G4HepRepFileSceneHandler::AddHepRepInstance called." << G4endl;
1399 #endif
1400 
1401   // Open the HepRep output file if it is not already open.
1402   CheckFileOpen();
1403 
1404   G4VPhysicalVolume* pCurrentPV = 0;
1405   G4LogicalVolume* pCurrentLV   = 0;
1406   G4int currentDepth            = 0;
1407   G4PhysicalVolumeModel* pPVModel =
1408     dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1409   if(pPVModel)
1410   {
1411     pCurrentPV   = pPVModel->GetCurrentPV();
1412     pCurrentLV   = pPVModel->GetCurrentLV();
1413     currentDepth = pPVModel->GetCurrentDepth();
1414   }
1415 
1416 #ifdef G4HEPREPFILEDEBUG
1417   G4cout << "pCurrentPV:" << pCurrentPV
1418          << ", readyForTransients:" << fReadyForTransients << G4endl;
1419 #endif
1420 
1421   if(drawingTraj || drawingHit)
1422   {
1423     // In this case, HepRep type, layer and instance were already created
1424     // in the AddCompound method.
1425   }
1426   else if(fReadyForTransients)
1427   {
1428     if(strcmp("Event Data", hepRepXMLWriter->prevTypeName[0]) != 0)
1429     {
1430       hepRepXMLWriter->addType("Event Data", 0);
1431       hepRepXMLWriter->addInstance();
1432     }
1433 
1434     // Applications have the option of either calling AddSolid(G4VTrajectory&)
1435     // and AddSolid(G4VHits&), or of just decomposing these into simpler
1436     // primitives. In the former case, drawing will be handled above and will
1437     // include setting of physics attributes. In the latter case, which is an
1438     // older style of working, we end up drawing the trajectories and hits here,
1439     // where we have no access to physics attributes. We receive primitives
1440     // here.  We can figure out that these are transients, but we have to guess
1441     // exactly what these transients represent. We assume the primitives are
1442     // being used as in G4VTrajectory, hence we assume: Lines are Trajectories
1443     // Squares that come after we've seen trajectories are Auxiliary Points
1444     // Circles that come after we've seen trajectories are Step Points
1445     // Other primitives are Hits
1446 
1447     int layer;
1448 
1449     if(strcmp("Text", primName) == 0)
1450     {
1451       hepRepXMLWriter->addType("EventID", 1);
1452     }
1453     else
1454     {
1455       if(strcmp("Line", primName) == 0)
1456       {
1457         hepRepXMLWriter->addType("TransientPolylines", 1);
1458         layer = 100;
1459       }
1460       else
1461       {
1462         if(strcmp(hepRepXMLWriter->prevTypeName[1], "TransientPolylines") ==
1463              0 &&
1464            strcmp("Square", primName) == 0)
1465         {
1466           hepRepXMLWriter->addType("AuxiliaryPoints", 2);
1467           layer = 110;
1468         }
1469         else
1470         {
1471           if(strcmp(hepRepXMLWriter->prevTypeName[1], "TransientPolylines") ==
1472                0 &&
1473              strcmp("Circle", primName) == 0)
1474           {
1475             hepRepXMLWriter->addType("StepPoints", 2);
1476             layer = 120;
1477           }
1478           else
1479           {
1480             hepRepXMLWriter->addType("Hits", 1);
1481             layer = 130;
1482           }
1483         }
1484       }
1485       hepRepXMLWriter->addAttValue("Layer", layer);
1486     }
1487 
1488     hepRepXMLWriter->addInstance();
1489 
1490     // Handle Type declaration for Axes, Ruler, etc.
1491   }
1492   else if(pCurrentPV == 0)
1493   {
1494     if(strcmp("AxesEtc", hepRepXMLWriter->prevTypeName[0]) != 0)
1495     {
1496       hepRepXMLWriter->addType("AxesEtc", 0);
1497       hepRepXMLWriter->addInstance();
1498     }
1499 
1500     int layer;
1501 
1502     if(strcmp("Text", primName) == 0)
1503     {
1504       hepRepXMLWriter->addType("Text", 1);
1505     }
1506     else
1507     {
1508       if(strcmp("Line", primName) == 0)
1509       {
1510         hepRepXMLWriter->addType("Polylines", 1);
1511         layer = 100;
1512       }
1513       else
1514       {
1515         hepRepXMLWriter->addType("Points", 1);
1516         layer = 130;
1517       }
1518       hepRepXMLWriter->addAttValue("Layer", layer);
1519     }
1520 
1521     hepRepXMLWriter->addInstance();
1522 
1523     // Handle Type declaration for Detector Geometry,
1524     // replacing G4's top geometry level name "worldPhysical" with the
1525     // name "Detector Geometry".
1526   }
1527   else
1528   {
1529     // G4cout << "CurrentDepth" << currentDepth << G4endl;
1530     // G4cout << "currentName" << pCurrentPV->GetName() << G4endl;
1531     if(strcmp("Detector Geometry", hepRepXMLWriter->prevTypeName[0]) != 0)
1532     {
1533       // G4cout << "Adding Det Geom type" << G4endl;
1534       hepRepXMLWriter->addType("Detector Geometry", 0);
1535       hepRepXMLWriter->addInstance();
1536     }
1537 
1538     // Re-insert any layers of the hierarchy that were removed by G4's culling
1539     // process. Don't bother checking if same type name as last instance.
1540     if(strcmp(hepRepXMLWriter->prevTypeName[currentDepth + 1],
1541               pCurrentPV->GetName()) != 0)
1542     {
1543       // G4cout << "Looking for mother of:" << pCurrentLV->GetName() << G4endl;
1544       typedef G4PhysicalVolumeModel::G4PhysicalVolumeNodeID PVNodeID;
1545       typedef std::vector<PVNodeID> PVPath;
1546       const PVPath& drawnPVPath         = pPVModel->GetDrawnPVPath();
1547       PVPath::const_reverse_iterator ri = ++drawnPVPath.rbegin();
1548       G4int drawnMotherDepth;
1549       if(ri != drawnPVPath.rend())
1550       {
1551         // This volume has a mother.
1552         drawnMotherDepth = ri->GetNonCulledDepth();
1553         // G4cout << "drawnMotherDepth" << drawnMotherDepth << G4endl;
1554       }
1555       else
1556       {
1557         // This volume has no mother.  Must be a top level volume.
1558         drawnMotherDepth = -1;
1559         // G4cout << "Mother must be very top" << G4endl;
1560       }
1561 
1562       while(drawnMotherDepth < (currentDepth - 1))
1563       {
1564         G4String culledParentName = "Culled parent of " + pCurrentPV->GetName();
1565         // G4cout << "Inserting culled layer " << culledParentName << " at
1566         // depth:" << drawnMotherDepth+2 << G4endl;
1567         hepRepXMLWriter->addType(culledParentName, drawnMotherDepth + 2);
1568         hepRepXMLWriter->addInstance();
1569         drawnMotherDepth++;
1570       }
1571     }
1572 
1573     // Add the HepRepType for the current volume.
1574     hepRepXMLWriter->addType(pCurrentPV->GetName(), currentDepth + 1);
1575     hepRepXMLWriter->addInstance();
1576 
1577     G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1578 
1579     if(fpVisAttribs && (fpVisAttribs->IsVisible() == 0) &&
1580        messenger->getCullInvisibles())
1581       return;
1582 
1583     // Additional attributes.
1584     hepRepXMLWriter->addAttValue("Layer", hepRepXMLWriter->typeDepth);
1585     hepRepXMLWriter->addAttValue("LVol", pCurrentLV->GetName());
1586     G4Region* region    = pCurrentLV->GetRegion();
1587     G4String regionName = region ? region->GetName() : G4String("No region");
1588     hepRepXMLWriter->addAttValue("Region", regionName);
1589     hepRepXMLWriter->addAttValue("RootRegion", pCurrentLV->IsRootRegion());
1590     hepRepXMLWriter->addAttValue("Solid", pCurrentLV->GetSolid()->GetName());
1591     hepRepXMLWriter->addAttValue("EType",
1592                                  pCurrentLV->GetSolid()->GetEntityType());
1593     G4Material* material = pPVModel->GetCurrentMaterial();
1594     G4String matName = material ? material->GetName() : G4String("No material");
1595     hepRepXMLWriter->addAttValue("Material", matName);
1596     G4double matDensity = material ? material->GetDensity() : 0.;
1597     hepRepXMLWriter->addAttValue("Density", matDensity * m3 / kg);
1598     G4State matState = material ? material->GetState() : kStateUndefined;
1599     hepRepXMLWriter->addAttValue("State", matState);
1600     G4double matRadlen = material ? material->GetRadlen() : 0.;
1601     hepRepXMLWriter->addAttValue("Radlen", matRadlen / m);
1602   }
1603 
1604   hepRepXMLWriter->addAttValue("DrawAs", primName);
1605 
1606   // Handle color and visibility attributes.
1607   float redness;
1608   float greenness;
1609   float blueness;
1610   G4bool isVisible;
1611 
1612   if(fpVisAttribs || haveVisible)
1613   {
1614     G4Colour colour;
1615 
1616     if(fpVisAttribs)
1617     {
1618       colour    = fpVisAttribs->GetColour();
1619       isVisible = fpVisAttribs->IsVisible();
1620     }
1621     else
1622     {
1623       colour = visible.GetVisAttributes()->GetColour();
1624       isVisible =
1625         fpViewer->GetApplicableVisAttributes(visible.GetVisAttributes())
1626           ->IsVisible();
1627     }
1628 
1629     redness   = colour.GetRed();
1630     greenness = colour.GetGreen();
1631     blueness  = colour.GetBlue();
1632 
1633     // Avoiding drawing anything black on black.
1634     if(redness == 0. && greenness == 0. && blueness == 0.)
1635     {
1636       redness   = 1.;
1637       greenness = 1.;
1638       blueness  = 1.;
1639     }
1640   }
1641   else
1642   {
1643 #ifdef G4HEPREPFILEDEBUG
1644     G4cout
1645       << "G4HepRepFileSceneHandler::AddHepRepInstance using default colour."
1646       << G4endl;
1647 #endif
1648     redness   = 1.;
1649     greenness = 1.;
1650     blueness  = 1.;
1651     isVisible = true;
1652   }
1653 
1654   if(strcmp(primName, "Point") == 0)
1655     hepRepXMLWriter->addAttValue("MarkColor", redness, greenness, blueness);
1656   else
1657     hepRepXMLWriter->addAttValue("LineColor", redness, greenness, blueness);
1658 
1659   hepRepXMLWriter->addAttValue("Visibility", isVisible);
1660 }
1661 
1662 void G4HepRepFileSceneHandler::CheckFileOpen()
1663 {
1664 #ifdef G4HEPREPFILEDEBUG
1665   G4cout << "G4HepRepFileSceneHandler::CheckFileOpen called." << G4endl;
1666 #endif
1667 
1668   if(!hepRepXMLWriter->isOpen)
1669   {
1670     G4String newFileSpec;
1671 
1672     G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1673 
1674     if(messenger->getOverwrite())
1675     {
1676       newFileSpec =
1677         messenger->getFileDir() + messenger->getFileName() + ".heprep";
1678     }
1679     else
1680     {
1681       newFileSpec = messenger->getFileDir() + messenger->getFileName() +
1682                     G4UIcommand::ConvertToString(fileCounter) + ".heprep";
1683     }
1684 
1685     G4cout << "HepRepFile writing to " << newFileSpec << G4endl;
1686 
1687     hepRepXMLWriter->open(newFileSpec);
1688 
1689     if(!messenger->getOverwrite())
1690       fileCounter++;
1691 
1692     hepRepXMLWriter->addAttDef("Generator", "HepRep Data Generator", "Physics",
1693                                "");
1694     G4String versionString = G4Version;
1695     versionString          = versionString.substr(1, versionString.size() - 2);
1696     versionString = " Geant4 version " + versionString + "   " + G4Date;
1697     hepRepXMLWriter->addAttValue("Generator", versionString);
1698 
1699     hepRepXMLWriter->addAttDef("LVol", "Logical Volume", "Physics", "");
1700     hepRepXMLWriter->addAttDef("Region", "Cuts Region", "Physics", "");
1701     hepRepXMLWriter->addAttDef("RootRegion", "Root Region", "Physics", "");
1702     hepRepXMLWriter->addAttDef("Solid", "Solid Name", "Physics", "");
1703     hepRepXMLWriter->addAttDef("EType", "Entity Type", "Physics", "");
1704     hepRepXMLWriter->addAttDef("Material", "Material Name", "Physics", "");
1705     hepRepXMLWriter->addAttDef("Density", "Material Density", "Physics",
1706                                "kg/m3");
1707     hepRepXMLWriter->addAttDef("State", "Material State", "Physics", "");
1708     hepRepXMLWriter->addAttDef("Radlen", "Material Radiation Length", "Physics",
1709                                "m");
1710   }
1711 }
1712 
1713 void G4HepRepFileSceneHandler::ClearTransientStore()
1714 {
1715   // This is typically called after an update and before drawing hits
1716   // of the next event.  To simulate the clearing of "transients"
1717   // (hits, etc.) the detector is redrawn...
1718   if(fpViewer)
1719   {
1720     fpViewer->SetView();
1721     fpViewer->ClearView();
1722     fpViewer->DrawView();
1723   }
1724 }
1725