Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4HepRepSceneHandler.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$
27 //
28 
33 #include <stdio.h>
34 
35 #include "globals.hh"
36 #include <vector>
37 #include <iostream>
38 // NOTE not available on Solaris 5.2 and Linux g++ 2.95.2
39 // #include <sstream>
40 #include <iomanip>
41 #include <fstream>
42 #include <cmath>
43 
44 //HepRep
45 #include "HEPREP/HepRep.h"
46 #include "G4HepRepMessenger.hh"
47 
48 //G4
49 #include "G4PhysicalConstants.hh"
50 #include "G4Vector3D.hh"
51 #include "G4Version.hh"
52 #include "G4Types.hh"
53 #include "G4Point3D.hh"
54 #include "G4Normal3D.hh"
55 #include "G4Polyline.hh"
56 #include "G4Polymarker.hh"
57 #include "G4Polyhedron.hh"
58 #include "G4Circle.hh"
59 #include "G4Square.hh"
60 #include "G4Text.hh"
61 #include "G4NURBS.hh"
62 #include "G4VPhysicalVolume.hh"
63 #include "G4VisAttributes.hh"
64 #include "G4VSolid.hh"
65 #include "G4VTrajectory.hh"
66 #include "G4VTrajectoryPoint.hh"
67 #include "G4VHit.hh"
68 #include "G4Scene.hh"
69 #include "G4Material.hh"
70 #include "G4AttDef.hh"
71 #include "G4AttValue.hh"
72 #include "G4AttCheck.hh"
73 
74 // CHepRep
76 
77 // This
78 #include "G4HepRep.hh"
79 #include "G4HepRepSceneHandler.hh"
80 #include "G4HepRepViewer.hh"
81 
82 
83 using namespace HEPREP;
84 using namespace cheprep;
85 using namespace std;
86 
87 G4int G4HepRepSceneHandler::sceneIdCount = 0;
88 
89 //#define LDEBUG 1
90 //#define SDEBUG 1
91 //#define PDEBUG 1
92 
94  : G4VSceneHandler (system, sceneIdCount++, name),
95  geometryLayer ("Geometry"),
96  eventLayer ("Event"),
97  calHitLayer ("CalHit"),
98  trajectoryLayer ("Trajectory"),
99  hitLayer ("Hit"),
100  rootVolumeName ("Geometry"),
101  baseName (""),
102  eventNumberPrefix (""),
103  eventNumberSuffix (""),
104  eventNumber (1),
105  eventNumberWidth (-1),
106  extension (""),
107  writeBinary (false),
108  writeZip (false),
109  writeGZ (false),
110  writeMultipleFiles (false),
111  currentHit (NULL),
112  currentTrack (NULL),
113  _heprep (NULL),
114  _heprepGeometry (NULL)
115 {
116 
117 #ifdef LDEBUG
118  cout << "G4HepRepSceneHandler::G4HepRepSceneHandler: " << system << endl;
119 #endif
120 
121  materialState[kStateSolid] = G4String("Solid");
122  materialState[kStateLiquid] = G4String("Liquid");
123  materialState[kStateGas] = G4String("Gas");
124  materialState[kStateUndefined] = G4String("Undefined");
125 
126  factory = new XMLHepRepFactory();
127  writer = NULL;
128 
129  // opening of file deferred to closeHepRep();
130  openHepRep();
131 }
132 
133 
135 #ifdef LDEBUG
136  cout << "G4HepRepSceneHandler::~G4HepRepSceneHandler() " << endl;
137 #endif
138  close();
139 
140  delete factory;
141  factory = NULL;
142 
143  dynamic_cast<G4HepRep*>(GetGraphicsSystem())->removeSceneHandler();
144 }
145 
146 
147 void G4HepRepSceneHandler::open(G4String name) {
148  if (writer != NULL) return;
149 
150  if (name == "stdout") {
151 #ifdef LDEBUG
152  cout << "G4HepRepSceneHandler::Open() stdout" << endl;
153 #endif
154  writer = factory->createHepRepWriter(&cout, false, false);
155  out = NULL;
156  baseName = name;
157  eventNumberPrefix = "";
158  eventNumberSuffix = "";
159  extension = "";
160  writeBinary = false;
161  writeZip = false;
162  writeGZ = false;
163  writeMultipleFiles = false;
164  eventNumber = 0;
165  eventNumberWidth = 0;
166  } else if (name == "stderr") {
167 #ifdef LDEBUG
168  cout << "G4HepRepSceneHandler::Open() stderr" << endl;
169 #endif
170  writer = factory->createHepRepWriter(&cerr, false, false);
171  out = NULL;
172  baseName = name;
173  eventNumberPrefix = "";
174  eventNumberSuffix = "";
175  extension = "";
176  writeBinary = false;
177  writeZip = false;
178  writeGZ = false;
179  writeMultipleFiles = false;
180  eventNumber = 0;
181  eventNumberWidth = 0;
182  } else {
183 #ifdef LDEBUG
184  cout << "G4HepRepSceneHandler::Open() " << name << endl;
185 #endif
186  if (eventNumberWidth < 0) {
187  // derive filename(s)
188  // check for extensions
189  const unsigned int numberOfExtensions = 8;
190  string ext[numberOfExtensions] = {".heprep", ".heprep.xml", ".heprep.zip", ".heprep.gz",
191  ".bheprep", ".bheprep.xml", ".bheprep.zip", ".bheprep.gz"};
192  unsigned int i=0;
193  while (i < numberOfExtensions) {
194  int dot = name.size() - ext[i].size();
195  if ((dot >= 0) &&
196  (name.substr(dot, ext[i].size()) == ext[i])) break;
197  i++;
198  }
199 
200  if (i != numberOfExtensions) {
201  extension = ext[i];
202  writeBinary = i >= (numberOfExtensions/2);
203  writeZip = (i == 2) || (i == 6);
204  writeGZ = (i == 3) || (i == 7);
205 
206  int dot = name.length() - extension.length();
207  baseName = (dot >= 0) ? name.substr(0, dot) : "";
208 
209 #ifndef G4LIB_USE_ZLIB
210  if (writeGZ) {
211  cerr << endl;
212  cerr << "WARNING: the .gz output file you are creating will be a plain file," << endl;
213  cerr << " since compression support (ZLIB) was not compiled into the Geant4." << endl;
214  cerr << " To avoid confusion with real gz files, the output filename has been" << endl;
215  cerr << " extended with the name '.no-gz'." << endl;
216  cerr << " A plain heprep or bheprep file can be fairly large." << endl;
217  }
218  if (writeZip) {
219  cerr << endl;
220  cerr << "WARNING: the .zip output file you are creating will not be compressed," << endl;
221  cerr << " since compression support (ZLIB) was not compiled into the Geant4." << endl;
222  cerr << " A zip file containing non-compressed heprep or bheprep files can" << endl;
223  cerr << " be fairly large." << endl;
224  }
225  if (writeGZ || writeZip) {
226  cerr << "SOLUTION: To add compression support using ZLIB, you need to:" << endl;
227  cerr << " 1. Define G4LIB_USE_ZLIB and recompile the visualization category." << endl;
228  cerr << " 2. Optionally define G4_LIB_BUILD_ZLIB if your system does not have" << endl;
229  cerr << " zlib installed (e.g. WIN32-VC)." << endl;
230  cerr << " 3. Relink your application code." << endl;
231  cerr << endl;
232  }
233  if (writeGZ) {
234  extension = extension + ".no-gz";
235  writeGZ = false;
236  }
237 #endif // G4LIB_USE_ZLIB
238  } else {
239  // Default for no extension
240  extension = ".heprep.zip";
241  writeBinary = false;
242  writeZip = true;
243  writeGZ = false;
244  baseName = name;
245  }
246 
247  writeMultipleFiles = false;
248  int startDigit = -1; int endDigit = -1;
249 
251 
252  string suffix = messenger->getEventNumberSuffix();
253  if (suffix != "") {
254  // look for 0000 pattern in suffix
255  endDigit = suffix.length()-1;
256  while (endDigit >= 0) {
257  if (isdigit(suffix.at(endDigit))) break;
258  endDigit--;
259  }
260  if (endDigit < 0) {
261  cerr << "/vis/heprep/appendEventNumberSuffix contains no digits" << endl;
262  } else {
263  writeMultipleFiles = true;
264  startDigit = endDigit;
265  while (startDigit >= 0) {
266  if (!isdigit(suffix.at(startDigit))) break;
267  startDigit--;
268  }
269  startDigit++;
270  }
271  }
272 
273  if (writeMultipleFiles) {
274  eventNumberPrefix = suffix.substr(0, startDigit);
275  eventNumber = atoi(suffix.substr(startDigit, endDigit).c_str());
276  eventNumberWidth = endDigit +1 - startDigit;
277  eventNumberSuffix = suffix.substr(endDigit+1);
278  } else {
279  // open single file here
280  openFile(baseName+extension);
281 
282  eventNumber = 1;
283  eventNumberWidth = 10;
284  eventNumberPrefix = "";
285  eventNumberSuffix = "";
286  }
287  }
288  }
289 }
290 
291 
293 #ifdef LDEBUG
294  cout << "G4HepRepSceneHandler::OpenHepRep() " << endl;
295 #endif
296 
297  if (_heprep != NULL) return;
298 
299  // all done on demand, once pointers are set to NULL
300  _heprepGeometry = NULL;
301  _geometryInstanceTree = NULL;
302  _geometryRootInstance = NULL;
303  _geometryInstance.clear();
304  _geometryTypeTree = NULL;
305  _geometryRootType = NULL;
306  _geometryTypeName.clear();
307  _geometryType.clear();
308  _eventInstanceTree = NULL;
309  _eventInstance = NULL;
310  _eventTypeTree = NULL;
311  _eventType = NULL;
312  _trajectoryType = NULL;
313  _hitType = NULL;
314  _calHitType = NULL;
315  _calHitFaceType = NULL;
316 }
317 
318 
323  if (_heprep == NULL) return true;
324 
325 #ifdef LDEBUG
326  cout << "G4HepRepSceneHandler::CloseHepRep() start" << endl;
327 #endif
328 
329  // if this is the final close, then there should not be any event pending to be written.
330  if (final) {
331  if (_eventInstanceTree != NULL) {
332  cerr << "WARNING: you probably used '/vis/viewer/endOfEventAction accumulate' and "
333  << "forgot to call /vis/viewer/update before exit. No event written." << endl;
334  }
335  } else {
336 
338 
339  // add geometry to the heprep if there is an event (separate geometries are written
340  // using DrawView() called from /vis/viewer/flush)
341  if (_eventInstanceTree != NULL) {
343 
344  // couple geometry
345 
346  if ( messenger->appendGeometry()) {
347  // couple geometry to event if geometry was written
348  if ((_geometryInstanceTree != NULL)) {
349  getEventInstanceTree()->addInstanceTree(getGeometryInstanceTree());
350  }
351  } else {
352  char name[128];
353  if (writeMultipleFiles) {
354  sprintf(name, "%s%s%s#%s", baseName.c_str(), "-geometry", extension.c_str(), "G4GeometryData");
355  } else {
356  sprintf(name, "%s%s#%s", "geometry", (writeBinary ? ".bheprep" : ".heprep"), "G4GeometryData");
357  }
358  getEventInstanceTree()->addInstanceTree(factory->createHepRepTreeID(name, "1.0"));
359  }
360  }
361 
362  // force inclusion of all subtypes of event
363  if (_eventInstanceTree != NULL) {
364  getEventType();
365  getTrajectoryType();
366  getHitType();
367  getCalHitType();
368  getCalHitFaceType();
369  }
370 
371  // Give this HepRep all of the layer order info for both geometry and event,
372  // since these will both end up in a single HepRep.
373  writeLayers(_heprepGeometry);
374  writeLayers(_heprep);
375 
376  // open heprep file
377  if (writer == NULL) {
378  open((GetScene() == NULL) ? G4String("G4HepRepOutput.heprep.zip") : GetScene()->GetName());
379  }
380 
381  // write out separate geometry
382  if (! messenger->appendGeometry() && (_heprepGeometry != NULL)) {
383  if (writeMultipleFiles) {
384  char fileName[128];
385  sprintf(fileName, "%s%s%s", baseName.c_str(), "-geometry", extension.c_str());
386  openFile(G4String(fileName));
387  }
388 
389  char name[128];
390  sprintf(name, "%s%s", "geometry", (writeBinary ? ".bheprep" : ".heprep"));
391  if (!writeMultipleFiles) {
392  writer->addProperty("RecordLoop.ignore", name);
393  }
394 
395  writer->write(_heprepGeometry, G4String(name));
396 
397  delete _heprepGeometry;
398  _heprepGeometry = NULL;
399 
400  if (writeMultipleFiles) closeFile();
401  }
402 
403  if (writeMultipleFiles) {
404 // NOTE: does not work on Solaris 5.2 and Linux 2.95.2
405 // stringstream fileName;
406 // fileName << baseName << eventNumberPrefix << setw(eventNumberWidth) << setfill('0') << eventNumber << eventNumberSuffix << extension;
407 // openFile(fileName.str());
408 // Use instead:
409  char fileName[128];
410  char fileFormat[128];
411  sprintf(fileFormat, "%s%d%s", "%s%s%0", eventNumberWidth, "d%s%s");
412  sprintf(fileName, fileFormat, baseName.c_str(), eventNumberPrefix.c_str(), eventNumber, eventNumberSuffix.c_str(), extension.c_str());
413  openFile(G4String(fileName));
414  }
415 
416  // write out the heprep
417 // NOTE: does not work on Solaris 5.2 and Linux 2.95.2
418 // stringstream eventName;
419 // eventName << "event-" << setw(eventNumberWidth) << setfill('0') << eventNumber << (writeBinary ? ".bheprep" : ".heprep");
420 // writer->write(_heprep, eventName.str());
421 // Use instead:
422  char eventName[128];
423  char eventFormat[128];
424  sprintf(eventFormat, "%s%d%s%s", "event-%0", eventNumberWidth, "d", (writeBinary ? ".bheprep" : ".heprep"));
425  sprintf(eventName, eventFormat, eventNumber);
426  writer->write(_heprep, G4String(eventName));
427 
428  eventNumber++;
429  }
430 
431  delete _heprep;
432  _heprep = NULL;
433 
434  if (writeMultipleFiles) closeFile();
435 
436  return true;
437 }
438 
439 
440 void G4HepRepSceneHandler::close() {
441 
442 #ifdef LDEBUG
443  cout << "G4HepRepSceneHandler::Close() " << endl;
444 #endif
445 
446  if (writer == NULL) return;
447 
448  if (!writeMultipleFiles) {
449  closeHepRep(true);
450  closeFile();
451  }
452 
453  G4HepRepViewer* viewer = dynamic_cast<G4HepRepViewer*>(GetCurrentViewer());
454  viewer->reset();
455 }
456 
458  out = new ofstream(name.c_str(), std::ios::out | std::ios::binary );
459  writer = factory->createHepRepWriter(out, writeZip, writeZip || writeGZ);
460 }
461 
463  writer->close();
464  delete writer;
465  writer = NULL;
466 
467  delete out;
468  out = NULL;
469 }
470 
471 void G4HepRepSceneHandler::writeLayers(HepRep* heprep) {
472  if (heprep == NULL) return;
473  heprep->addLayer(geometryLayer);
474  heprep->addLayer(eventLayer);
475  heprep->addLayer(calHitLayer);
476  heprep->addLayer(trajectoryLayer);
477  heprep->addLayer(hitLayer);
478 }
479 
481 #ifdef SDEBUG
482  cout << "G4HepRepSceneHandler::BeginModeling() " << endl;
483 #endif
485 }
486 
487 
489 #ifdef SDEBUG
490  cout << "G4HepRepSceneHandler::EndModeling() " << endl;
491 #endif
493 }
494 
496 #ifdef SDEBUG
497  cout << "G4HepRepSceneHandler::AddSolid(const G4Box& box)" << endl;
498 #endif
499 
500  if (dontWrite()) return;
501 
503 
504  if (! messenger->useSolids()) {
506  return;
507  }
508 
509  G4double dx = box.GetXHalfLength();
510  G4double dy = box.GetYHalfLength();
511  G4double dz = box.GetZHalfLength();
512 
513  G4Point3D vertex1(G4Point3D( dx, dy,-dz));
514  G4Point3D vertex2(G4Point3D( dx,-dy,-dz));
515  G4Point3D vertex3(G4Point3D(-dx,-dy,-dz));
516  G4Point3D vertex4(G4Point3D(-dx, dy,-dz));
517  G4Point3D vertex5(G4Point3D( dx, dy, dz));
518  G4Point3D vertex6(G4Point3D( dx,-dy, dz));
519  G4Point3D vertex7(G4Point3D(-dx,-dy, dz));
520  G4Point3D vertex8(G4Point3D(-dx, dy, dz));
521 
522  vertex1 = (transform) * vertex1;
523  vertex2 = (transform) * vertex2;
524  vertex3 = (transform) * vertex3;
525  vertex4 = (transform) * vertex4;
526  vertex5 = (transform) * vertex5;
527  vertex6 = (transform) * vertex6;
528  vertex7 = (transform) * vertex7;
529  vertex8 = (transform) * vertex8;
530 
531  HepRepInstance* instance = getGeometryOrEventInstance(getCalHitType());
532  addAttributes(instance, getCalHitType());
533 
534  setAttribute(instance, "DrawAs", G4String("Prism"));
535 
536  setVisibility(instance, box);
537  setLine(instance, box);
538  setColor(instance, getColorFor(box));
539 
540  factory->createHepRepPoint(instance, vertex1.x(), vertex1.y(), vertex1.z());
541  factory->createHepRepPoint(instance, vertex2.x(), vertex2.y(), vertex2.z());
542  factory->createHepRepPoint(instance, vertex3.x(), vertex3.y(), vertex3.z());
543  factory->createHepRepPoint(instance, vertex4.x(), vertex4.y(), vertex4.z());
544  factory->createHepRepPoint(instance, vertex5.x(), vertex5.y(), vertex5.z());
545  factory->createHepRepPoint(instance, vertex6.x(), vertex6.y(), vertex6.z());
546  factory->createHepRepPoint(instance, vertex7.x(), vertex7.y(), vertex7.z());
547  factory->createHepRepPoint(instance, vertex8.x(), vertex8.y(), vertex8.z());
548 }
549 
550 
552 #ifdef SDEBUG
553  cout << "G4HepRepSceneHandler::AddSolid(const G4Cons& cons)" << endl;
554 #endif
555 
556  if (dontWrite()) return;
557 
559 
560  if (! messenger->useSolids() || (cons.GetDeltaPhiAngle() < twopi)) {
562  return;
563  }
564 
565  G4PhysicalVolumeModel* pPVModel =
566  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
567  if (!pPVModel) {
569  return;
570  }
571 
572  G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
573  G4int currentDepth = pPVModel->GetCurrentDepth();
574  G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
575 
576  G4Point3D vertex1(G4Point3D( 0., 0., cons.GetZHalfLength()));
577  G4Point3D vertex2(G4Point3D( 0., 0.,-cons.GetZHalfLength()));
578 
579  vertex1 = (transform) * vertex1;
580  vertex2 = (transform) * vertex2;
581 
582  HepRepInstance* instance = getGeometryInstance(pCurrentLV, pCurrentMaterial, currentDepth);
583  setAttribute(instance, "DrawAs", G4String("Cylinder"));
584 
585  setVisibility(instance, cons);
586  setLine(instance, cons);
587  setColor(instance, getColorFor(cons));
588 
589  HepRepType* type = getGeometryType(pCurrentLV->GetName(), currentDepth);
590 
591  // Outer cylinder.
592  HepRepInstance* outer = factory->createHepRepInstance(instance, type);
593  outer->addAttValue("pickParent",true);
594  outer->addAttValue("showParentAttributes",true);
595 
596  HepRepPoint* op1 = factory->createHepRepPoint(outer, vertex1.x(), vertex1.y(), vertex1.z());
597  op1->addAttValue("Radius",cons.GetOuterRadiusPlusZ());
598 
599  HepRepPoint* op2 = factory->createHepRepPoint(outer, vertex2.x(), vertex2.y(), vertex2.z());
600  op2->addAttValue("Radius",cons.GetOuterRadiusMinusZ());
601 
602  // Inner cylinder.
603  HepRepInstance* inner = factory->createHepRepInstance(instance, type);
604  inner->addAttValue("pickParent",true);
605  inner->addAttValue("showParentAttributes",true);
606 
607  HepRepPoint* ip1 = factory->createHepRepPoint(inner, vertex1.x(), vertex1.y(), vertex1.z());
608  ip1->addAttValue("Radius",cons.GetInnerRadiusPlusZ());
609 
610  HepRepPoint* ip2 = factory->createHepRepPoint(inner, vertex2.x(), vertex2.y(), vertex2.z());
611  ip2->addAttValue("Radius",cons.GetInnerRadiusMinusZ());
612 }
613 
614 
616 #ifdef SDEBUG
617  cout << "G4HepRepSceneHandler::AddSolid(const G4Tubs& tubs)" << endl;
618 #endif
619 
620  if (dontWrite()) return;
621 
623 
624  if (! messenger->useSolids() || (tubs.GetDeltaPhiAngle() < twopi)) {
626  return;
627  }
628 
629  G4PhysicalVolumeModel* pPVModel =
630  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
631  if (!pPVModel) {
633  return;
634  }
635 
636  G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
637  G4int currentDepth = pPVModel->GetCurrentDepth();
638  G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
639 
640  G4Point3D vertex1(G4Point3D( 0., 0., tubs.GetZHalfLength()));
641  G4Point3D vertex2(G4Point3D( 0., 0.,-tubs.GetZHalfLength()));
642 
643  vertex1 = (transform) * vertex1;
644  vertex2 = (transform) * vertex2;
645 
646  HepRepInstance* instance = getGeometryInstance(pCurrentLV, pCurrentMaterial, currentDepth);
647  setAttribute(instance, "DrawAs", G4String("Cylinder"));
648 
649  setVisibility(instance, tubs);
650  setLine(instance, tubs);
651  setColor(instance, getColorFor(tubs));
652 
653  HepRepType* type = getGeometryType(pCurrentLV->GetName(), currentDepth);
654 
655  // Outer cylinder.
656  HepRepInstance* outer = factory->createHepRepInstance(instance, type);
657  outer->addAttValue("Radius",tubs.GetOuterRadius());
658  outer->addAttValue("pickParent",true);
659  outer->addAttValue("showParentAttributes",true);
660  factory->createHepRepPoint(outer, vertex1.x(), vertex1.y(), vertex1.z());
661  factory->createHepRepPoint(outer, vertex2.x(), vertex2.y(), vertex2.z());
662 
663  // Inner cylinder.
664  if (tubs.GetInnerRadius() > 0.) {
665  HepRepInstance* inner = factory->createHepRepInstance(instance, type);
666  inner->addAttValue("Radius",tubs.GetInnerRadius());
667  inner->addAttValue("pickParent",true);
668  inner->addAttValue("showParentAttributes",true);
669  factory->createHepRepPoint(inner, vertex1.x(), vertex1.y(), vertex1.z());
670  factory->createHepRepPoint(inner, vertex2.x(), vertex2.y(), vertex2.z());
671  }
672 }
673 
674 
676 #ifdef SDEBUG
677  cout << "G4HepRepSceneHandler::AddSolid(const G4Trd& trd)" << endl;
678 #endif
679  if (dontWrite()) return;
680 
682 
683  if (! messenger->useSolids()) {
685  return;
686  }
687 
688  G4double dx1 = trd.GetXHalfLength1();
689  G4double dy1 = trd.GetYHalfLength1();
690  G4double dx2 = trd.GetXHalfLength2();
691  G4double dy2 = trd.GetYHalfLength2();
692  G4double dz = trd.GetZHalfLength();
693 
694  G4Point3D vertex1(G4Point3D( dx1, dy1,-dz));
695  G4Point3D vertex2(G4Point3D( dx1,-dy1,-dz));
696  G4Point3D vertex3(G4Point3D(-dx1,-dy1,-dz));
697  G4Point3D vertex4(G4Point3D(-dx1, dy1,-dz));
698  G4Point3D vertex5(G4Point3D( dx2, dy2, dz));
699  G4Point3D vertex6(G4Point3D( dx2,-dy2, dz));
700  G4Point3D vertex7(G4Point3D(-dx2,-dy2, dz));
701  G4Point3D vertex8(G4Point3D(-dx2, dy2, dz));
702 
703  vertex1 = (transform) * vertex1;
704  vertex2 = (transform) * vertex2;
705  vertex3 = (transform) * vertex3;
706  vertex4 = (transform) * vertex4;
707  vertex5 = (transform) * vertex5;
708  vertex6 = (transform) * vertex6;
709  vertex7 = (transform) * vertex7;
710  vertex8 = (transform) * vertex8;
711 
712  HepRepInstance* instance = getGeometryOrEventInstance(getCalHitType());
713 
714  addAttributes(instance, getCalHitType());
715 
716  setAttribute(instance, "DrawAs", G4String("Prism"));
717 
718  setVisibility(instance, trd);
719  setLine(instance, trd);
720  setColor(instance, getColorFor(trd));
721 
722  factory->createHepRepPoint(instance, vertex1.x(), vertex1.y(), vertex1.z());
723  factory->createHepRepPoint(instance, vertex2.x(), vertex2.y(), vertex2.z());
724  factory->createHepRepPoint(instance, vertex3.x(), vertex3.y(), vertex3.z());
725  factory->createHepRepPoint(instance, vertex4.x(), vertex4.y(), vertex4.z());
726  factory->createHepRepPoint(instance, vertex5.x(), vertex5.y(), vertex5.z());
727  factory->createHepRepPoint(instance, vertex6.x(), vertex6.y(), vertex6.z());
728  factory->createHepRepPoint(instance, vertex7.x(), vertex7.y(), vertex7.z());
729  factory->createHepRepPoint(instance, vertex8.x(), vertex8.y(), vertex8.z());
730 }
731 
733  if (dontWrite()) return;
735 }
736 
738  if (dontWrite()) return;
739  G4VSceneHandler::AddSolid (sphere);
740 }
741 
743  if (dontWrite()) return;
745 }
746 
748  if (dontWrite()) return;
749  G4VSceneHandler::AddSolid (torus);
750 }
751 
753  if (dontWrite()) return;
754  G4VSceneHandler::AddSolid (polycone);
755 }
756 
757 void G4HepRepSceneHandler::AddSolid (const G4Polyhedra& polyhedra) {
758  if (dontWrite()) return;
759  G4VSceneHandler::AddSolid (polyhedra);
760 }
761 
763  if (dontWrite()) return;
765 }
766 
767 
769 
770 #ifdef PDEBUG
771  cout << "G4HepRepSceneHandler::AddPrimitive(G4Polyline&) " << line.size() << endl;
772 #endif
773  if (dontWrite()) return;
774 
775  if (fProcessing2D) {
776  static G4bool warned = false;
777  if (!warned) {
778  warned = true;
780  ("G4HepRepSceneHandler::AddPrimitive (const G4Polyline&)",
781  "vis-HepRep1001", JustWarning,
782  "2D polylines not implemented. Ignored.");
783  }
784  return;
785  }
786 
787  HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getTrajectoryType());
788 
789  addAttributes(instance, getTrajectoryType());
790 
791  setColor(instance, GetColor(line));
792 
793  setVisibility(instance, line);
794 
795  setLine(instance, line);
796 
797  for (size_t i=0; i < line.size(); i++) {
798  G4Point3D vertex = transform * line[i];
799  factory->createHepRepPoint(instance, vertex.x(), vertex.y(), vertex.z());
800  }
801 }
802 
803 
805 
806 #ifdef PDEBUG
807  cout << "G4HepRepSceneHandler::AddPrimitive(G4Polymarker&) " << line.size() << endl;
808 #endif
809  if (dontWrite()) return;
810 
811  if (fProcessing2D) {
812  static G4bool warned = false;
813  if (!warned) {
814  warned = true;
816  ("G4HepRepSceneHandler::AddPrimitive (const G4Polymarker&)",
817  "vis-HepRep1002", JustWarning,
818  "2D polymarkers not implemented. Ignored.");
819  }
820  return;
821  }
822 
823  HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getHitType());
824 
825  addAttributes(instance, getHitType());
826 
827  setColor(instance, GetColor(line));
828 
829  setVisibility(instance, line);
830 
831  setMarker(instance, line);
832 
833  // Default MarkName is set to Circle for this Type.
834  int mtype = line.GetMarkerType();
835 
836  // Cannot be case statement since line.xxx is not a constant
837  if (mtype == line.dots) {
838  setAttribute(instance, "Fill", true);
839  setColor(instance, GetColor(line), G4String("FillColor"));
840  } else if (mtype == line.circles) {
841  } else if (line.squares) {
842  setAttribute(instance, "MarkName", G4String("Box"));
843  } else {
844  // line.line + default
845  setAttribute(instance, "MarkName", G4String("Plus"));
846  }
847 
848  for (size_t i=0; i < line.size(); i++) {
849  G4Point3D vertex = transform * line[i];
850  factory->createHepRepPoint(instance, vertex.x(), vertex.y(), vertex.z());
851  }
852 }
853 
854 
856 #ifdef PDEBUG
857  cout << "G4HepRepSceneHandler::AddPrimitive(G4Circle&) " << endl;
858 #endif
859  if (dontWrite()) return;
860 
861  if (fProcessing2D) {
862  static G4bool warned = false;
863  if (!warned) {
864  warned = true;
866  ("G4HepRepSceneHandler::AddPrimitive (const G4Circle&)",
867  "vis-HepRep1003", JustWarning,
868  "2D circles not implemented. Ignored.");
869  }
870  return;
871  }
872 
873  HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getHitType());
874 
875  addAttributes(instance, getHitType());
876 
877  G4Point3D center = transform * circle.GetPosition();
878 
879  setColor (instance, GetColor(circle));
880 
881  setVisibility(instance, circle);
882 
883  setMarker(instance, circle);
884 
885  factory->createHepRepPoint(instance, center.x(), center.y(), center.z());
886 }
887 
888 
890 
891 #ifdef PDEBUG
892  cout << "G4HepRepSceneHandler::AddPrimitive(G4Polyhedron&) " << endl;
893 #endif
894  if (dontWrite()) return;
895 
896  if (fProcessing2D) {
897  static G4bool warned = false;
898  if (!warned) {
899  warned = true;
901  ("G4HepRepSceneHandler::AddPrimitive (const G4Polyhedron&)",
902  "vis-HepRep1004", JustWarning,
903  "2D polyhedra not implemented. Ignored.");
904  }
905  return;
906  }
907 
908  G4Normal3D surfaceNormal;
909  G4Point3D vertex;
910 
911  if (polyhedron.GetNoFacets()==0) return;
912 
913  HepRepInstance* instance = getGeometryOrEventInstance(getCalHitType());
914 
915  addAttributes(instance, getCalHitType());
916 
917  setVisibility(instance, polyhedron);
918 
919  G4int currentDepth = 0;
920  G4PhysicalVolumeModel* pPVModel =
921  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
922  if (pPVModel) currentDepth = pPVModel->GetCurrentDepth();
923 
924  G4bool notLastFace;
925  do {
926  HepRepInstance* face;
927  if (isEventData()) {
928  face = factory->createHepRepInstance(instance, getCalHitFaceType());
929  } else {
930  face = getGeometryInstance("*Face", currentDepth+1);
931  setAttribute(face, "PickParent", true);
932  setAttribute(face, "DrawAs", G4String("Polygon"));
933  }
934 
935  setLine(face, polyhedron);
936  setColor(face, GetColor(polyhedron));
937  if (isEventData()) setColor(face, GetColor(polyhedron), G4String("FillColor"));
938 
939  notLastFace = polyhedron.GetNextNormal (surfaceNormal);
940 
941  G4int edgeFlag = 1;
942  G4bool notLastEdge;
943  do {
944  notLastEdge = polyhedron.GetNextVertex (vertex, edgeFlag);
945  vertex = transform * vertex;
946  factory->createHepRepPoint(face, vertex.x(), vertex.y(), vertex.z());
947  } while (notLastEdge);
948  } while (notLastFace);
949 }
950 
951 
953 #ifdef PDEBUG
954  cout << "G4HepRepSceneHandler::AddPrimitive(G4Text&) " << endl;
955 #endif
956  if (dontWrite()) return;
957 
958  /*** You may need this
959  if (fProcessing2D) {
960  static G4bool warned = false;
961  if (!warned) {
962  warned = true;
963  G4Exception
964  ("G4HepRepSceneHandler::AddPrimitive (const G4Text&)",
965  "vis-HepRep1005", JustWarning,
966  "2D text not implemented. Ignored.");
967  }
968  return;
969  }
970  ***/
971 
972  cout << "G4HepRepSceneHandler::AddPrimitive G4Text : not yet implemented. " << endl;
973 }
974 
975 
977 #ifdef PDEBUG
978  cout << "G4HepRepSceneHandler::AddPrimitive(G4Square&) " << endl;
979 #endif
980  if (dontWrite()) return;
981 
982  if (fProcessing2D) {
983  static G4bool warned = false;
984  if (!warned) {
985  warned = true;
987  ("G4HepRepSceneHandler::AddPrimitive (const G4Square&)",
988  "vis-HepRep1006", JustWarning,
989  "2D squares not implemented. Ignored.");
990  }
991  return;
992  }
993 
994  HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getHitType());
995 
996  addAttributes(instance, getHitType());
997 
998  G4Point3D center = transform * square.GetPosition();
999 
1000  setColor (instance, getColorFor(square));
1001 
1002  setVisibility(instance, square);
1003 
1004  setMarker(instance, square);
1005 
1006  factory->createHepRepPoint(instance, center.x(), center.y(), center.z());
1007 }
1008 
1010 #ifdef PDEBUG
1011  cout << "G4HepRepSceneHandler::AddPrimitive(G4NURBS&) " << endl;
1012 #endif
1013  if (dontWrite()) return;
1014 
1015  /*** You may need this
1016  if (fProcessing2D) {
1017  static G4bool warned = false;
1018  if (!warned) {
1019  warned = true;
1020  G4Exception
1021  ("G4HepRepSceneHandler::AddPrimitive (const G4NURBS&)",
1022  "vis-HepRep1007", JustWarning,
1023  "2D NURBS not implemented. Ignored.");
1024  }
1025  return;
1026  }
1027  ***/
1028 
1029  cout << "G4HepRepSceneHandler::AddPrimitive G4NURBS : not yet implemented. " << endl;
1030 }
1031 
1033  if (dontWrite()) return;
1035 }
1036 
1038 #ifdef PDEBUG
1039  cout << "G4HepRepSceneHandler::AddCompound(G4VTrajectory&) " << endl;
1040 #endif
1041  if (dontWrite()) return;
1042 
1043  currentTrack = &trajectory;
1044  G4VSceneHandler::AddCompound(trajectory);
1045  currentTrack = NULL;
1046 }
1047 
1048 
1050 #ifdef PDEBUG
1051  cout << "G4HepRepSceneHandler::AddCompound(G4VHit&) " << endl;
1052 #endif
1053  if (dontWrite()) return;
1054 
1055  currentHit = &hit;
1057  currentHit = NULL;
1058 }
1059 
1060 void G4HepRepSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
1061  const G4VisAttributes& visAttribs) {
1062 
1063  G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
1064 
1065  transform = objectTransformation;
1066 #ifdef SDEBUG
1067  cout << "G4HepRepSceneHandler::PreAddSolid(G4Transform3D&, G4VisAttributes&)" << endl;
1068 #endif
1069 }
1070 
1071 
1073 #ifdef SDEBUG
1074  cout << "G4HepRepSceneHandler::PostAddSolid()" << endl;
1075 #endif
1077 }
1078 
1079 
1080 void G4HepRepSceneHandler::BeginPrimitives (const G4Transform3D& objectTransformation) {
1081 #ifdef SDEBUG
1082  cout << "G4HepRepSceneHandler::BeginPrimitives(G4Transform3D&)" << endl;
1083 #endif
1084 
1085  G4VSceneHandler::BeginPrimitives (objectTransformation);
1086  transform = objectTransformation;
1087 }
1088 
1089 
1091 #ifdef SDEBUG
1092  cout << "G4HepRepSceneHandler::EndPrimitives" << endl;
1093 #endif
1095 }
1096 
1097 
1098 G4bool G4HepRepSceneHandler::dontWrite() {
1100  return !( messenger->writeInvisibles() || (fpVisAttribs ? (bool)fpVisAttribs->IsVisible() : true));
1101 }
1102 
1103 void G4HepRepSceneHandler::setColor (HepRepAttribute *attribute,
1104  const G4Color& color,
1105  const G4String& key) {
1106 #ifdef CDEBUG
1107  cout << "G4HepRepSceneHandler::setColor : red : " << color.GetRed () <<
1108  " green : " << color.GetGreen () <<
1109  " blue : " << color.GetBlue () << endl;
1110 #endif
1111 
1112  setAttribute(attribute, key, color.GetRed(), color.GetGreen(), color.GetBlue(), color.GetAlpha());
1113 }
1114 
1115 G4Color G4HepRepSceneHandler::getColorFor (const G4VSolid& /* solid */) {
1116  return fpVisAttribs ? fpVisAttribs->GetColor() : GetColor(NULL);
1117 }
1118 
1119 G4Color G4HepRepSceneHandler::getColorFor (const G4Visible& visible) {
1120  return GetColor(visible);
1121 }
1122 
1123 void G4HepRepSceneHandler::setVisibility (HepRepAttribute *attribute, const G4VSolid& /* solid */) {
1124  setAttribute(attribute, "Visibility", (fpVisAttribs ? (bool)fpVisAttribs->IsVisible() : true));
1125 }
1126 
1127 void G4HepRepSceneHandler::setVisibility ( HepRepAttribute *attribute, const G4Visible& visible) {
1128  const G4VisAttributes* atts = visible.GetVisAttributes();
1129 
1130  setAttribute(attribute, "Visibility", (atts && (atts->IsVisible()==0)) ? false : true);
1131 }
1132 
1133 void G4HepRepSceneHandler::setLine (HepRepAttribute *attribute, const G4VSolid& /* solid*/) {
1134  setAttribute(attribute, "LineWidth", 1.0);
1135 }
1136 
1137 void G4HepRepSceneHandler::setLine (HepRepAttribute *attribute, const G4Visible& visible) {
1138  const G4VisAttributes* atts = visible.GetVisAttributes();
1139 
1140  setAttribute(attribute, "LineWidth", (atts != NULL) ? atts->GetLineWidth() : 1.0);
1141 
1142  if (atts != NULL) {
1143  switch (atts->GetLineStyle()) {
1145  setAttribute(attribute, "LineStyle", G4String("Dotted"));
1146  break;
1148  setAttribute(attribute, "LineStyle", G4String("Dashed"));
1149  break;
1151  default:
1152  break;
1153  }
1154  }
1155 }
1156 
1157 void G4HepRepSceneHandler::setMarker (HepRepAttribute *attribute, const G4VMarker& marker) {
1158  MarkerSizeType markerType;
1159  G4double size = GetMarkerRadius( marker , markerType );
1160 
1161  setAttribute(attribute, "MarkSize", size);
1162 
1163  if (markerType == screen) setAttribute(attribute, "MarkType", G4String("Symbol"));
1164  if (marker.GetFillStyle() == G4VMarker::noFill) {
1165  setAttribute(attribute, "Fill", false);
1166  } else {
1167  setColor(attribute, GetColor(marker), G4String("FillColor"));
1168  }
1169 }
1170 
1171 void G4HepRepSceneHandler::addAttributes(HepRepInstance* instance, HepRepType* type) {
1172  if (currentHit) {
1173  vector<G4AttValue>* hitAttValues = currentHit->CreateAttValues();
1174  const map<G4String,G4AttDef>* hitAttDefs = currentHit->GetAttDefs();
1175 
1176  addAttDefs(getHitType(), hitAttDefs);
1177 
1178  // these attValues are non-standard, so can only be added when we have the attDef.
1179  type->addAttValue("LVol", G4String(""));
1180  type->addAttValue("HitType", G4String(""));
1181  type->addAttValue("ID", -1);
1182  type->addAttValue("Column", -1);
1183  type->addAttValue("Row", -1);
1184  type->addAttValue("Energy", 0.0);
1185  type->addAttValue("Pos", G4String(""));
1186 
1187  addAttVals(instance, hitAttDefs, hitAttValues);
1188 
1189  delete hitAttValues;
1190 
1191  } else if (currentTrack) {
1192  vector<G4AttValue>* trajectoryAttValues = currentTrack->CreateAttValues();
1193  const map<G4String,G4AttDef>* trajectoryAttDefs = currentTrack->GetAttDefs();
1194 
1195  addAttDefs(type, trajectoryAttDefs);
1196 
1197  // these attValues are non-standard, so can only be added when we have the attDef.
1198  type->addAttValue("Ch", 0.0);
1199  type->addAttValue("Color", 1.0, 1.0, 1.0, 1.0);
1200  type->addAttValue("ID", -1);
1201  type->addAttValue("IMom", G4String(""));
1202  type->addAttValue("IMag", 0.0);
1203  type->addAttValue("PDG", -1);
1204  type->addAttValue("PN", G4String(""));
1205  type->addAttValue("PID", -1);
1206 
1207  addAttVals(instance, trajectoryAttDefs, trajectoryAttValues);
1208 
1209  delete trajectoryAttValues;
1210 
1211  }
1212 }
1213 
1214 void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, G4String value) {
1215  HepRepAttValue* attValue = attribute->getAttValue(name);
1216  if ((attValue == NULL) || (attValue->getString() != value)) {
1217  HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1218  if (point != NULL) {
1219  if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1220  attribute = point->getInstance();
1221  }
1222  }
1223 
1224  HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1225  if (instance != NULL) {
1226  // look for definition on type (node only)
1227  if (instance->getType()->getAttValueFromNode(name) == NULL) {
1228  attribute = instance->getType();
1229  }
1230  }
1231 
1232  attribute->addAttValue(name, value);
1233  }
1234 }
1235 
1236 void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, bool value) {
1237  HepRepAttValue* attValue = attribute->getAttValue(name);
1238  if ((attValue == NULL) || (attValue->getBoolean() != value)) {
1239  HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1240  if (point != NULL) {
1241  if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1242  attribute = point->getInstance();
1243  }
1244  }
1245 
1246  HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1247  if (instance != NULL) {
1248  // look for definition on type (node only)
1249  if (instance->getType()->getAttValueFromNode(name) == NULL) {
1250  attribute = instance->getType();
1251  }
1252  }
1253 
1254  attribute->addAttValue(name, value);
1255  }
1256 }
1257 
1258 void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, double value) {
1259  HepRepAttValue* attValue = attribute->getAttValue(name);
1260  if ((attValue == NULL) || (attValue->getDouble() != value)) {
1261  HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1262  if (point != NULL) {
1263  if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1264  attribute = point->getInstance();
1265  }
1266  }
1267 
1268  HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1269  if (instance != NULL) {
1270  // look for definition on type (node only)
1271  if (instance->getType()->getAttValueFromNode(name) == NULL) {
1272  attribute = instance->getType();
1273  }
1274  }
1275 
1276  attribute->addAttValue(name, value);
1277  }
1278 }
1279 
1280 void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, int value) {
1281  HepRepAttValue* attValue = attribute->getAttValue(name);
1282  if ((attValue == NULL) || (attValue->getInteger() != value)) {
1283  HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1284  if (point != NULL) {
1285  if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1286  attribute = point->getInstance();
1287  }
1288  }
1289 
1290  HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1291  if (instance != NULL) {
1292  // look for definition on type (node only)
1293  if (instance->getType()->getAttValueFromNode(name) == NULL) {
1294  attribute = instance->getType();
1295  }
1296  }
1297 
1298  attribute->addAttValue(name, value);
1299  }
1300 }
1301 
1302 void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, double red, double green, double blue, double alpha) {
1303  HepRepAttValue* attValue = attribute->getAttValue(name);
1304  vector<double> color;
1305  if (attValue != NULL) color = attValue->getColor();
1306  if ((color.size() == 0) ||
1307  (color[0] != red) ||
1308  (color[1] != green) ||
1309  (color[2] != blue) ||
1310  ((color.size() > 3) && (color[3] != alpha))) {
1311 
1312  HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1313  if (point != NULL) {
1314  if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1315  attribute = point->getInstance();
1316  }
1317  }
1318 
1319  HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1320  if (instance != NULL) {
1321  // look for definition on type (node only)
1322  if (instance->getType()->getAttValueFromNode(name) == NULL) {
1323  attribute = instance->getType();
1324  }
1325  }
1326 
1327  attribute->addAttValue(name, red, green, blue, alpha);
1328  }
1329 }
1330 
1331 void G4HepRepSceneHandler::addAttDefs(HepRepDefinition* definition, const map<G4String,G4AttDef>* attDefs) {
1332  if (attDefs == NULL) return;
1333 
1334  // Specify additional attribute definitions.
1335  map<G4String,G4AttDef>::const_iterator attDefIterator = attDefs->begin();
1336  while (attDefIterator != attDefs->end()) {
1337  definition->addAttDef(attDefIterator->first, attDefIterator->second.GetDesc(),
1338  attDefIterator->second.GetCategory(), attDefIterator->second.GetExtra());
1339  attDefIterator++;
1340  }
1341 }
1342 
1343 void G4HepRepSceneHandler::addAttVals(HepRepAttribute* attribute, const map<G4String,G4AttDef>* attDefs, vector<G4AttValue>* attValues) {
1344  if (attValues == NULL) return;
1345 
1346  // Copy the instance's G4AttValues to HepRepAttValues.
1347  for (vector<G4AttValue>::iterator attValIterator = attValues->begin(); attValIterator != attValues->end(); attValIterator++) {
1348  G4String name = attValIterator->GetName();
1349 
1350  HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1351  if ((name == "Pos") && (point != NULL)) {
1352  G4String pos = attValIterator->GetValue();
1353 // cout << "Pos* " << pos << endl;
1354  int is = 0;
1355  int in = 0;
1356  int im = 0;
1357  G4String unit;
1358  for (unsigned int i=0; i<pos.length(); i++) {
1359  if (pos[i] == ' ') {
1360  if (in == 0) {
1361  // first coordinate
1362  double factor = atof(pos.substr(is, i-is).c_str())/point->getX();
1363  im = (int)(std::log10(factor)+((factor < 1) ? -0.5 : 0.5));
1364 // cout << factor << ", " << im << endl;
1365  } else if (in == 3) {
1366  // unit
1367  unit = pos.substr(is, i-is);
1368  if (unit == G4String("mum")) {
1369  im += -6;
1370  } else if (unit == G4String("mm")) {
1371  im += -3;
1372  } else if (unit == G4String("cm")) {
1373  im += -2;
1374  } else if (unit == G4String("m")) {
1375  im += 0;
1376  } else if (unit == G4String("km")) {
1377  im += 3;
1378  } else {
1379  cerr << "HepRepSceneHandler: Unrecognized Unit: '" << unit << "'" << endl;
1380  }
1381  }
1382  is = i+1;
1383  in++;
1384  }
1385  }
1386  switch(im) {
1387  case -6:
1388  unit = G4String("mum");
1389  break;
1390  case -3:
1391  unit = G4String("mm");
1392  break;
1393  case -2:
1394  unit = G4String("cm");
1395  break;
1396  case 0:
1397  unit = G4String("m");
1398  break;
1399  case 3:
1400  unit = G4String("km");
1401  break;
1402  default:
1403  cerr << "HepRepSceneHandler: No valid unit found for im: " << im << endl;
1404  unit = G4String("*im");
1405  break;
1406  }
1407 // cout << "U: " << unit << endl;
1408  setAttribute(attribute, G4String("PointUnit"), unit);
1409  continue;
1410  }
1411 
1412  // NTP already in points being written
1413  if (name == "NTP") continue;
1414 
1415  // find type of attribute using def
1416  const map<G4String,G4AttDef>::const_iterator attDefIterator = attDefs->find(name);
1417  G4String type = attDefIterator->second.GetValueType();
1418 
1419  // set based on type
1420  if ((type == "G4double") || (type == "double")) {
1421  setAttribute(attribute, attValIterator->GetName(), atof(attValIterator->GetValue()));
1422  } else if ((type == "G4int") || (type == "int")) {
1423  setAttribute(attribute, attValIterator->GetName(), atoi(attValIterator->GetValue()));
1424  } else { // G4String, string and others
1425  setAttribute(attribute, attValIterator->GetName(), attValIterator->GetValue());
1426  }
1427  }
1428 }
1429 
1430 
1431 bool G4HepRepSceneHandler::isEventData () {
1432  G4PhysicalVolumeModel* pPVModel =
1433  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1434  return !pPVModel || fReadyForTransients || currentHit || currentTrack;
1435 }
1436 
1437 void G4HepRepSceneHandler::addTopLevelAttributes(HepRepType* type) {
1438 
1439  // Some non-standard attributes
1440  type->addAttDef( "Generator", "Generator of the file", "General", "");
1441  type->addAttValue("Generator", G4String("Geant4"));
1442 
1443  type->addAttDef( "GeneratorVersion", "Version of the Generator", "General", "");
1444  G4String versionString = G4Version;
1445  versionString = versionString.substr(1,versionString.size()-2);
1446  versionString = " Geant4 version " + versionString + " " + G4Date;
1447  type->addAttValue("GeneratorVersion", versionString);
1448 
1449  const G4ViewParameters parameters = GetCurrentViewer()->GetViewParameters();
1450  const G4Vector3D& viewPointDirection = parameters.GetViewpointDirection();
1451  type->addAttDef( "ViewTheta", "Theta of initial suggested viewpoint", "Draw", "rad");
1452  type->addAttValue("ViewTheta", viewPointDirection.theta());
1453 
1454  type->addAttDef( "ViewPhi", "Phi of initial suggested viewpoint", "Draw", "rad");
1455  type->addAttValue("ViewPhi", viewPointDirection.phi());
1456 
1457  type->addAttDef( "ViewScale", "Scale of initial suggested viewpoint", "Draw", "");
1458  type->addAttValue("ViewScale", parameters.GetZoomFactor());
1459 
1460 // FIXME, no way to set these
1461  type->addAttDef( "ViewTranslateX", "Translate in X of initial suggested viewpoint", "Draw", "");
1462  type->addAttValue("ViewTranslateX", 0.0);
1463 
1464  type->addAttDef( "ViewTranslateY", "Translate in Y of initial suggested viewpoint", "Draw", "");
1465  type->addAttValue("ViewTranslateY", 0.0);
1466 
1467  type->addAttDef( "ViewTranslateZ", "Translate in Z of initial suggested viewpoint", "Draw", "");
1468  type->addAttValue("ViewTranslateZ", 0.0);
1469 
1470  type->addAttDef( "PointUnit", "Length", "Physics", "");
1471  type->addAttValue("PointUnit", G4String("m"));
1472 
1474 
1475  type->addAttDef( "UseSolids", "Use HepRep Solids rather than Geant4 Primitives", "Draw", "");
1476  type->addAttValue("UseSolids", messenger->useSolids());
1477 
1478  type->addAttDef( "WriteInvisibles", "Write Invisible Objects", "Draw", "");
1479  type->addAttValue("WriteInvisibles", messenger->writeInvisibles());
1480 }
1481 
1482 
1483 HepRepInstance* G4HepRepSceneHandler::getGeometryOrEventInstance(HepRepType* type) {
1484  if (isEventData()) {
1485  return factory->createHepRepInstance(getEventInstance(), type);
1486  } else {
1487  G4PhysicalVolumeModel* pPVModel =
1488  dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1489  G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
1490  G4int currentDepth = pPVModel->GetCurrentDepth();
1491  G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
1492  return getGeometryInstance(pCurrentLV, pCurrentMaterial, currentDepth);
1493  }
1494 }
1495 
1496 HepRep* G4HepRepSceneHandler::getHepRep() {
1497  if (_heprep == NULL) {
1498  // Create the HepRep that holds the Trees.
1499  _heprep = factory->createHepRep();
1500  }
1501  return _heprep;
1502 }
1503 
1504 HepRep* G4HepRepSceneHandler::getHepRepGeometry() {
1505  if (_heprepGeometry == NULL) {
1506  // Create the HepRep that holds the Trees.
1507  _heprepGeometry = factory->createHepRep();
1508  }
1509  return _heprepGeometry;
1510 }
1511 
1512 HepRepInstanceTree* G4HepRepSceneHandler::getGeometryInstanceTree() {
1513  if (_geometryInstanceTree == NULL) {
1514  // Create the Geometry InstanceTree.
1515  _geometryInstanceTree = factory->createHepRepInstanceTree("G4GeometryData", "1.0", getGeometryTypeTree());
1516 
1518  if ( messenger->appendGeometry()) {
1519  getHepRep()->addInstanceTree(_geometryInstanceTree);
1520  } else {
1521  getHepRepGeometry()->addInstanceTree(_geometryInstanceTree);
1522  }
1523  }
1524  return _geometryInstanceTree;
1525 }
1526 
1527 HepRepInstance* G4HepRepSceneHandler::getGeometryRootInstance() {
1528  if (_geometryRootInstance == NULL) {
1529  // Create the top level Geometry Instance.
1530  _geometryRootInstance = factory->createHepRepInstance(getGeometryInstanceTree(), getGeometryRootType());
1531  }
1532  return _geometryRootInstance;
1533 }
1534 
1535 HepRepInstance* G4HepRepSceneHandler::getGeometryInstance(G4LogicalVolume* volume, G4Material* material, int depth) {
1536  HepRepInstance* instance = getGeometryInstance(volume->GetName(), depth);
1537 
1538  setAttribute(instance, "LVol", volume->GetName());
1539  G4Region* region = volume->GetRegion();
1540  G4String regionName = region? region->GetName(): G4String("No region");
1541  setAttribute(instance, "Region", regionName);
1542  setAttribute(instance, "RootRegion", volume->IsRootRegion());
1543  setAttribute(instance, "Solid", volume->GetSolid()->GetName());
1544  setAttribute(instance, "EType", volume->GetSolid()->GetEntityType());
1545  G4String matName = material? material->GetName(): G4String("No material");
1546  setAttribute(instance, "Material", matName );
1547  G4double matDensity = material? material->GetDensity(): 0.;
1548  setAttribute(instance, "Density", matDensity);
1549  G4double matRadlen = material? material->GetRadlen(): 0.;
1550  setAttribute(instance, "Radlen", matRadlen);
1551 
1552  G4State matState = material? material->GetState(): kStateUndefined;
1553  G4String state = materialState[matState];
1554  setAttribute(instance, "State", state);
1555 
1556  return instance;
1557 }
1558 
1559 HepRepInstance* G4HepRepSceneHandler::getGeometryInstance(G4String volumeName, int depth) {
1560  // no extra checks since these are done in the geometryType already
1561 
1562  // adjust depth, also pop the current instance
1563  while ((int)_geometryInstance.size() > depth) {
1564  _geometryInstance.pop_back();
1565  }
1566 
1567  // get parent
1568  HepRepInstance* parent = (_geometryInstance.empty()) ? getGeometryRootInstance() : _geometryInstance.back();
1569 
1570  // get type
1571  HepRepType* type = getGeometryType(volumeName, depth);
1572 
1573  // create instance
1574  HepRepInstance* instance = factory->createHepRepInstance(parent, type);
1575  _geometryInstance.push_back(instance);
1576 
1577  return instance;
1578 }
1579 
1580 HepRepTypeTree* G4HepRepSceneHandler::getGeometryTypeTree() {
1581  if (_geometryTypeTree == NULL) {
1582  // Create the Geometry TypeTree.
1583  HepRepTreeID* geometryTreeID = factory->createHepRepTreeID("G4GeometryTypes", "1.0");
1584  _geometryTypeTree = factory->createHepRepTypeTree(geometryTreeID);
1585 
1587  if ( messenger->appendGeometry()) {
1588  getHepRep()->addTypeTree(_geometryTypeTree);
1589  } else {
1590  getHepRepGeometry()->addTypeTree(_geometryTypeTree);
1591  }
1592  }
1593  return _geometryTypeTree;
1594 }
1595 
1596 HepRepType* G4HepRepSceneHandler::getGeometryRootType() {
1597  if (_geometryRootType == NULL) {
1598  // Create the top level Geometry Type.
1599  _geometryRootType = factory->createHepRepType(getGeometryTypeTree(), rootVolumeName);
1600  _geometryRootType->addAttValue("Layer", geometryLayer);
1601 
1602  // Add attdefs used by all geometry types.
1603  _geometryRootType->addAttDef ("LVol", "Logical Volume", "Physics","");
1604  _geometryRootType->addAttValue("LVol", G4String(""));
1605  _geometryRootType->addAttDef ("Region", "Cuts Region", "Physics","");
1606  _geometryRootType->addAttValue("Region", G4String(""));
1607  _geometryRootType->addAttDef ("RootRegion", "Root Region", "Physics","");
1608  _geometryRootType->addAttValue("RootRegion", false);
1609  _geometryRootType->addAttDef ("Solid", "Solid Name", "Physics","");
1610  _geometryRootType->addAttValue("Solid", G4String(""));
1611  _geometryRootType->addAttDef ("EType", "Entity Type", "Physics","");
1612  _geometryRootType->addAttValue("EType", G4String("G4Box"));
1613  _geometryRootType->addAttDef ("Material", "Material Name", "Physics","");
1614  _geometryRootType->addAttValue("Material", G4String("Air"));
1615  _geometryRootType->addAttDef ("Density", "Material Density", "Physics","");
1616  _geometryRootType->addAttValue("Density", 0.0);
1617  _geometryRootType->addAttDef ("State", "Material State", "Physics","");
1618  _geometryRootType->addAttValue("State", G4String("Gas"));
1619  _geometryRootType->addAttDef ("Radlen", "Material Radiation Length", "Physics","");
1620  _geometryRootType->addAttValue("Radlen", 0.0);
1621 
1622  // add defaults for Geometry
1623  _geometryRootType->addAttValue("Color", 0.8, 0.8, 0.8, 1.0);
1624  _geometryRootType->addAttValue("Visibility", true);
1625  _geometryRootType->addAttValue("FillColor", 0.8, 0.8, 0.8, 1.0);
1626  _geometryRootType->addAttValue("LineWidth", 1.0);
1627  _geometryRootType->addAttValue("DrawAs", G4String("Polygon"));
1628  _geometryRootType->addAttValue("PickParent", false);
1629  _geometryRootType->addAttValue("ShowParentAttributes", true);
1630 
1631  _geometryRootType->addAttValue("MarkSizeMultiplier", 4.0);
1632  _geometryRootType->addAttValue("LineWidthMultiplier", 1.0);
1633 
1634  addTopLevelAttributes(_geometryRootType);
1635 
1636  _geometryType["/"+_geometryRootType->getName()] = _geometryRootType;
1637  }
1638  return _geometryRootType;
1639 }
1640 
1641 HepRepType* G4HepRepSceneHandler::getGeometryType(G4String volumeName, int depth) {
1642  // make sure we have a root
1643  getGeometryRootType();
1644 
1645  // construct the full name for this volume
1646  G4String name = getFullTypeName(volumeName, depth);
1647 
1648  // lookup type and create if necessary
1649  HepRepType* type = _geometryType[name];
1650  if (type == NULL) {
1651  G4String parentName = getParentTypeName(depth);
1652  HepRepType* parentType = _geometryType[parentName];
1653  // HepRep uses hierarchical names
1654  type = factory->createHepRepType(parentType, volumeName);
1655  _geometryType[name] = type;
1656  }
1657  return type;
1658 }
1659 
1660 G4String G4HepRepSceneHandler::getFullTypeName(G4String volumeName, int depth) {
1661  // check for name depth
1662  if (depth > (int)_geometryTypeName.size()) {
1663  // there is a problem, book this type under problems
1664  G4String problem = "HierarchyProblem";
1665  if (_geometryType["/"+problem] == NULL) {
1666  // HepRep uses hierarchical names
1667  HepRepType* type = factory->createHepRepType(getGeometryRootType(), problem);
1668  _geometryType["/"+problem] = type;
1669  }
1670  return "/" + problem + "/" + volumeName;
1671  }
1672 
1673  // adjust name depth, also pop the current volumeName
1674  while ((int)_geometryTypeName.size() > depth) {
1675  _geometryTypeName.pop_back();
1676  }
1677 
1678  // construct full name and push it
1679  G4String name = (_geometryTypeName.empty()) ? G4String("/"+rootVolumeName) : _geometryTypeName.back();
1680  name = name + "/" + volumeName;
1681  _geometryTypeName.push_back(name);
1682  return name;
1683 }
1684 
1685 G4String G4HepRepSceneHandler::getParentTypeName(int depth) {
1686  return (depth >= 1) ? _geometryTypeName[depth-1] : G4String("/"+rootVolumeName);
1687 }
1688 
1689 HepRepInstanceTree* G4HepRepSceneHandler::getEventInstanceTree() {
1690  if (_eventInstanceTree == NULL) {
1691  // Create the Event InstanceTree.
1692  _eventInstanceTree = factory->createHepRepInstanceTree("G4EventData", "1.0", getEventTypeTree());
1693  getHepRep()->addInstanceTree(_eventInstanceTree);
1694  }
1695  return _eventInstanceTree;
1696 }
1697 
1698 HepRepInstance* G4HepRepSceneHandler::getEventInstance() {
1699  if (_eventInstance == NULL) {
1700  // Create the top level Event Instance.
1701  _eventInstance = factory->createHepRepInstance(getEventInstanceTree(), getEventType());
1702  }
1703  return _eventInstance;
1704 }
1705 
1706 HepRepTypeTree* G4HepRepSceneHandler::getEventTypeTree() {
1707  if (_eventTypeTree == NULL) {
1708  // Create the Event TypeTree.
1709  HepRepTreeID* eventTreeID = factory->createHepRepTreeID("G4EventTypes", "1.0");
1710  _eventTypeTree = factory->createHepRepTypeTree(eventTreeID);
1711  getHepRep()->addTypeTree(_eventTypeTree);
1712  }
1713 
1714  return _eventTypeTree;
1715 }
1716 
1717 HepRepType* G4HepRepSceneHandler::getEventType() {
1718  if (_eventType == NULL) {
1719  // Create the top level Event Type.
1720  _eventType = factory->createHepRepType(getEventTypeTree(), "Event");
1721  _eventType->addAttValue("Layer", eventLayer);
1722 
1723  // add defaults for Events
1724  _eventType->addAttValue("Visibility", true);
1725  _eventType->addAttValue("Color", 1.0, 1.0, 1.0, 1.0);
1726  _eventType->addAttValue("FillColor", 1.0, 1.0, 1.0, 1.0);
1727  _eventType->addAttValue("LineWidth", 1.0);
1728  _eventType->addAttValue("HasFrame", true);
1729  _eventType->addAttValue("PickParent", false);
1730  _eventType->addAttValue("ShowParentAttributes", false);
1731 
1732  _eventType->addAttValue("MarkSizeMultiplier", 4.0);
1733  _eventType->addAttValue("LineWidthMultiplier", 1.0);
1734 
1735  addTopLevelAttributes(_eventType);
1736  }
1737 
1738  return _eventType;
1739 }
1740 
1741 HepRepType* G4HepRepSceneHandler::getTrajectoryType() {
1742  if (_trajectoryType == NULL) {
1743  _trajectoryType = factory->createHepRepType(getEventType(), "Trajectory");
1744 
1745  _trajectoryType->addAttValue("Layer", trajectoryLayer);
1746  _trajectoryType->addAttValue("DrawAs", G4String("Line"));
1747 
1748  _trajectoryType->addAttValue("LineWidthMultiplier", 2.0);
1749 
1750  // attributes to draw the points of a track as markers.
1751  _trajectoryType->addAttValue("MarkName", G4String("Box"));
1752  _trajectoryType->addAttValue("MarkSize", 4);
1753  _trajectoryType->addAttValue("MarkType", G4String("Symbol"));
1754  _trajectoryType->addAttValue("Fill", true);
1755  }
1756  return _trajectoryType;
1757 }
1758 
1759 HepRepType* G4HepRepSceneHandler::getHitType() {
1760  if (_hitType == NULL) {
1761  _hitType = factory->createHepRepType(getEventType(), "Hit");
1762  _hitType->addAttValue("Layer", hitLayer);
1763  _hitType->addAttValue("DrawAs", G4String("Point"));
1764  _hitType->addAttValue("MarkName", G4String("Box"));
1765  _hitType->addAttValue("MarkSize", 4.0);
1766  _hitType->addAttValue("MarkType", G4String("Symbol"));
1767  _hitType->addAttValue("Fill", true);
1768  }
1769  return _hitType;
1770 }
1771 
1772 HepRepType* G4HepRepSceneHandler::getCalHitType() {
1773  if (_calHitType == NULL) {
1774  _calHitType = factory->createHepRepType(getEventType(), "CalHit");
1775  _calHitType->addAttValue("Layer", calHitLayer);
1776  _calHitType->addAttValue("Fill", true);
1777  _calHitType->addAttValue("DrawAs", G4String("Polygon"));
1778  }
1779  return _calHitType;
1780 }
1781 
1782 HepRepType* G4HepRepSceneHandler::getCalHitFaceType() {
1783  if (_calHitFaceType == NULL) {
1784  _calHitFaceType = factory->createHepRepType(getCalHitType(), "CalHitFace");
1785  _calHitFaceType->addAttValue("PickParent", true);
1786  }
1787  return _calHitFaceType;
1788 }
1789