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

Last change on this file was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 60.7 KB
RevLine 
[834]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//
[1228]26// $Id: G4HepRepSceneHandler.cc,v 1.102 2009/11/23 05:42:28 perl Exp $
[1347]27// GEANT4 tag $Name: geant4-09-04-ref-00 $
[834]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"
[1228]47#include "G4HepRepMessenger.hh"
[834]48
49//G4
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
75#include "cheprep/XMLHepRepFactory.h"
76
77// This
78#include "G4HepRep.hh"
79#include "G4HepRepSceneHandler.hh"
80#include "G4HepRepViewer.hh"
81
82
83using namespace HEPREP;
84using namespace cheprep;
85using namespace std;
86
87G4int G4HepRepSceneHandler::sceneIdCount = 0;
88
89//#define LDEBUG 1
90//#define SDEBUG 1
91//#define PDEBUG 1
92
[1228]93G4HepRepSceneHandler::G4HepRepSceneHandler (G4VGraphicsSystem& system, const G4String& name)
[834]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
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;
[1228]249
250                        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
251
252            string suffix = messenger->getEventNumberSuffix();
[834]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
292void G4HepRepSceneHandler::openHepRep() {
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
319/**
320 * Returns true if the HepRep was (already) closed, false if the HepRep is still open
321 */
322bool G4HepRepSceneHandler::closeHepRep(bool final) {
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 {
[1228]336               
337                G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
[834]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) {
342            GetCurrentViewer()->DrawView();
343       
344            // couple geometry
[1228]345
346            if ( messenger->appendGeometry()) {
[834]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
[1228]382        if (! messenger->appendGeometry() && (_heprepGeometry != NULL)) {
[834]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
440void 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
457void G4HepRepSceneHandler::openFile(G4String name) {
458    out = new ofstream(name.c_str(), std::ios::out | std::ios::binary );
459    writer = factory->createHepRepWriter(out, writeZip, writeZip || writeGZ);
460}
461
462void G4HepRepSceneHandler::closeFile() {
463    writer->close();
464    delete writer;
465    writer = NULL;
466
467    delete out;
468    out = NULL;
469}
470
471void 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
480void G4HepRepSceneHandler::BeginModeling() {
481#ifdef SDEBUG
482    cout << "G4HepRepSceneHandler::BeginModeling() " << endl;
483#endif
484    G4VSceneHandler::BeginModeling();
485}
486
487
488void G4HepRepSceneHandler::EndModeling() {
489#ifdef SDEBUG
490    cout << "G4HepRepSceneHandler::EndModeling() " << endl;
491#endif
492    G4VSceneHandler::EndModeling();
493}
494
495void G4HepRepSceneHandler::AddSolid(const G4Box& box) {
496#ifdef SDEBUG
497    cout << "G4HepRepSceneHandler::AddSolid(const G4Box& box)" << endl;
498#endif
499
500    if (dontWrite()) return;
[1228]501       
502        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
[834]503
[1228]504    if (! messenger->useSolids()) {
[834]505        G4VSceneHandler::AddSolid(box);
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
551void G4HepRepSceneHandler::AddSolid(const G4Cons& cons) {
552#ifdef SDEBUG
553    cout << "G4HepRepSceneHandler::AddSolid(const G4Cons& cons)" << endl;
554#endif
555
556    if (dontWrite()) return;
[1228]557       
558        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
[834]559
[1228]560    if (! messenger->useSolids() || (cons.GetDeltaPhiAngle() < twopi)) {
[834]561        G4VSceneHandler::AddSolid(cons);
562        return;
563    }
564
565    G4PhysicalVolumeModel* pPVModel =
566      dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
567    if (!pPVModel) {
568      G4VSceneHandler::AddSolid(cons);
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
615void G4HepRepSceneHandler::AddSolid(const G4Tubs& tubs) {
616#ifdef SDEBUG
617    cout << "G4HepRepSceneHandler::AddSolid(const G4Tubs& tubs)" << endl;
618#endif
619
620    if (dontWrite()) return;
[1228]621       
622        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
[834]623
[1228]624    if (! messenger->useSolids() || (tubs.GetDeltaPhiAngle() < twopi)) {
[834]625        G4VSceneHandler::AddSolid(tubs);
626        return;
627    }
628   
629    G4PhysicalVolumeModel* pPVModel =
630      dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
631    if (!pPVModel) {
632      G4VSceneHandler::AddSolid(tubs);
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
675void G4HepRepSceneHandler::AddSolid(const G4Trd& trd) {
676#ifdef SDEBUG
677    cout << "G4HepRepSceneHandler::AddSolid(const G4Trd& trd)" << endl;
678#endif
679    if (dontWrite()) return;
[1228]680       
681        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
[834]682
[1228]683    if (! messenger->useSolids()) {
[834]684        G4VSceneHandler::AddSolid(trd);
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
732void G4HepRepSceneHandler::AddSolid (const G4Trap& trap) { 
733    if (dontWrite()) return;
734    G4VSceneHandler::AddSolid (trap); 
735}
736
737void G4HepRepSceneHandler::AddSolid (const G4Sphere& sphere) { 
738    if (dontWrite()) return;
739    G4VSceneHandler::AddSolid (sphere); 
740}
741
742void G4HepRepSceneHandler::AddSolid (const G4Para& para) { 
743    if (dontWrite()) return;
744    G4VSceneHandler::AddSolid (para); 
745}
746
747void G4HepRepSceneHandler::AddSolid (const G4Torus& torus) { 
748    if (dontWrite()) return;
749    G4VSceneHandler::AddSolid (torus); 
750}
751
752void G4HepRepSceneHandler::AddSolid (const G4Polycone& polycone) { 
753    if (dontWrite()) return;
754    G4VSceneHandler::AddSolid (polycone); 
755}
756
757void G4HepRepSceneHandler::AddSolid (const G4Polyhedra& polyhedra) { 
758    if (dontWrite()) return;
759    G4VSceneHandler::AddSolid (polyhedra); 
760}
761
762void G4HepRepSceneHandler::AddSolid (const G4VSolid& solid) { 
763    if (dontWrite()) return;
764    G4VSceneHandler::AddSolid(solid); 
765}
766
767
768void G4HepRepSceneHandler::AddPrimitive (const G4Polyline& line) {
769
770#ifdef PDEBUG
771    cout << "G4HepRepSceneHandler::AddPrimitive(G4Polyline&) " << line.size() << endl;
772#endif
773    if (dontWrite()) return;
774
775    HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getTrajectoryType());
776
777    addAttributes(instance, getTrajectoryType());
778
779    setColor(instance, GetColor(line));
780
781    setVisibility(instance, line);
782
783    setLine(instance, line);
784
785    for (size_t i=0; i < line.size(); i++) {
786        G4Point3D vertex = transform * line[i];
787        factory->createHepRepPoint(instance, vertex.x(), vertex.y(), vertex.z());
788    }
789}
790
791
792void G4HepRepSceneHandler::AddPrimitive (const G4Polymarker& line) {
793
794#ifdef PDEBUG
795    cout << "G4HepRepSceneHandler::AddPrimitive(G4Polymarker&) " << line.size() << endl;
796#endif
797    if (dontWrite()) return;
798
799    HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getHitType());
800
801    addAttributes(instance, getHitType());
802
803    setColor(instance, GetColor(line));
804
805    setVisibility(instance, line);
806
807    setMarker(instance, line);
808
809    // Default MarkName is set to Circle for this Type.
810    int mtype = line.GetMarkerType();
811   
812    // Cannot be case statement since line.xxx is not a constant
813    if (mtype == line.dots) {
814        setAttribute(instance, "Fill", true);
815        setColor(instance, GetColor(line), G4String("FillColor"));
816    } else if (mtype == line.circles) {
817    } else if (line.squares) {
818        setAttribute(instance, "MarkName", G4String("Box"));
819    } else {
820        // line.line + default
821        setAttribute(instance, "MarkName", G4String("Plus"));
822    }
823
824    for (size_t i=0; i < line.size(); i++) {
825        G4Point3D vertex = transform * line[i];
826        factory->createHepRepPoint(instance, vertex.x(), vertex.y(), vertex.z());
827    }
828}
829
830
831void G4HepRepSceneHandler::AddPrimitive (const G4Circle& circle) {
832#ifdef PDEBUG
833    cout << "G4HepRepSceneHandler::AddPrimitive(G4Circle&) " << endl;
834#endif
835    if (dontWrite()) return;
836
837    HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getHitType());
838
839    addAttributes(instance, getHitType());
840
841    G4Point3D center = transform * circle.GetPosition();
842
843    setColor (instance, GetColor(circle));
844
845    setVisibility(instance, circle);
846
847    setMarker(instance, circle);
848
849    factory->createHepRepPoint(instance, center.x(), center.y(), center.z());
850}
851
852
853void G4HepRepSceneHandler::AddPrimitive (const G4Polyhedron& polyhedron) {
854
855#ifdef PDEBUG
856    cout << "G4HepRepSceneHandler::AddPrimitive(G4Polyhedron&) " << endl;
857#endif
858    if (dontWrite()) return;
859
860    G4Normal3D surfaceNormal;
861    G4Point3D vertex;
862
863    if (polyhedron.GetNoFacets()==0) return;
864
865    HepRepInstance* instance = getGeometryOrEventInstance(getCalHitType());
866       
867    addAttributes(instance, getCalHitType());
868
869    setVisibility(instance, polyhedron);
870       
871    G4int currentDepth = 0;
872    G4PhysicalVolumeModel* pPVModel =
873      dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
874    if (pPVModel) currentDepth = pPVModel->GetCurrentDepth();
875
876    G4bool notLastFace;
877    do {
878        HepRepInstance* face;
879        if (isEventData()) {
880            face = factory->createHepRepInstance(instance, getCalHitFaceType());
881        } else {
882            face = getGeometryInstance("*Face", currentDepth+1);
883            setAttribute(face, "PickParent", true);
884            setAttribute(face, "DrawAs", G4String("Polygon"));
885        }
886       
887        setLine(face, polyhedron);
888        setColor(face, GetColor(polyhedron));
889        if (isEventData()) setColor(face, GetColor(polyhedron), G4String("FillColor"));
890
891        notLastFace = polyhedron.GetNextNormal (surfaceNormal);
892
893        G4int edgeFlag = 1;
894        G4bool notLastEdge;
895        do {
896                notLastEdge = polyhedron.GetNextVertex (vertex, edgeFlag);
897            vertex = transform * vertex;
898            factory->createHepRepPoint(face, vertex.x(), vertex.y(), vertex.z());
899        } while (notLastEdge);
900    } while (notLastFace);
901}
902
903
904void G4HepRepSceneHandler::AddPrimitive (const G4Text&) {
905#ifdef PDEBUG
906    cout << "G4HepRepSceneHandler::AddPrimitive(G4Text&) " << endl;
907#endif
908    if (dontWrite()) return;
909
910    cout << "G4HepRepSceneHandler::AddPrimitive G4Text : not yet implemented. " << endl;
911}
912
913
914void G4HepRepSceneHandler::AddPrimitive (const G4Square& square) {
915#ifdef PDEBUG
916    cout << "G4HepRepSceneHandler::AddPrimitive(G4Square&) " << endl;
917#endif
918    if (dontWrite()) return;
919
920    HepRepInstance* instance = factory->createHepRepInstance(getEventInstance(), getHitType());
921
922    addAttributes(instance, getHitType());
923
924    G4Point3D center = transform * square.GetPosition();
925
926    setColor (instance, getColorFor(square));
927
928    setVisibility(instance, square);
929
930    setMarker(instance, square);
931
932    factory->createHepRepPoint(instance, center.x(), center.y(), center.z());
933}
934
935void G4HepRepSceneHandler::AddPrimitive (const G4NURBS&) {
936#ifdef PDEBUG
937    cout << "G4HepRepSceneHandler::AddPrimitive(G4NURBS&) " << endl;
938#endif
939    if (dontWrite()) return;
940
941    cout << "G4HepRepSceneHandler::AddPrimitive G4NURBS : not yet implemented. " << endl;
942}
943
944void G4HepRepSceneHandler::AddPrimitive (const G4Scale& scale) { 
945    if (dontWrite()) return;
946    G4VSceneHandler::AddPrimitive(scale); 
947}
948
949void G4HepRepSceneHandler::AddCompound (const G4VTrajectory& trajectory) {
950#ifdef PDEBUG
951    cout << "G4HepRepSceneHandler::AddCompound(G4VTrajectory&) " << endl;
952#endif
953    if (dontWrite()) return;
954   
955    currentTrack = &trajectory;
956    G4VSceneHandler::AddCompound(trajectory); 
957    currentTrack = NULL;   
958}
959
960
961void G4HepRepSceneHandler::AddCompound (const G4VHit& hit) { 
962#ifdef PDEBUG
963    cout << "G4HepRepSceneHandler::AddCompound(G4VHit&) " << endl;
964#endif
965    if (dontWrite()) return;
966   
967    currentHit = &hit;
968    G4VSceneHandler::AddCompound(hit); 
969    currentHit = NULL;
970}
971
972void G4HepRepSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
973                                  const G4VisAttributes& visAttribs) {
974
975    G4VSceneHandler::PreAddSolid (objectTransformation, visAttribs);
976
977    transform = objectTransformation;
978#ifdef SDEBUG
979    cout << "G4HepRepSceneHandler::PreAddSolid(G4Transform3D&, G4VisAttributes&)" << endl;
980#endif
981}
982
983
984void G4HepRepSceneHandler::PostAddSolid () {
985#ifdef SDEBUG
986    cout << "G4HepRepSceneHandler::PostAddSolid()" << endl;
987#endif
988    G4VSceneHandler::PostAddSolid();
989}
990
991
992void G4HepRepSceneHandler::BeginPrimitives (const G4Transform3D& objectTransformation) {
993#ifdef SDEBUG
994    cout << "G4HepRepSceneHandler::BeginPrimitives(G4Transform3D&)" << endl;
995#endif
996
997    G4VSceneHandler::BeginPrimitives (objectTransformation);
998    transform = objectTransformation;
999}
1000
1001
1002void G4HepRepSceneHandler::EndPrimitives () {
1003#ifdef SDEBUG
1004    cout << "G4HepRepSceneHandler::EndPrimitives" << endl;
1005#endif
1006    G4VSceneHandler::EndPrimitives ();
1007}
1008
1009
[1228]1010G4bool G4HepRepSceneHandler::dontWrite() {     
1011        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1012    return !( messenger->writeInvisibles() || (fpVisAttribs ? (bool)fpVisAttribs->IsVisible() : true));
[834]1013}
1014
1015void G4HepRepSceneHandler::setColor (HepRepAttribute *attribute,
1016                                      const G4Color& color,
1017                                      const G4String& key) {
1018#ifdef CDEBUG
1019    cout << "G4HepRepSceneHandler::setColor : red : " << color.GetRed ()   <<
1020                                  " green : " << color.GetGreen () <<
1021                                  " blue : " << color.GetBlue ()   << endl;
1022#endif
1023
1024    setAttribute(attribute, key, color.GetRed(), color.GetGreen(), color.GetBlue(), color.GetAlpha());
1025}
1026
1027G4Color G4HepRepSceneHandler::getColorFor (const G4VSolid& /* solid */) {
1028    return fpVisAttribs ? fpVisAttribs->GetColor() : GetColor(NULL);   
1029}
1030
1031G4Color G4HepRepSceneHandler::getColorFor (const G4Visible& visible) {
1032    return GetColor(visible);
1033}
1034
1035void G4HepRepSceneHandler::setVisibility (HepRepAttribute *attribute, const G4VSolid& /* solid */) {
1036    setAttribute(attribute, "Visibility", (fpVisAttribs ? (bool)fpVisAttribs->IsVisible() : true));
1037}
1038
1039void G4HepRepSceneHandler::setVisibility ( HepRepAttribute *attribute, const G4Visible& visible) {
1040    const G4VisAttributes* atts = visible.GetVisAttributes();
1041
1042    setAttribute(attribute, "Visibility", (atts && (atts->IsVisible()==0)) ? false : true);
1043}
1044
1045void G4HepRepSceneHandler::setLine (HepRepAttribute *attribute, const G4VSolid& /* solid*/) {
1046    setAttribute(attribute, "LineWidth", 1.0);
1047}
1048
1049void G4HepRepSceneHandler::setLine (HepRepAttribute *attribute, const G4Visible& visible) {
1050    const G4VisAttributes* atts = visible.GetVisAttributes();
1051   
1052    setAttribute(attribute, "LineWidth", (atts != NULL) ? atts->GetLineWidth() : 1.0);
1053   
1054    if (atts != NULL) {
1055        switch (atts->GetLineStyle()) {
1056            case G4VisAttributes::dotted:
1057                setAttribute(attribute, "LineStyle", G4String("Dotted"));
1058                break;
1059            case G4VisAttributes::dashed:
1060                setAttribute(attribute, "LineStyle", G4String("Dashed"));
1061                break;
1062            case G4VisAttributes::unbroken:
1063            default:
1064                break;
1065        }
1066    }
1067}
1068
1069void G4HepRepSceneHandler::setMarker (HepRepAttribute *attribute, const G4VMarker& marker) {
1070    MarkerSizeType markerType;
1071    G4double size = GetMarkerRadius( marker , markerType );
1072
1073    setAttribute(attribute, "MarkSize", size);
1074
1075    if (markerType == screen) setAttribute(attribute, "MarkType", G4String("Symbol"));
1076    if (marker.GetFillStyle() == G4VMarker::noFill) {
1077        setAttribute(attribute, "Fill", false);
1078    } else {
1079        setColor(attribute, GetColor(marker), G4String("FillColor"));
1080    }
1081}
1082
1083void G4HepRepSceneHandler::addAttributes(HepRepInstance* instance, HepRepType* type) {
1084    if (currentHit) {
1085        vector<G4AttValue>* hitAttValues = currentHit->CreateAttValues();
1086        const map<G4String,G4AttDef>* hitAttDefs = currentHit->GetAttDefs();
1087
1088        addAttDefs(getHitType(), hitAttDefs);
1089
1090        // these attValues are non-standard, so can only be added when we have the attDef.
1091        type->addAttValue("LVol", G4String(""));
1092        type->addAttValue("HitType", G4String(""));
1093        type->addAttValue("ID", -1);
1094        type->addAttValue("Column", -1);
1095        type->addAttValue("Row", -1);
1096        type->addAttValue("Energy", 0.0);
1097        type->addAttValue("Pos", G4String(""));
1098
1099        addAttVals(instance, hitAttDefs, hitAttValues);
1100   
1101        delete hitAttValues;
1102   
1103    } else if (currentTrack) {
1104        vector<G4AttValue>* trajectoryAttValues = currentTrack->CreateAttValues();
1105        const map<G4String,G4AttDef>* trajectoryAttDefs = currentTrack->GetAttDefs();
1106
1107        addAttDefs(type, trajectoryAttDefs);
1108   
1109        // these attValues are non-standard, so can only be added when we have the attDef.
1110        type->addAttValue("Ch", 0.0);
1111        type->addAttValue("Color", 1.0, 1.0, 1.0, 1.0);
1112        type->addAttValue("ID", -1);
1113        type->addAttValue("IMom", G4String(""));
1114        type->addAttValue("IMag", 0.0);
1115        type->addAttValue("PDG", -1);
1116        type->addAttValue("PN", G4String(""));
1117        type->addAttValue("PID", -1);
1118
1119        addAttVals(instance, trajectoryAttDefs, trajectoryAttValues);
1120   
1121        delete trajectoryAttValues;
1122       
1123    }
1124}
1125
1126void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, G4String value) {
1127    HepRepAttValue* attValue = attribute->getAttValue(name);
1128    if ((attValue == NULL) || (attValue->getString() != value)) {
1129        HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1130        if (point != NULL) {
1131            if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1132                attribute = point->getInstance();
1133            }
1134        }
1135       
1136        HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1137        if (instance != NULL) {
1138            // look for definition on type (node only)
1139            if (instance->getType()->getAttValueFromNode(name) == NULL) {
1140                attribute = instance->getType();
1141            }   
1142        }
1143       
1144        attribute->addAttValue(name, value);
1145    }
1146}
1147
1148void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, bool value) {
1149    HepRepAttValue* attValue = attribute->getAttValue(name);
1150    if ((attValue == NULL) || (attValue->getBoolean() != value)) {
1151        HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1152        if (point != NULL) {
1153            if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1154                attribute = point->getInstance();
1155            }
1156        }
1157       
1158        HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1159        if (instance != NULL) {
1160            // look for definition on type (node only)
1161            if (instance->getType()->getAttValueFromNode(name) == NULL) {
1162                attribute = instance->getType();
1163            }   
1164        }
1165       
1166        attribute->addAttValue(name, value);
1167    }
1168}
1169
1170void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, double value) {
1171    HepRepAttValue* attValue = attribute->getAttValue(name);
1172    if ((attValue == NULL) || (attValue->getDouble() != value)) {
1173        HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1174        if (point != NULL) {
1175            if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1176                attribute = point->getInstance();
1177            }
1178        }
1179       
1180        HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1181        if (instance != NULL) {
1182            // look for definition on type (node only)
1183            if (instance->getType()->getAttValueFromNode(name) == NULL) {
1184                attribute = instance->getType();
1185            }   
1186        }
1187       
1188        attribute->addAttValue(name, value);
1189    }
1190}
1191
1192void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, int value) {
1193    HepRepAttValue* attValue = attribute->getAttValue(name);
1194    if ((attValue == NULL) || (attValue->getInteger() != value)) {
1195        HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1196        if (point != NULL) {
1197            if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1198                attribute = point->getInstance();
1199            }
1200        }
1201       
1202        HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1203        if (instance != NULL) {
1204            // look for definition on type (node only)
1205            if (instance->getType()->getAttValueFromNode(name) == NULL) {
1206                attribute = instance->getType();
1207            }   
1208        }
1209       
1210        attribute->addAttValue(name, value);
1211    }
1212}
1213
1214void G4HepRepSceneHandler::setAttribute(HepRepAttribute* attribute, G4String name, double red, double green, double blue, double alpha) {
1215    HepRepAttValue* attValue = attribute->getAttValue(name);
1216    vector<double> color;
1217    if (attValue != NULL) color = attValue->getColor();
1218    if ((color.size() == 0) ||
1219        (color[0] != red) ||
1220        (color[1] != green) ||
1221        (color[2] != blue) ||
1222        ((color.size() > 3) && (color[3] != alpha))) {
1223       
1224        HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1225        if (point != NULL) {
1226            if (point->getInstance()->getAttValueFromNode(name) == NULL) {
1227                attribute = point->getInstance();
1228            }
1229        }
1230       
1231        HepRepInstance* instance = dynamic_cast<HepRepInstance*>(attribute);
1232        if (instance != NULL) {
1233            // look for definition on type (node only)
1234            if (instance->getType()->getAttValueFromNode(name) == NULL) {
1235                attribute = instance->getType();
1236            }   
1237        }
1238       
1239        attribute->addAttValue(name, red, green, blue, alpha);
1240    }
1241}
1242
1243void G4HepRepSceneHandler::addAttDefs(HepRepDefinition* definition, const map<G4String,G4AttDef>* attDefs) {
1244    if (attDefs == NULL) return;
1245
1246    // Specify additional attribute definitions.
1247    map<G4String,G4AttDef>::const_iterator attDefIterator = attDefs->begin();
1248    while (attDefIterator != attDefs->end()) {
1249        definition->addAttDef(attDefIterator->first, attDefIterator->second.GetDesc(),
1250                        attDefIterator->second.GetCategory(), attDefIterator->second.GetExtra());
1251        attDefIterator++;
1252    }
1253}
1254
1255void G4HepRepSceneHandler::addAttVals(HepRepAttribute* attribute, const map<G4String,G4AttDef>* attDefs, vector<G4AttValue>* attValues) {
1256    if (attValues == NULL) return;
1257
1258    // Copy the instance's G4AttValues to HepRepAttValues.
1259    for (vector<G4AttValue>::iterator attValIterator = attValues->begin(); attValIterator != attValues->end(); attValIterator++) {       
1260        G4String name = attValIterator->GetName();
1261
1262        HepRepPoint* point = dynamic_cast<HepRepPoint*>(attribute);
1263        if ((name == "Pos") && (point != NULL)) {
1264            G4String pos = attValIterator->GetValue();
1265//            cout << "Pos* " << pos << endl;
1266            int s = 0;
1267            int n = 0;
1268            int m = 0;
1269            G4String unit;
1270            for (unsigned int i=0; i<pos.length(); i++) {
1271                if (pos[i] == ' ') {
1272                    if (n == 0) {
1273                        // first coordinate
1274                        double factor = atof(pos.substr(s, i-s).c_str())/point->getX();
1275                        m = (int)(std::log10(factor)+((factor < 1) ? -0.5 : 0.5));
1276//                        cout << factor << ", " << m << endl;
1277                    } else if (n == 3) {
1278                        // unit
1279                        unit = pos.substr(s, i-s);
1280                        if (unit == G4String("mum")) {
1281                            m += -6;
1282                        } else if (unit == G4String("mm")) {
1283                            m += -3;
1284                        } else if (unit == G4String("cm")) {
1285                            m += -2;
1286                        } else if (unit == G4String("m")) {
1287                            m += 0;
1288                        } else if (unit == G4String("km")) {
1289                            m += 3;
1290                        } else {
1291                            cerr << "HepRepSceneHandler: Unrecognized Unit: '" << unit << "'" << endl;
1292                        }
1293                    }
1294                    s = i+1;
1295                    n++;
1296                }
1297            }
1298            switch(m) {
1299                case -6:
1300                    unit = G4String("mum");
1301                    break;
1302                case -3:
1303                    unit = G4String("mm");
1304                    break;
1305                case -2:
1306                    unit = G4String("cm");
1307                    break;
1308                case 0:
1309                    unit = G4String("m");
1310                    break;
1311                case 3:
1312                    unit = G4String("km");
1313                    break;
1314                default:
1315                    cerr << "HepRepSceneHandler: No valid unit found for m: " << m << endl;
1316                    unit = G4String("*m");                       
1317                    break;
1318            }
1319//            cout << "U: " << unit << endl;
1320            setAttribute(attribute, G4String("PointUnit"), unit);
1321            continue;
1322        }
1323               
1324        // NTP already in points being written
1325        if (name == "NTP") continue;
1326       
1327        // find type of attribute using def
1328        const map<G4String,G4AttDef>::const_iterator attDefIterator = attDefs->find(name);
1329        G4String type = attDefIterator->second.GetValueType();
1330               
1331        // set based on type
1332        if ((type == "G4double") || (type == "double")) {   
1333            setAttribute(attribute, attValIterator->GetName(), atof(attValIterator->GetValue()));           
1334        } else if ((type == "G4int") || (type == "int")) { 
1335            setAttribute(attribute, attValIterator->GetName(), atoi(attValIterator->GetValue()));           
1336        } else { // G4String, string and others
1337            setAttribute(attribute, attValIterator->GetName(), attValIterator->GetValue());
1338        }
1339    }
1340}
1341
1342
1343bool G4HepRepSceneHandler::isEventData () {
1344    G4PhysicalVolumeModel* pPVModel =
1345      dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1346    return !pPVModel || fReadyForTransients || currentHit || currentTrack;
1347}
1348
1349void G4HepRepSceneHandler::addTopLevelAttributes(HepRepType* type) {
1350   
1351    // Some non-standard attributes
1352    type->addAttDef(  "Generator", "Generator of the file", "General", "");
1353    type->addAttValue("Generator", G4String("Geant4"));
1354
1355    type->addAttDef(  "GeneratorVersion", "Version of the Generator", "General", "");
1356    G4String versionString = G4Version;
1357    versionString = versionString.substr(1,versionString.size()-2);
1358    versionString = " Geant4 version " + versionString + "   " + G4Date;
1359    type->addAttValue("GeneratorVersion", versionString);
1360   
1361    const G4ViewParameters parameters = GetCurrentViewer()->GetViewParameters();
1362    const G4Vector3D& viewPointDirection = parameters.GetViewpointDirection();   
1363    type->addAttDef(  "ViewTheta", "Theta of initial suggested viewpoint", "Draw", "rad");
1364    type->addAttValue("ViewTheta", viewPointDirection.theta());   
1365   
1366    type->addAttDef(  "ViewPhi", "Phi of initial suggested viewpoint", "Draw", "rad");
1367    type->addAttValue("ViewPhi", viewPointDirection.phi());   
1368   
1369    type->addAttDef(  "ViewScale", "Scale of initial suggested viewpoint", "Draw", "");
1370    type->addAttValue("ViewScale", parameters.GetZoomFactor());   
1371   
1372// FIXME, no way to set these
1373    type->addAttDef(  "ViewTranslateX", "Translate in X of initial suggested viewpoint", "Draw", "");
1374    type->addAttValue("ViewTranslateX", 0.0);   
1375   
1376    type->addAttDef(  "ViewTranslateY", "Translate in Y of initial suggested viewpoint", "Draw", "");
1377    type->addAttValue("ViewTranslateY", 0.0);   
1378   
1379    type->addAttDef(  "ViewTranslateZ", "Translate in Z of initial suggested viewpoint", "Draw", "");
1380    type->addAttValue("ViewTranslateZ", 0.0);   
1381   
1382    type->addAttDef(  "PointUnit", "Length", "Physics", "");
1383    type->addAttValue("PointUnit", G4String("m"));
[1228]1384       
1385        G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
[834]1386
1387    type->addAttDef(  "UseSolids", "Use HepRep Solids rather than Geant4 Primitives", "Draw", "");
[1228]1388    type->addAttValue("UseSolids", messenger->useSolids());
[834]1389
1390    type->addAttDef(  "WriteInvisibles", "Write Invisible Objects", "Draw", "");
[1228]1391    type->addAttValue("WriteInvisibles", messenger->writeInvisibles());
[834]1392}           
1393
1394
1395HepRepInstance* G4HepRepSceneHandler::getGeometryOrEventInstance(HepRepType* type) {
1396    if (isEventData()) {
1397      return factory->createHepRepInstance(getEventInstance(), type);
1398    } else {
1399      G4PhysicalVolumeModel* pPVModel =
1400        dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
1401      G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
1402      G4int currentDepth = pPVModel->GetCurrentDepth();
1403      G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
1404      return getGeometryInstance(pCurrentLV, pCurrentMaterial, currentDepth);
1405    }
1406}
1407
1408HepRep* G4HepRepSceneHandler::getHepRep() {
1409    if (_heprep == NULL) {
1410        // Create the HepRep that holds the Trees.
1411        _heprep = factory->createHepRep();
1412    }
1413    return _heprep;
1414}   
1415
1416HepRep* G4HepRepSceneHandler::getHepRepGeometry() {
1417    if (_heprepGeometry == NULL) {
1418        // Create the HepRep that holds the Trees.
1419        _heprepGeometry = factory->createHepRep();
1420    }
1421    return _heprepGeometry;
1422}   
1423
1424HepRepInstanceTree* G4HepRepSceneHandler::getGeometryInstanceTree() {
1425    if (_geometryInstanceTree == NULL) {
1426        // Create the Geometry InstanceTree.
1427        _geometryInstanceTree = factory->createHepRepInstanceTree("G4GeometryData", "1.0", getGeometryTypeTree());
[1228]1428               
1429                G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1430        if ( messenger->appendGeometry()) {
[834]1431            getHepRep()->addInstanceTree(_geometryInstanceTree);
1432        } else {
1433            getHepRepGeometry()->addInstanceTree(_geometryInstanceTree);
1434        }
1435    }
1436    return _geometryInstanceTree;
1437}
1438
1439HepRepInstance* G4HepRepSceneHandler::getGeometryRootInstance() {
1440    if (_geometryRootInstance == NULL) {
1441        // Create the top level Geometry Instance.
1442        _geometryRootInstance = factory->createHepRepInstance(getGeometryInstanceTree(), getGeometryRootType());
1443    }
1444    return _geometryRootInstance;
1445}
1446
1447HepRepInstance* G4HepRepSceneHandler::getGeometryInstance(G4LogicalVolume* volume, G4Material* material, int depth) {
1448    HepRepInstance* instance = getGeometryInstance(volume->GetName(), depth);
1449
1450    setAttribute(instance, "LVol",       volume->GetName());
1451    G4Region* region = volume->GetRegion();
1452    G4String regionName = region? region->GetName(): G4String("No region");
1453    setAttribute(instance, "Region",     regionName);
1454    setAttribute(instance, "RootRegion", volume->IsRootRegion());
1455    setAttribute(instance, "Solid",      volume->GetSolid()->GetName());
1456    setAttribute(instance, "EType",      volume->GetSolid()->GetEntityType());
1457    G4String matName = material? material->GetName(): G4String("No material");
1458    setAttribute(instance, "Material",   matName );
1459    G4double matDensity = material? material->GetDensity(): 0.;
1460    setAttribute(instance, "Density",    matDensity);
1461    G4double matRadlen = material? material->GetRadlen(): 0.;
1462    setAttribute(instance, "Radlen",     matRadlen);
1463   
1464    G4State matState = material? material->GetState(): kStateUndefined;
1465    G4String state = materialState[matState];
1466    setAttribute(instance, "State", state);
1467   
1468    return instance;
1469}
1470
1471HepRepInstance* G4HepRepSceneHandler::getGeometryInstance(G4String volumeName, int depth) {
1472    // no extra checks since these are done in the geometryType already
1473
1474    // adjust depth, also pop the current instance
1475    while ((int)_geometryInstance.size() > depth) {
1476        _geometryInstance.pop_back();
1477    }
1478   
1479    // get parent
1480    HepRepInstance* parent = (_geometryInstance.empty()) ? getGeometryRootInstance() : _geometryInstance.back();
1481   
1482    // get type
1483    HepRepType* type = getGeometryType(volumeName, depth);
1484   
1485    // create instance
1486    HepRepInstance* instance = factory->createHepRepInstance(parent, type);
1487    _geometryInstance.push_back(instance);
1488
1489    return instance;
1490}
1491
1492HepRepTypeTree* G4HepRepSceneHandler::getGeometryTypeTree() {
1493    if (_geometryTypeTree == NULL) {
1494        // Create the Geometry TypeTree.
1495        HepRepTreeID* geometryTreeID = factory->createHepRepTreeID("G4GeometryTypes", "1.0");
1496        _geometryTypeTree = factory->createHepRepTypeTree(geometryTreeID);
[1228]1497               
1498                G4HepRepMessenger* messenger = G4HepRepMessenger::GetInstance();
1499        if ( messenger->appendGeometry()) {
[834]1500            getHepRep()->addTypeTree(_geometryTypeTree);
1501        } else {
1502            getHepRepGeometry()->addTypeTree(_geometryTypeTree);
1503        }
1504    }
1505    return _geometryTypeTree;
1506}
1507
1508HepRepType* G4HepRepSceneHandler::getGeometryRootType() {
1509    if (_geometryRootType == NULL) {
1510        // Create the top level Geometry Type.
1511        _geometryRootType = factory->createHepRepType(getGeometryTypeTree(), rootVolumeName);
1512        _geometryRootType->addAttValue("Layer", geometryLayer);
1513   
1514        // Add attdefs used by all geometry types.
1515        _geometryRootType->addAttDef  ("LVol", "Logical Volume", "Physics","");
1516        _geometryRootType->addAttValue("LVol", G4String(""));
1517        _geometryRootType->addAttDef  ("Region", "Cuts Region", "Physics","");
1518        _geometryRootType->addAttValue("Region", G4String(""));
1519        _geometryRootType->addAttDef  ("RootRegion", "Root Region", "Physics","");
1520        _geometryRootType->addAttValue("RootRegion", false);
1521        _geometryRootType->addAttDef  ("Solid", "Solid Name", "Physics","");
1522        _geometryRootType->addAttValue("Solid", G4String(""));
1523        _geometryRootType->addAttDef  ("EType", "Entity Type", "Physics","");
1524        _geometryRootType->addAttValue("EType", G4String("G4Box"));
1525        _geometryRootType->addAttDef  ("Material", "Material Name", "Physics","");
1526        _geometryRootType->addAttValue("Material", G4String("Air"));
1527        _geometryRootType->addAttDef  ("Density", "Material Density", "Physics","");
1528        _geometryRootType->addAttValue("Density", 0.0);
1529        _geometryRootType->addAttDef  ("State", "Material State", "Physics","");
1530        _geometryRootType->addAttValue("State", G4String("Gas"));
1531        _geometryRootType->addAttDef  ("Radlen", "Material Radiation Length", "Physics","");
1532        _geometryRootType->addAttValue("Radlen", 0.0);
1533       
1534        // add defaults for Geometry
1535        _geometryRootType->addAttValue("Color", 0.8, 0.8, 0.8, 1.0);
1536        _geometryRootType->addAttValue("Visibility", true);
1537        _geometryRootType->addAttValue("FillColor", 0.8, 0.8, 0.8, 1.0);
1538        _geometryRootType->addAttValue("LineWidth", 1.0);
1539        _geometryRootType->addAttValue("DrawAs", G4String("Polygon"));
1540        _geometryRootType->addAttValue("PickParent", false);
1541        _geometryRootType->addAttValue("ShowParentAttributes", true);
1542       
1543        _geometryRootType->addAttValue("MarkSizeMultiplier", 4.0);
1544        _geometryRootType->addAttValue("LineWidthMultiplier", 1.0);
1545
1546        addTopLevelAttributes(_geometryRootType);       
1547               
1548        _geometryType["/"+_geometryRootType->getName()] = _geometryRootType;
1549    }
1550    return _geometryRootType;
1551}
1552
1553HepRepType* G4HepRepSceneHandler::getGeometryType(G4String volumeName, int depth) {
1554    // make sure we have a root
1555    getGeometryRootType();   
1556
1557    // construct the full name for this volume
1558    G4String name = getFullTypeName(volumeName, depth);
1559   
1560    // lookup type and create if necessary
1561    HepRepType* type = _geometryType[name];   
1562    if (type == NULL) {
1563        G4String parentName = getParentTypeName(depth);
1564        HepRepType* parentType = _geometryType[parentName];
1565        // HepRep uses hierarchical names
1566        type = factory->createHepRepType(parentType, volumeName);
1567        _geometryType[name] = type;
1568    }
1569    return type;   
1570}
1571
1572G4String G4HepRepSceneHandler::getFullTypeName(G4String volumeName, int depth) {
1573    // check for name depth
1574    if (depth > (int)_geometryTypeName.size()) {
1575        // there is a problem, book this type under problems
1576        G4String problem = "HierarchyProblem";
1577        if (_geometryType["/"+problem] == NULL) {
1578            // HepRep uses hierarchical names
1579            HepRepType* type = factory->createHepRepType(getGeometryRootType(), problem);
1580            _geometryType["/"+problem] = type;
1581        }
1582        return "/" + problem + "/" + volumeName;
1583    }
1584   
1585    // adjust name depth, also pop the current volumeName
1586    while ((int)_geometryTypeName.size() > depth) {
1587        _geometryTypeName.pop_back();
1588    }
1589   
1590    // construct full name and push it
1591    G4String name = (_geometryTypeName.empty()) ? G4String("/"+rootVolumeName) : _geometryTypeName.back();
1592    name = name + "/" + volumeName;
1593    _geometryTypeName.push_back(name);
1594    return name;
1595}
1596
1597G4String G4HepRepSceneHandler::getParentTypeName(int depth) {
1598    return (depth >= 1) ? _geometryTypeName[depth-1] : G4String("/"+rootVolumeName);
1599}       
1600
1601HepRepInstanceTree* G4HepRepSceneHandler::getEventInstanceTree() {
1602    if (_eventInstanceTree == NULL) {
1603        // Create the Event InstanceTree.
1604        _eventInstanceTree = factory->createHepRepInstanceTree("G4EventData", "1.0", getEventTypeTree());
1605        getHepRep()->addInstanceTree(_eventInstanceTree);
1606    }
1607    return _eventInstanceTree;
1608}
1609
1610HepRepInstance* G4HepRepSceneHandler::getEventInstance() {
1611    if (_eventInstance == NULL) {
1612        // Create the top level Event Instance.
1613        _eventInstance = factory->createHepRepInstance(getEventInstanceTree(), getEventType());
1614    }
1615    return _eventInstance;
1616}
1617
1618HepRepTypeTree* G4HepRepSceneHandler::getEventTypeTree() {
1619    if (_eventTypeTree == NULL) {
1620        // Create the Event TypeTree.
1621        HepRepTreeID* eventTreeID = factory->createHepRepTreeID("G4EventTypes", "1.0");
1622        _eventTypeTree = factory->createHepRepTypeTree(eventTreeID);
1623        getHepRep()->addTypeTree(_eventTypeTree);   
1624    }
1625   
1626    return _eventTypeTree;
1627}
1628
1629HepRepType* G4HepRepSceneHandler::getEventType() {
1630    if (_eventType == NULL) {
1631        // Create the top level Event Type.
1632        _eventType = factory->createHepRepType(getEventTypeTree(), "Event");
1633        _eventType->addAttValue("Layer", eventLayer);
1634
1635        // add defaults for Events
1636        _eventType->addAttValue("Visibility", true);
1637        _eventType->addAttValue("Color", 1.0, 1.0, 1.0, 1.0);
1638        _eventType->addAttValue("FillColor", 1.0, 1.0, 1.0, 1.0);
1639        _eventType->addAttValue("LineWidth", 1.0);
1640        _eventType->addAttValue("HasFrame", true);
1641        _eventType->addAttValue("PickParent", false);
1642        _eventType->addAttValue("ShowParentAttributes", false);
1643
1644        _eventType->addAttValue("MarkSizeMultiplier", 4.0);
1645        _eventType->addAttValue("LineWidthMultiplier", 1.0);
1646
1647        addTopLevelAttributes(_eventType);
1648    }
1649   
1650    return _eventType;
1651}
1652
1653HepRepType* G4HepRepSceneHandler::getTrajectoryType() {
1654    if (_trajectoryType == NULL) {
1655        _trajectoryType = factory->createHepRepType(getEventType(), "Trajectory");
1656       
1657        _trajectoryType->addAttValue("Layer", trajectoryLayer);
1658        _trajectoryType->addAttValue("DrawAs", G4String("Line"));
1659
1660        _trajectoryType->addAttValue("LineWidthMultiplier", 2.0);
1661       
1662        // attributes to draw the points of a track as markers.
1663        _trajectoryType->addAttValue("MarkName", G4String("Box"));
1664        _trajectoryType->addAttValue("MarkSize", 4);
1665        _trajectoryType->addAttValue("MarkType", G4String("Symbol"));
1666        _trajectoryType->addAttValue("Fill", true);
1667    }
1668    return _trajectoryType;
1669}
1670
1671HepRepType* G4HepRepSceneHandler::getHitType() {
1672    if (_hitType == NULL) {
1673        _hitType = factory->createHepRepType(getEventType(), "Hit");
1674        _hitType->addAttValue("Layer", hitLayer);
1675        _hitType->addAttValue("DrawAs", G4String("Point"));
1676        _hitType->addAttValue("MarkName", G4String("Box"));
1677        _hitType->addAttValue("MarkSize", 4.0);
1678        _hitType->addAttValue("MarkType", G4String("Symbol"));
1679        _hitType->addAttValue("Fill", true);
1680    }
1681    return _hitType;
1682}
1683
1684HepRepType* G4HepRepSceneHandler::getCalHitType() {
1685    if (_calHitType == NULL) {
1686        _calHitType = factory->createHepRepType(getEventType(), "CalHit");
1687        _calHitType->addAttValue("Layer", calHitLayer);
1688        _calHitType->addAttValue("Fill", true);
1689        _calHitType->addAttValue("DrawAs", G4String("Polygon"));
1690    }
1691    return _calHitType;
1692}
1693
1694HepRepType* G4HepRepSceneHandler::getCalHitFaceType() {
1695    if (_calHitFaceType == NULL) {
1696        _calHitFaceType = factory->createHepRepType(getCalHitType(), "CalHitFace");
1697        _calHitFaceType->addAttValue("PickParent", true);
1698    }
1699    return _calHitFaceType;
1700}
1701
Note: See TracBrowser for help on using the repository browser.