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