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