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