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

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

geant4 tag 9.4

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