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

Last change on this file since 1215 was 944, checked in by garnier, 17 years ago

mise a jour des tags

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