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

Last change on this file since 1277 was 1228, checked in by garnier, 14 years ago

update geant4.9.3 tag

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