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