source: trunk/source/visualization/HepRep/src/G4HepRepSceneHandler.cc @ 1098

Last change on this file since 1098 was 944, checked in by garnier, 15 years ago

mise a jour des tags

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