source: trunk/geant4/visualization/management/src/G4VSceneHandler.cc@ 597

Last change on this file since 597 was 593, checked in by garnier, 18 years ago

r627@mac-90108: laurentgarnier | 2007-11-09 07:57:42 +0100
modif dans les includes directives

  • Property svn:mime-type set to text/cpp
File size: 27.9 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//
27// $Id: G4VSceneHandler.cc,v 1.82 2007/05/16 15:47:44 allison Exp $
28// GEANT4 tag $Name: geant4-09-00-ref-01 $
29//
30//
31// John Allison 19th July 1996
32// Abstract interface class for graphics scenes.
33
34#include "G4VSceneHandler.hh"
35
36#include "G4ios.hh"
37#include <sstream>
38
39#include "G4VisManager.hh"
40#include "G4VGraphicsSystem.hh"
41#include "G4VViewer.hh"
42#include "G4VSolid.hh"
43#include "G4RotationMatrix.hh"
44#include "G4ThreeVector.hh"
45#include "G4VPhysicalVolume.hh"
46#include "G4Material.hh"
47#include "G4Polyline.hh"
48#include "G4Scale.hh"
49#include "G4Text.hh"
50#include "G4Circle.hh"
51#include "G4Square.hh"
52#include "G4Polymarker.hh"
53#include "G4Polyhedron.hh"
54#include "G4NURBS.hh"
55#include "G4Visible.hh"
56#include "G4VisAttributes.hh"
57#include "G4VModel.hh"
58#include "G4TrajectoriesModel.hh"
59#include "G4Box.hh"
60#include "G4Cons.hh"
61#include "G4Tubs.hh"
62#include "G4Trd.hh"
63#include "G4Trap.hh"
64#include "G4Sphere.hh"
65#include "G4Para.hh"
66#include "G4Torus.hh"
67#include "G4Polycone.hh"
68#include "G4Polyhedra.hh"
69#include "G4LogicalVolume.hh"
70#include "G4PhysicalVolumeModel.hh"
71#include "G4ModelingParameters.hh"
72#include "G4VTrajectory.hh"
73#include "G4VTrajectoryPoint.hh"
74#include "G4HitsModel.hh"
75#include "G4VHit.hh"
76#include "Randomize.hh"
77#include "G4StateManager.hh"
78#include "G4RunManager.hh"
79#include "G4Run.hh"
80#include "G4Transform3D.hh"
81#include "G4AttHolder.hh"
82#include "G4AttDef.hh"
83
84G4VSceneHandler::G4VSceneHandler (G4VGraphicsSystem& system, G4int id, const G4String& name):
85 fSystem (system),
86 fSceneHandlerId (id),
87 fViewCount (0),
88 fpViewer (0),
89 fpScene (0),
90 fMarkForClearingTransientStore (true), // Ready for first
91 // ClearTransientStoreIfMarked(),
92 // e.g., at end of run (see
93 // G4VisManager.cc).
94 fReadyForTransients (true), // Only false while processing scene.
95 fProcessingSolid (false),
96 fSecondPassRequested (false),
97 fSecondPass (false),
98 fpModel (0),
99 fpObjectTransformation (0),
100 fNestingDepth (0),
101 fpVisAttribs (0)
102{
103 G4VisManager* pVMan = G4VisManager::GetInstance ();
104 fpScene = pVMan -> GetCurrentScene ();
105 if (name == "") {
106 std::ostringstream ost;
107 ost << fSystem.GetName () << '-' << fSceneHandlerId;
108 fName = ost.str();
109 }
110 else {
111 fName = name;
112 }
113 fTransientsDrawnThisEvent = pVMan->GetTransientsDrawnThisEvent();
114 fTransientsDrawnThisRun = pVMan->GetTransientsDrawnThisRun();
115}
116
117G4VSceneHandler::~G4VSceneHandler () {
118 G4ViewerListIterator i;
119 for (i = fViewerList.begin(); i != fViewerList.end(); ++i) {
120 delete *i;
121 }
122}
123
124void G4VSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
125 const G4VisAttributes& visAttribs) {
126 fpObjectTransformation = &objectTransformation;
127 fpVisAttribs = &visAttribs;
128 fProcessingSolid = true;
129}
130
131void G4VSceneHandler::PostAddSolid () {
132 fpObjectTransformation = 0;
133 fpVisAttribs = 0;
134 fProcessingSolid = false;
135 if (fReadyForTransients) {
136 fTransientsDrawnThisEvent = true;
137 fTransientsDrawnThisRun = true;
138 }
139}
140
141void G4VSceneHandler::BeginPrimitives
142(const G4Transform3D& objectTransformation) {
143 fNestingDepth++;
144 if (fNestingDepth > 1)
145 G4Exception("G4VSceneHandler::BeginPrimitives: Nesting detected."
146 "\n It is illegal to nest Begin/EndPrimitives.");
147 fpObjectTransformation = &objectTransformation;
148}
149
150void G4VSceneHandler::EndPrimitives () {
151 if (fNestingDepth <= 0)
152 G4Exception("G4VSceneHandler::EndPrimitives: Nesting error");
153 fNestingDepth--;
154 fpObjectTransformation = 0;
155 if (fReadyForTransients) {
156 fTransientsDrawnThisEvent = true;
157 fTransientsDrawnThisRun = true;
158 }
159}
160
161void G4VSceneHandler::BeginPrimitives2D () {
162 fNestingDepth++;
163 if (fNestingDepth > 1)
164 G4Exception("G4VSceneHandler::BeginPrimitives2D: Nesting detected."
165 "\n It is illegal to nest Begin/EndPrimitives.");
166 // Not actually required for 2D operations but some drivers do an
167 // initial transformation...
168 fpObjectTransformation = &fIdentityTransformation;
169}
170
171void G4VSceneHandler::EndPrimitives2D () {
172 if (fNestingDepth <= 0)
173 G4Exception("G4VSceneHandler::EndPrimitives2D: Nesting error");
174 fNestingDepth--;
175 fpObjectTransformation = 0;
176 if (fReadyForTransients) {
177 fTransientsDrawnThisEvent = true;
178 fTransientsDrawnThisRun = true;
179 }
180}
181
182void G4VSceneHandler::BeginModeling () {
183}
184
185void G4VSceneHandler::EndModeling ()
186{
187 fpModel = 0;
188}
189
190void G4VSceneHandler::ClearStore () {
191 // if (fpViewer) fpViewer -> NeedKernelVisit (true);
192 // ?? Viewer is supposed to be smart enough to know when to visit
193 // kernel, but a problem in OpenGL Stored seems to require a forced
194 // kernel visit triggered by the above code. John Allison Aug 2001
195 // Feb 2005 - commented out. Let's fix OpenGL if necessary.
196}
197
198void G4VSceneHandler::ClearTransientStore () {
199}
200
201void G4VSceneHandler::AddSolid (const G4Box& box) {
202 RequestPrimitives (box);
203// If your graphics system is sophisticated enough to handle a
204// particular solid shape as a primitive, in your derived class write a
205// function to override this. (Note: some compilers warn that your
206// function "hides" this one. That's OK.)
207// Your function might look like this...
208// void G4MyScene::AddSolid (const G4Box& box) {
209// Get parameters of appropriate object, e.g.:
210// G4double dx = box.GetXHalfLength ();
211// G4double dy = box.GetYHalfLength ();
212// G4double dz = box.GetZHalfLength ();
213// and Draw or Store in your display List.
214}
215
216void G4VSceneHandler::AddSolid (const G4Tubs& tubs) {
217 RequestPrimitives (tubs);
218}
219
220void G4VSceneHandler::AddSolid (const G4Cons& cons) {
221 RequestPrimitives (cons);
222}
223
224void G4VSceneHandler::AddSolid (const G4Trd& trd) {
225 RequestPrimitives (trd);
226}
227
228void G4VSceneHandler::AddSolid (const G4Trap& trap) {
229 RequestPrimitives (trap);
230}
231
232void G4VSceneHandler::AddSolid (const G4Sphere& sphere) {
233 RequestPrimitives (sphere );
234}
235
236void G4VSceneHandler::AddSolid (const G4Para& para) {
237 RequestPrimitives (para);
238}
239
240void G4VSceneHandler::AddSolid (const G4Torus& torus) {
241 RequestPrimitives (torus);
242}
243
244void G4VSceneHandler::AddSolid (const G4Polycone& polycone) {
245 RequestPrimitives (polycone);
246}
247
248void G4VSceneHandler::AddSolid (const G4Polyhedra& polyhedra) {
249 RequestPrimitives (polyhedra);
250}
251
252void G4VSceneHandler::AddSolid (const G4VSolid& solid) {
253 RequestPrimitives (solid);
254}
255
256void G4VSceneHandler::AddCompound (const G4VTrajectory& traj) {
257 G4TrajectoriesModel* pTrModel =
258 dynamic_cast<G4TrajectoriesModel*>(fpModel);
259 if (!pTrModel) G4Exception
260 ("G4VSceneHandler::AddCompound(const G4VTrajectory&): Not a G4TrajectoriesModel.");
261 traj.DrawTrajectory(pTrModel->GetDrawingMode());
262}
263
264void G4VSceneHandler::AddCompound (const G4VHit& hit) {
265 ((G4VHit&)hit).Draw(); // Cast to non-const because Draw is non-const!!!!
266}
267
268void G4VSceneHandler::AddViewerToList (G4VViewer* pViewer) {
269 fViewerList.push_back (pViewer);
270}
271
272void G4VSceneHandler::AddPrimitive (const G4Scale& scale) {
273
274 const G4double margin(0.01);
275 // Fractional margin - ensures scale is comfortably inside viewing
276 // volume.
277 const G4double oneMinusMargin (1. - margin);
278
279 const G4VisExtent& sceneExtent = fpScene->GetExtent();
280
281 // Useful constants...
282 const G4double length(scale.GetLength());
283 const G4double halfLength(length / 2.);
284 const G4double tickLength(length / 20.);
285 const G4double piBy2(halfpi);
286
287 // Get size of scene...
288 const G4double xmin = sceneExtent.GetXmin();
289 const G4double xmax = sceneExtent.GetXmax();
290 const G4double ymin = sceneExtent.GetYmin();
291 const G4double ymax = sceneExtent.GetYmax();
292 const G4double zmin = sceneExtent.GetZmin();
293 const G4double zmax = sceneExtent.GetZmax();
294
295 // Create (empty) polylines having the same vis attributes...
296 G4Polyline scaleLine, tick11, tick12, tick21, tick22;
297 G4VisAttributes visAtts(*scale.GetVisAttributes()); // Long enough life.
298 scaleLine.SetVisAttributes(&visAtts);
299 tick11.SetVisAttributes(&visAtts);
300 tick12.SetVisAttributes(&visAtts);
301 tick21.SetVisAttributes(&visAtts);
302 tick22.SetVisAttributes(&visAtts);
303
304 // Add points to the polylines to represent an scale parallel to the
305 // x-axis centred on the origin...
306 G4Point3D r1(G4Point3D(-halfLength, 0., 0.));
307 G4Point3D r2(G4Point3D( halfLength, 0., 0.));
308 scaleLine.push_back(r1);
309 scaleLine.push_back(r2);
310 G4Point3D ticky(0., tickLength, 0.);
311 G4Point3D tickz(0., 0., tickLength);
312 tick11.push_back(r1 + ticky);
313 tick11.push_back(r1 - ticky);
314 tick12.push_back(r1 + tickz);
315 tick12.push_back(r1 - tickz);
316 tick21.push_back(r2 + ticky);
317 tick21.push_back(r2 - ticky);
318 tick22.push_back(r2 + tickz);
319 tick22.push_back(r2 - tickz);
320 G4Point3D textPosition(0., tickLength, 0.);
321
322 // Transform appropriately...
323
324 G4Transform3D transformation;
325 if (scale.GetAutoPlacing()) {
326 G4Transform3D rotation;
327 switch (scale.GetDirection()) {
328 case G4Scale::x:
329 break;
330 case G4Scale::y:
331 rotation = G4RotateZ3D(piBy2);
332 break;
333 case G4Scale::z:
334 rotation = G4RotateY3D(piBy2);
335 break;
336 }
337 G4double sxmid(scale.GetXmid());
338 G4double symid(scale.GetYmid());
339 G4double szmid(scale.GetZmid());
340 sxmid = xmin + oneMinusMargin * (xmax - xmin);
341 symid = ymin + margin * (ymax - ymin);
342 szmid = zmin + oneMinusMargin * (zmax - zmin);
343 switch (scale.GetDirection()) {
344 case G4Scale::x:
345 sxmid -= halfLength;
346 break;
347 case G4Scale::y:
348 symid += halfLength;
349 break;
350 case G4Scale::z:
351 szmid -= halfLength;
352 break;
353 }
354 G4Translate3D translation(sxmid, symid, szmid);
355 transformation = translation * rotation;
356 } else {
357 if (fpModel) transformation = fpModel->GetTransformation();
358 }
359
360 // Draw...
361 // We would like to call BeginPrimitives(transformation) here but
362 // calling BeginPrimitives from within an AddPrimitive is not
363 // allowed! So we have to do our own transformation...
364 AddPrimitive(scaleLine.transform(transformation));
365 AddPrimitive(tick11.transform(transformation));
366 AddPrimitive(tick12.transform(transformation));
367 AddPrimitive(tick21.transform(transformation));
368 AddPrimitive(tick22.transform(transformation));
369 G4Text text(scale.GetAnnotation(),textPosition.transform(transformation));
370 text.SetScreenSize(12.);
371 AddPrimitive(text);
372}
373
374void G4VSceneHandler::AddPrimitive (const G4Polymarker& polymarker) {
375 switch (polymarker.GetMarkerType()) {
376 default:
377 case G4Polymarker::dots:
378 {
379 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
380 G4Circle dot (polymarker);
381 dot.SetPosition (polymarker[iPoint]);
382 dot.SetWorldSize (0.);
383 dot.SetScreenSize (0.1); // Very small circle.
384 AddPrimitive (dot);
385 }
386 }
387 break;
388 case G4Polymarker::circles:
389 {
390 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
391 G4Circle circle (polymarker);
392 circle.SetPosition (polymarker[iPoint]);
393 AddPrimitive (circle);
394 }
395 }
396 break;
397 case G4Polymarker::squares:
398 {
399 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
400 G4Square square (polymarker);
401 square.SetPosition (polymarker[iPoint]);
402 AddPrimitive (square);
403 }
404 }
405 break;
406 }
407}
408
409void G4VSceneHandler::RemoveViewerFromList (G4VViewer* pViewer) {
410 fViewerList.remove(pViewer);
411}
412
413void G4VSceneHandler::SetScene (G4Scene* pScene) {
414 fpScene = pScene;
415 // Notify all viewers that a kernel visit is required.
416 G4ViewerListIterator i;
417 for (i = fViewerList.begin(); i != fViewerList.end(); i++) {
418 (*i) -> SetNeedKernelVisit (true);
419 }
420}
421
422void G4VSceneHandler::RequestPrimitives (const G4VSolid& solid) {
423 BeginPrimitives (*fpObjectTransformation);
424 G4NURBS* pNURBS = 0;
425 G4Polyhedron* pPolyhedron = 0;
426 switch (fpViewer -> GetViewParameters () . GetRepStyle ()) {
427 case G4ViewParameters::nurbs:
428 pNURBS = solid.CreateNURBS ();
429 if (pNURBS) {
430 pNURBS -> SetVisAttributes (fpVisAttribs);
431 AddPrimitive (*pNURBS);
432 delete pNURBS;
433 break;
434 }
435 else {
436 G4VisManager::Verbosity verbosity =
437 G4VisManager::GetInstance()->GetVerbosity();
438 if (verbosity >= G4VisManager::errors) {
439 G4cout <<
440 "ERROR: G4VSceneHandler::RequestPrimitives"
441 "\n NURBS not available for "
442 << solid.GetName () << G4endl;
443 G4cout << "Trying polyhedron." << G4endl;
444 }
445 }
446 // Dropping through to polyhedron...
447 case G4ViewParameters::polyhedron:
448 default:
449 G4Polyhedron::SetNumberOfRotationSteps (GetNoOfSides (fpVisAttribs));
450 pPolyhedron = solid.GetPolyhedron ();
451 G4Polyhedron::ResetNumberOfRotationSteps ();
452 if (pPolyhedron) {
453 pPolyhedron -> SetVisAttributes (fpVisAttribs);
454 AddPrimitive (*pPolyhedron);
455 }
456 else {
457 G4VisManager::Verbosity verbosity =
458 G4VisManager::GetInstance()->GetVerbosity();
459 if (verbosity >= G4VisManager::errors) {
460 G4cout <<
461 "ERROR: G4VSceneHandler::RequestPrimitives"
462 "\n Polyhedron not available for " << solid.GetName () <<
463 ".\n This means it cannot be visualized on most systems."
464 "\n Contact the Visualization Coordinator." << G4endl;
465 }
466 }
467 break;
468 }
469 EndPrimitives ();
470}
471
472void G4VSceneHandler::ProcessScene (G4VViewer&) {
473
474 if (!fpScene) return;
475
476 G4VisManager* visManager = G4VisManager::GetInstance();
477
478 if (!visManager->GetConcreteInstance()) return;
479
480 G4VisManager::Verbosity verbosity = visManager->GetVerbosity();
481
482 fReadyForTransients = false;
483
484 // Clear stored scene, if any, i.e., display lists, scene graphs.
485 ClearStore ();
486
487 // Reset fMarkForClearingTransientStore. No need to clear transient
488 // store since it has just been cleared above. (Leaving
489 // fMarkForClearingTransientStore true causes problems with
490 // recomputing transients below.) Restore it again at end...
491 G4bool tmpMarkForClearingTransientStore = fMarkForClearingTransientStore;
492 fMarkForClearingTransientStore = false;
493
494 // Traverse geometry tree and send drawing primitives to window(s).
495
496 const std::vector<G4VModel*>& runDurationModelList =
497 fpScene -> GetRunDurationModelList ();
498
499 if (runDurationModelList.size ()) {
500 if (verbosity >= G4VisManager::confirmations) {
501 G4cout << "Traversing scene data..." << G4endl;
502 }
503
504 BeginModeling ();
505
506 // Create modeling parameters from view parameters...
507 G4ModelingParameters* pMP = CreateModelingParameters ();
508
509 for (size_t i = 0; i < runDurationModelList.size (); i++) {
510 G4VModel* pModel = runDurationModelList[i];
511 // Note: this is not the place to take action on
512 // pModel->GetTransformation(). The model must take care of
513 // this in pModel->DescribeYourselfTo(*this). See, for example,
514 // G4PhysicalVolumeModel and /vis/scene/add/logo.
515 pModel -> SetModelingParameters (pMP);
516 SetModel (pModel); // Store for use by derived class.
517 pModel -> DescribeYourselfTo (*this);
518 pModel -> SetModelingParameters (0);
519 }
520
521 // Repeat if required...
522 if (fSecondPassRequested) {
523 fSecondPass = true;
524 for (size_t i = 0; i < runDurationModelList.size (); i++) {
525 G4VModel* pModel = runDurationModelList[i];
526 pModel -> SetModelingParameters (pMP);
527 SetModel (pModel); // Store for use by derived class.
528 pModel -> DescribeYourselfTo (*this);
529 pModel -> SetModelingParameters (0);
530 }
531 fSecondPass = false;
532 fSecondPassRequested = false;
533 }
534
535 delete pMP;
536 EndModeling ();
537
538 }
539
540 fpViewer->FinishView(); // Flush streams and/or swap buffers.
541
542 fReadyForTransients = true;
543
544 // Refresh event from end-of-event model list. Allow only in Idle state...
545 G4StateManager* stateManager = G4StateManager::GetStateManager();
546 G4ApplicationState state = stateManager->GetCurrentState();
547 if (state == G4State_Idle) {
548
549 visManager->SetEventRefreshing(true);
550
551 if (visManager->GetRequestedEvent()) {
552 DrawEvent(visManager->GetRequestedEvent());
553
554 } else {
555
556 G4RunManager* runManager = G4RunManager::GetRunManager();
557 if (runManager) {
558 const G4Run* run = runManager->GetCurrentRun();
559 const std::vector<const G4Event*>* events =
560 run? run->GetEventVector(): 0;
561 size_t nKeptEvents = 0;
562 if (events) nKeptEvents = events->size();
563 if (nKeptEvents) {
564
565 if (fpScene->GetRefreshAtEndOfEvent()) {
566
567 if (verbosity >= G4VisManager::confirmations) {
568 G4cout << "Refreshing event..." << G4endl;
569 }
570 const G4Event* event = 0;
571 if (events && events->size()) event = events->back();
572 if (event) DrawEvent(event);
573
574 } else { // Accumulating events.
575
576 if (verbosity >= G4VisManager::confirmations) {
577 G4cout << "Refreshing events in run..." << G4endl;
578 }
579 for (size_t i = 0; i < nKeptEvents; ++i) {
580 const G4Event* event = (*events)[i];
581 if (event) DrawEvent(event);
582 }
583
584 if (!fpScene->GetRefreshAtEndOfRun()) {
585 if (verbosity >= G4VisManager::warnings) {
586 G4cout <<
587 "WARNING: Cannot refresh events accumulated over more"
588 "\n than one runs. Refreshed just the last run."
589 << G4endl;
590 }
591 }
592 }
593 }
594 }
595 }
596 visManager->SetEventRefreshing(false);
597 }
598
599 fMarkForClearingTransientStore = tmpMarkForClearingTransientStore;
600}
601
602void G4VSceneHandler::DrawEvent(const G4Event* event)
603{
604 const std::vector<G4VModel*>& EOEModelList =
605 fpScene -> GetEndOfEventModelList ();
606 size_t nModels = EOEModelList.size();
607 if (nModels) {
608 G4ModelingParameters* pMP = CreateModelingParameters();
609 pMP->SetEvent(event);
610 for (size_t i = 0; i < nModels; i++) {
611 G4VModel* pModel = EOEModelList [i];
612 pModel -> SetModelingParameters(pMP);
613 SetModel (pModel);
614 pModel -> DescribeYourselfTo (*this);
615 pModel -> SetModelingParameters(0);
616 }
617 delete pMP;
618 SetModel (0);
619 }
620}
621
622G4ModelingParameters* G4VSceneHandler::CreateModelingParameters ()
623{
624 // Create modeling parameters from View Parameters...
625 const G4ViewParameters& vp = fpViewer -> GetViewParameters ();
626
627 // Convert drawing styles...
628 G4ModelingParameters::DrawingStyle modelDrawingStyle =
629 G4ModelingParameters::wf;
630 switch (vp.GetDrawingStyle ()) {
631 default:
632 case G4ViewParameters::wireframe:
633 modelDrawingStyle = G4ModelingParameters::wf;
634 break;
635 case G4ViewParameters::hlr:
636 modelDrawingStyle = G4ModelingParameters::hlr;
637 break;
638 case G4ViewParameters::hsr:
639 modelDrawingStyle = G4ModelingParameters::hsr;
640 break;
641 case G4ViewParameters::hlhsr:
642 modelDrawingStyle = G4ModelingParameters::hlhsr;
643 break;
644 }
645
646 // Decide if covered daughters are really to be culled...
647 G4bool reallyCullCovered =
648 vp.IsCullingCovered() // Culling daughters depends also on...
649 && !vp.IsSection () // Sections (DCUT) not requested.
650 && !vp.IsCutaway () // Cutaways not requested.
651 ;
652
653 G4ModelingParameters* pModelingParams = new G4ModelingParameters
654 (vp.GetDefaultVisAttributes (),
655 modelDrawingStyle,
656 vp.IsCulling (),
657 vp.IsCullingInvisible (),
658 vp.IsDensityCulling (),
659 vp.GetVisibleDensity (),
660 reallyCullCovered,
661 vp.GetNoOfSides ()
662 );
663
664 pModelingParams->SetWarning
665 (G4VisManager::GetInstance()->GetVerbosity() >= G4VisManager::warnings);
666
667 pModelingParams->SetExplodeFactor(vp.GetExplodeFactor());
668 pModelingParams->SetExplodeCentre(vp.GetExplodeCentre());
669
670 pModelingParams->SetSectionPolyhedron(CreateSectionPolyhedron());
671 pModelingParams->SetCutawayPolyhedron(CreateCutawayPolyhedron());
672 // The polyhedron objects are deleted in the modeling parameters destructor.
673
674 return pModelingParams;
675}
676
677const G4Polyhedron* G4VSceneHandler::CreateSectionPolyhedron()
678{
679 /* Disable for now. Boolean processor not up to it.
680 const G4ViewParameters& vp = fpViewer->GetViewParameters();
681 if (vp.IsSection () ) {
682 G4double radius = fpScene->GetExtent().GetExtentRadius();
683 G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
684 G4Box sectionBox("clipper",
685 safe, safe, 1.e-5 * radius); // Thin in z-plane.
686 G4Polyhedron* sectioner = sectionBox.CreatePolyhedron();
687 const G4Plane3D& s = vp.GetSectionPlane ();
688 G4double a = s.a();
689 G4double b = s.b();
690 G4double c = s.c();
691 G4double d = s.d();
692 G4Transform3D transform = G4TranslateZ3D(-d);
693 const G4Normal3D normal(a,b,c);
694 if (normal != G4Normal3D(0,0,1)) {
695 const G4double angle = std::acos(normal.dot(G4Normal3D(0,0,1)));
696 const G4Vector3D axis = G4Normal3D(0,0,1).cross(normal);
697 transform = G4Rotate3D(angle, axis) * transform;
698 }
699 sectioner->Transform(transform);
700 return sectioner;
701 } else {
702 return 0;
703 }
704 */
705 return 0;
706}
707
708const G4Polyhedron* G4VSceneHandler::CreateCutawayPolyhedron()
709{
710 return 0;
711}
712
713void G4VSceneHandler::LoadAtts(const G4Visible& visible, G4AttHolder* holder)
714{
715 // Load G4Atts from G4VisAttributes, if any...
716 const G4VisAttributes* va = visible.GetVisAttributes();
717 if (va) {
718 const std::map<G4String,G4AttDef>* vaDefs =
719 va->GetAttDefs();
720 if (vaDefs) {
721 holder->AddAtts(visible.GetVisAttributes()->CreateAttValues(), vaDefs);
722 }
723 }
724
725 G4PhysicalVolumeModel* pPVModel =
726 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
727 if (pPVModel) {
728 // Load G4Atts from G4PhysicalVolumeModel...
729 const std::map<G4String,G4AttDef>* defs = pPVModel->GetAttDefs();
730 if (defs) {
731 holder->AddAtts(pPVModel->CreateCurrentAttValues(), defs);
732 }
733 }
734
735 G4TrajectoriesModel* trajModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
736 if (trajModel) {
737 // Load G4Atts from trajectory...
738 const G4VTrajectory* traj = trajModel->GetCurrentTrajectory();
739 const std::map<G4String,G4AttDef>* defs = traj->GetAttDefs();
740 if (defs) {
741 holder->AddAtts(traj->CreateAttValues(), defs);
742 }
743 G4int nPoints = traj->GetPointEntries();
744 for (G4int i = 0; i < nPoints; ++i) {
745 G4VTrajectoryPoint* trajPoint = traj->GetPoint(i);
746 const std::map<G4String,G4AttDef>* defs = trajPoint->GetAttDefs();
747 if (defs) {
748 holder->AddAtts(trajPoint->CreateAttValues(), defs);
749 }
750 }
751 }
752
753 G4HitsModel* hitsModel = dynamic_cast<G4HitsModel*>(fpModel);
754 if (hitsModel) {
755 // Load G4Atts from hit...
756 const G4VHit* hit = hitsModel->GetCurrentHit();
757 const std::map<G4String,G4AttDef>* defs = hit->GetAttDefs();
758 if (defs) {
759 holder->AddAtts(hit->CreateAttValues(), defs);
760 }
761 }
762}
763
764const G4Colour& G4VSceneHandler::GetColour (const G4Visible& visible) {
765 // Colour is determined by the applicable vis attributes.
766 const G4Colour& colour = fpViewer ->
767 GetApplicableVisAttributes (visible.GetVisAttributes ()) -> GetColour ();
768 return colour;
769}
770
771const G4Colour& G4VSceneHandler::GetTextColour (const G4Text& text) {
772 const G4VisAttributes* pVA = text.GetVisAttributes ();
773 if (!pVA) {
774 pVA = fpViewer -> GetViewParameters (). GetDefaultTextVisAttributes ();
775 }
776 const G4Colour& colour = pVA -> GetColour ();
777 return colour;
778}
779
780G4double G4VSceneHandler::GetLineWidth(const G4VisAttributes* pVisAttribs)
781{
782 G4double lineWidth = pVisAttribs->GetLineWidth();
783 if (lineWidth < 1.) lineWidth = 1.;
784 lineWidth *= fpViewer -> GetViewParameters().GetGlobalLineWidthScale();
785 if (lineWidth < 1.) lineWidth = 1.;
786 return lineWidth;
787}
788
789G4ViewParameters::DrawingStyle G4VSceneHandler::GetDrawingStyle
790(const G4VisAttributes* pVisAttribs) {
791 // Drawing style is normally determined by the view parameters, but
792 // it can be overriddden by the ForceDrawingStyle flag in the vis
793 // attributes.
794 G4ViewParameters::DrawingStyle style =
795 fpViewer->GetViewParameters().GetDrawingStyle();
796 if (pVisAttribs -> IsForceDrawingStyle ()) {
797 G4VisAttributes::ForcedDrawingStyle forcedStyle =
798 pVisAttribs -> GetForcedDrawingStyle ();
799 // This is complicated because if hidden line and surface removal
800 // has been requested we wish to preserve this sometimes.
801 switch (forcedStyle) {
802 case (G4VisAttributes::solid):
803 switch (style) {
804 case (G4ViewParameters::hlr):
805 style = G4ViewParameters::hlhsr;
806 break;
807 case (G4ViewParameters::wireframe):
808 style = G4ViewParameters::hsr;
809 break;
810 case (G4ViewParameters::hlhsr):
811 case (G4ViewParameters::hsr):
812 default:
813 break;
814 }
815 break;
816 case (G4VisAttributes::wireframe):
817 default:
818 // But if forced style is wireframe, do it, because one of its
819 // main uses is in displaying the consituent solids of a Boolean
820 // solid and their surfaces overlap with the resulting Booean
821 // solid, making a mess if hlr is specified.
822 style = G4ViewParameters::wireframe;
823 break;
824 }
825 }
826 return style;
827}
828
829G4bool G4VSceneHandler::GetAuxEdgeVisible (const G4VisAttributes* pVisAttribs) {
830 G4bool isAuxEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible ();
831 if (pVisAttribs -> IsForceAuxEdgeVisible()) isAuxEdgeVisible = true;
832 return isAuxEdgeVisible;
833}
834
835G4double G4VSceneHandler::GetMarkerSize
836(const G4VMarker& marker,
837 G4VSceneHandler::MarkerSizeType& markerSizeType)
838{
839 G4bool userSpecified = marker.GetWorldSize() || marker.GetScreenSize();
840 const G4VMarker& defaultMarker =
841 fpViewer -> GetViewParameters().GetDefaultMarker();
842 G4double size = userSpecified ?
843 marker.GetWorldSize() : defaultMarker.GetWorldSize();
844 if (size) {
845 // Draw in world coordinates.
846 markerSizeType = world;
847 }
848 else {
849 size = userSpecified ?
850 marker.GetScreenSize() : defaultMarker.GetScreenSize();
851 // Draw in screen coordinates.
852 markerSizeType = screen;
853 }
854 if (size <= 1.) size = 1.;
855 size *= fpViewer -> GetViewParameters().GetGlobalMarkerScale();
856 if (size <= 1.) size = 1.;
857 return size;
858}
859
860G4int G4VSceneHandler::GetNoOfSides(const G4VisAttributes* pVisAttribs)
861{
862 // No. of sides (lines segments per circle) is normally determined
863 // by the view parameters, but it can be overriddden by the
864 // ForceLineSegmentsPerCircle in the vis attributes.
865 G4int lineSegmentsPerCircle = fpViewer->GetViewParameters().GetNoOfSides();
866 if (pVisAttribs) {
867 if (pVisAttribs->IsForceLineSegmentsPerCircle())
868 lineSegmentsPerCircle = pVisAttribs->GetForcedLineSegmentsPerCircle();
869 const G4int nSegmentsMin = 12;
870 if (lineSegmentsPerCircle < nSegmentsMin) {
871 lineSegmentsPerCircle = nSegmentsMin;
872 G4cout <<
873 "G4VSceneHandler::GetNoOfSides: attempt to set the"
874 "\nnumber of line segements per circle < " << nSegmentsMin
875 << "; forced to " << lineSegmentsPerCircle << G4endl;
876 }
877 }
878 return lineSegmentsPerCircle;
879}
880
881std::ostream& operator << (std::ostream& os, const G4VSceneHandler& s) {
882
883 os << "Scene handler " << s.fName << " has "
884 << s.fViewerList.size () << " viewer(s):";
885 for (size_t i = 0; i < s.fViewerList.size (); i++) {
886 os << "\n " << *(s.fViewerList [i]);
887 }
888
889 if (s.fpScene) {
890 os << "\n " << *s.fpScene;
891 }
892 else {
893 os << "\n This scene handler currently has no scene.";
894 }
895
896 return os;
897}
Note: See TracBrowser for help on using the repository browser.