source: trunk/source/visualization/HepRep/src/G4HepRepFileSceneHandler.cc @ 834

Last change on this file since 834 was 834, checked in by garnier, 16 years ago

import all except CVS

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