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

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

r660@mac-90108: laurentgarnier | 2007-06-25 16:10:12 +0200
ajout de fichiers NON modifies

  • Property svn:mime-type set to text/cpp
File size: 26.3 KB
RevLine 
[531]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.77 2006/11/21 14:23:53 allison Exp $
28// GEANT4 tag $Name: geant4-08-02-patch-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 "G4VHit.hh"
74#include "Randomize.hh"
75#include "G4StateManager.hh"
76#include "G4RunManager.hh"
77#include "G4Run.hh"
78#include "G4Transform3D.hh"
79
80G4VSceneHandler::G4VSceneHandler (G4VGraphicsSystem& system, G4int id, const G4String& name):
81 fSystem (system),
82 fSceneHandlerId (id),
83 fViewCount (0),
84 fpViewer (0),
85 fpScene (0),
86 fMarkForClearingTransientStore (true), // Ready for first
87 // ClearTransientStoreIfMarked(),
88 // e.g., at end of run (see
89 // G4VisManager.cc).
90 fReadyForTransients (true), // Only false while processing scene.
91 fProcessingSolid (false),
92 fSecondPassRequested (false),
93 fSecondPass (false),
94 fpModel (0),
95 fpObjectTransformation (0),
96 fNestingDepth (0),
97 fpVisAttribs (0),
98 fRequestedEvent (0)
99{
100 G4VisManager* pVMan = G4VisManager::GetInstance ();
101 fpScene = pVMan -> GetCurrentScene ();
102 if (name == "") {
103 std::ostringstream ost;
104 ost << fSystem.GetName () << '-' << fSceneHandlerId;
105 fName = ost.str();
106 }
107 else {
108 fName = name;
109 }
110 fTransientsDrawnThisEvent = pVMan->GetTransientsDrawnThisEvent();
111 fTransientsDrawnThisRun = pVMan->GetTransientsDrawnThisRun();
112}
113
114G4VSceneHandler::~G4VSceneHandler () {
115 G4ViewerListIterator i;
116 for (i = fViewerList.begin(); i != fViewerList.end(); ++i) {
117 delete *i;
118 }
119}
120
121void G4VSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
122 const G4VisAttributes& visAttribs) {
123 fpObjectTransformation = &objectTransformation;
124 fpVisAttribs = &visAttribs;
125 fProcessingSolid = true;
126}
127
128void G4VSceneHandler::PostAddSolid () {
129 fpObjectTransformation = 0;
130 fpVisAttribs = 0;
131 fProcessingSolid = false;
132 if (fReadyForTransients) {
133 fTransientsDrawnThisEvent = true;
134 fTransientsDrawnThisRun = true;
135 }
136}
137
138void G4VSceneHandler::BeginPrimitives
139(const G4Transform3D& objectTransformation) {
140 fNestingDepth++;
141 if (fNestingDepth > 1)
142 G4Exception("G4VSceneHandler::BeginPrimitives: Nesting detected."
143 "\n It is illegal to nest Begin/EndPrimitives.");
144 fpObjectTransformation = &objectTransformation;
145}
146
147void G4VSceneHandler::EndPrimitives () {
148 if (fNestingDepth <= 0)
149 G4Exception("G4VSceneHandler::EndPrimitives: Nesting error");
150 fNestingDepth--;
151 fpObjectTransformation = 0;
152 if (fReadyForTransients) {
153 fTransientsDrawnThisEvent = true;
154 fTransientsDrawnThisRun = true;
155 }
156}
157
158void G4VSceneHandler::BeginPrimitives2D () {
159 fNestingDepth++;
160 if (fNestingDepth > 1)
161 G4Exception("G4VSceneHandler::BeginPrimitives2D: Nesting detected."
162 "\n It is illegal to nest Begin/EndPrimitives.");
163 // Not actually required for 2D operations but some drivers do an
164 // initial transformation...
165 fpObjectTransformation = &fIdentityTransformation;
166}
167
168void G4VSceneHandler::EndPrimitives2D () {
169 if (fNestingDepth <= 0)
170 G4Exception("G4VSceneHandler::EndPrimitives2D: Nesting error");
171 fNestingDepth--;
172 fpObjectTransformation = 0;
173 if (fReadyForTransients) {
174 fTransientsDrawnThisEvent = true;
175 fTransientsDrawnThisRun = true;
176 }
177}
178
179void G4VSceneHandler::BeginModeling () {
180}
181
182void G4VSceneHandler::EndModeling ()
183{
184 fpModel = 0;
185}
186
187void G4VSceneHandler::ClearStore () {
188 // if (fpViewer) fpViewer -> NeedKernelVisit (true);
189 // ?? Viewer is supposed to be smart enough to know when to visit
190 // kernel, but a problem in OpenGL Stored seems to require a forced
191 // kernel visit triggered by the above code. John Allison Aug 2001
192 // Feb 2005 - commented out. Let's fix OpenGL if necessary.
193}
194
195void G4VSceneHandler::ClearTransientStore () {
196}
197
198void G4VSceneHandler::AddSolid (const G4Box& box) {
199 RequestPrimitives (box);
200// If your graphics system is sophisticated enough to handle a
201// particular solid shape as a primitive, in your derived class write a
202// function to override this. (Note: some compilers warn that your
203// function "hides" this one. That's OK.)
204// Your function might look like this...
205// void G4MyScene::AddSolid (const G4Box& box) {
206// Get parameters of appropriate object, e.g.:
207// G4double dx = box.GetXHalfLength ();
208// G4double dy = box.GetYHalfLength ();
209// G4double dz = box.GetZHalfLength ();
210// and Draw or Store in your display List.
211}
212
213void G4VSceneHandler::AddSolid (const G4Tubs& tubs) {
214 RequestPrimitives (tubs);
215}
216
217void G4VSceneHandler::AddSolid (const G4Cons& cons) {
218 RequestPrimitives (cons);
219}
220
221void G4VSceneHandler::AddSolid (const G4Trd& trd) {
222 RequestPrimitives (trd);
223}
224
225void G4VSceneHandler::AddSolid (const G4Trap& trap) {
226 RequestPrimitives (trap);
227}
228
229void G4VSceneHandler::AddSolid (const G4Sphere& sphere) {
230 RequestPrimitives (sphere );
231}
232
233void G4VSceneHandler::AddSolid (const G4Para& para) {
234 RequestPrimitives (para);
235}
236
237void G4VSceneHandler::AddSolid (const G4Torus& torus) {
238 RequestPrimitives (torus);
239}
240
241void G4VSceneHandler::AddSolid (const G4Polycone& polycone) {
242 RequestPrimitives (polycone);
243}
244
245void G4VSceneHandler::AddSolid (const G4Polyhedra& polyhedra) {
246 RequestPrimitives (polyhedra);
247}
248
249void G4VSceneHandler::AddSolid (const G4VSolid& solid) {
250 RequestPrimitives (solid);
251}
252
253void G4VSceneHandler::AddCompound (const G4VTrajectory& traj) {
254 G4TrajectoriesModel* pTrModel =
255 dynamic_cast<G4TrajectoriesModel*>(fpModel);
256 if (!pTrModel) G4Exception
257 ("G4VSceneHandler::AddCompound(const G4VTrajectory&): Not a G4TrajectoriesModel.");
258 traj.DrawTrajectory(pTrModel->GetDrawingMode());
259}
260
261void G4VSceneHandler::AddCompound (const G4VHit& hit) {
262 ((G4VHit&)hit).Draw(); // Cast to non-const because Draw is non-const!!!!
263}
264
265void G4VSceneHandler::AddViewerToList (G4VViewer* pViewer) {
266 fViewerList.push_back (pViewer);
267}
268
269void G4VSceneHandler::AddPrimitive (const G4Scale& scale) {
270
271 const G4double margin(0.01);
272 // Fractional margin - ensures scale is comfortably inside viewing
273 // volume.
274 const G4double oneMinusMargin (1. - margin);
275
276 const G4VisExtent& sceneExtent = fpScene->GetExtent();
277
278 // Useful constants...
279 const G4double length(scale.GetLength());
280 const G4double halfLength(length / 2.);
281 const G4double tickLength(length / 20.);
282 const G4double piBy2(halfpi);
283
284 // Get size of scene...
285 const G4double xmin = sceneExtent.GetXmin();
286 const G4double xmax = sceneExtent.GetXmax();
287 const G4double ymin = sceneExtent.GetYmin();
288 const G4double ymax = sceneExtent.GetYmax();
289 const G4double zmin = sceneExtent.GetZmin();
290 const G4double zmax = sceneExtent.GetZmax();
291
292 // Create (empty) polylines having the same vis attributes...
293 G4Polyline scaleLine, tick11, tick12, tick21, tick22;
294 G4VisAttributes visAtts(*scale.GetVisAttributes()); // Long enough life.
295 scaleLine.SetVisAttributes(&visAtts);
296 tick11.SetVisAttributes(&visAtts);
297 tick12.SetVisAttributes(&visAtts);
298 tick21.SetVisAttributes(&visAtts);
299 tick22.SetVisAttributes(&visAtts);
300
301 // Add points to the polylines to represent an scale parallel to the
302 // x-axis centred on the origin...
303 G4Point3D r1(G4Point3D(-halfLength, 0., 0.));
304 G4Point3D r2(G4Point3D( halfLength, 0., 0.));
305 scaleLine.push_back(r1);
306 scaleLine.push_back(r2);
307 G4Point3D ticky(0., tickLength, 0.);
308 G4Point3D tickz(0., 0., tickLength);
309 tick11.push_back(r1 + ticky);
310 tick11.push_back(r1 - ticky);
311 tick12.push_back(r1 + tickz);
312 tick12.push_back(r1 - tickz);
313 tick21.push_back(r2 + ticky);
314 tick21.push_back(r2 - ticky);
315 tick22.push_back(r2 + tickz);
316 tick22.push_back(r2 - tickz);
317 G4Point3D textPosition(0., tickLength, 0.);
318
319 // Transform appropriately...
320
321 G4Transform3D transformation;
322 if (scale.GetAutoPlacing()) {
323 G4Transform3D rotation;
324 switch (scale.GetDirection()) {
325 case G4Scale::x:
326 break;
327 case G4Scale::y:
328 rotation = G4RotateZ3D(piBy2);
329 break;
330 case G4Scale::z:
331 rotation = G4RotateY3D(piBy2);
332 break;
333 }
334 G4double sxmid(scale.GetXmid());
335 G4double symid(scale.GetYmid());
336 G4double szmid(scale.GetZmid());
337 sxmid = xmin + oneMinusMargin * (xmax - xmin);
338 symid = ymin + margin * (ymax - ymin);
339 szmid = zmin + oneMinusMargin * (zmax - zmin);
340 switch (scale.GetDirection()) {
341 case G4Scale::x:
342 sxmid -= halfLength;
343 break;
344 case G4Scale::y:
345 symid += halfLength;
346 break;
347 case G4Scale::z:
348 szmid -= halfLength;
349 break;
350 }
351 G4Translate3D translation(sxmid, symid, szmid);
352 transformation = translation * rotation;
353 } else {
354 if (fpModel) transformation = fpModel->GetTransformation();
355 }
356
357 // Draw...
358 // We would like to call BeginPrimitives(transformation) here but
359 // calling BeginPrimitives from within an AddPrimitive is not
360 // allowed! So we have to do our own transformation...
361 AddPrimitive(scaleLine.transform(transformation));
362 AddPrimitive(tick11.transform(transformation));
363 AddPrimitive(tick12.transform(transformation));
364 AddPrimitive(tick21.transform(transformation));
365 AddPrimitive(tick22.transform(transformation));
366 G4Text text(scale.GetAnnotation(),textPosition.transform(transformation));
367 text.SetScreenSize(12.);
368 AddPrimitive(text);
369}
370
371void G4VSceneHandler::AddPrimitive (const G4Polymarker& polymarker) {
372 switch (polymarker.GetMarkerType()) {
373 default:
374 case G4Polymarker::dots:
375 {
376 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
377 G4Circle dot (polymarker);
378 dot.SetPosition (polymarker[iPoint]);
379 dot.SetWorldSize (0.);
380 dot.SetScreenSize (0.1); // Very small circle.
381 AddPrimitive (dot);
382 }
383 }
384 break;
385 case G4Polymarker::circles:
386 {
387 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
388 G4Circle circle (polymarker);
389 circle.SetPosition (polymarker[iPoint]);
390 AddPrimitive (circle);
391 }
392 }
393 break;
394 case G4Polymarker::squares:
395 {
396 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
397 G4Square square (polymarker);
398 square.SetPosition (polymarker[iPoint]);
399 AddPrimitive (square);
400 }
401 }
402 break;
403 }
404}
405
406void G4VSceneHandler::RemoveViewerFromList (G4VViewer* pViewer) {
407 fViewerList.remove(pViewer);
408}
409
410void G4VSceneHandler::SetScene (G4Scene* pScene) {
411 fpScene = pScene;
412 // Notify all viewers that a kernel visit is required.
413 G4ViewerListIterator i;
414 for (i = fViewerList.begin(); i != fViewerList.end(); i++) {
415 (*i) -> SetNeedKernelVisit (true);
416 }
417}
418
419void G4VSceneHandler::RequestPrimitives (const G4VSolid& solid) {
420 BeginPrimitives (*fpObjectTransformation);
421 G4NURBS* pNURBS = 0;
422 G4Polyhedron* pPolyhedron = 0;
423 const G4VisAttributes* pVisAttribs =
424 fpViewer -> GetApplicableVisAttributes (fpVisAttribs);
425 switch (fpViewer -> GetViewParameters () . GetRepStyle ()) {
426 case G4ViewParameters::nurbs:
427 pNURBS = solid.CreateNURBS ();
428 if (pNURBS) {
429 pNURBS -> SetVisAttributes
430 (fpViewer -> GetApplicableVisAttributes (pVisAttribs));
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 (pVisAttribs));
450 pPolyhedron = solid.GetPolyhedron ();
451 G4Polyhedron::ResetNumberOfRotationSteps ();
452 if (pPolyhedron) {
453 pPolyhedron -> SetVisAttributes (pVisAttribs);
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 (fRequestedEvent) {
552 DrawEvent(fRequestedEvent);
553
554 } else {
555
556 G4RunManager* runManager = G4RunManager::GetRunManager();
557 const G4Run* run = runManager->GetCurrentRun();
558 const std::vector<const G4Event*>* events =
559 run? run->GetEventVector(): 0;
560 size_t nKeptEvents = 0;
561 if (events) nKeptEvents = events->size();
562 if (runManager) {
563 if (fpScene->GetRefreshAtEndOfEvent()) {
564
565 if (verbosity >= G4VisManager::confirmations) {
566 G4cout << "Refreshing event..." << G4endl;
567 }
568 const G4Event* event = 0;
569 if (events && events->size()) event = events->back();
570 if (event) DrawEvent(event);
571
572 } else { // Accumulating events.
573
574 if (verbosity >= G4VisManager::confirmations) {
575 G4cout << "Refreshing events in run..." << G4endl;
576 }
577 for (size_t i = 0; i < nKeptEvents; ++i) {
578 const G4Event* event = (*events)[i];
579 if (event) DrawEvent(event);
580 }
581
582 if (!fpScene->GetRefreshAtEndOfRun()) {
583 if (verbosity >= G4VisManager::warnings) {
584 G4cout <<
585 "WARNING: Cannot refresh events accumulated over more"
586 "\n than one runs. Refreshed just the last run..."
587 << G4endl;
588 }
589 }
590 }
591 }
592 }
593 visManager->SetEventRefreshing(false);
594 }
595
596 fMarkForClearingTransientStore = tmpMarkForClearingTransientStore;
597}
598
599void G4VSceneHandler::DrawEvent(const G4Event* event)
600{
601 const std::vector<G4VModel*>& EOEModelList =
602 fpScene -> GetEndOfEventModelList ();
603 size_t nModels = EOEModelList.size();
604 if (nModels) {
605 G4ModelingParameters* pMP = CreateModelingParameters();
606 pMP->SetEvent(event);
607 for (size_t i = 0; i < nModels; i++) {
608 G4VModel* pModel = EOEModelList [i];
609 pModel -> SetModelingParameters(pMP);
610 SetModel (pModel);
611 pModel -> DescribeYourselfTo (*this);
612 pModel -> SetModelingParameters(0);
613 }
614 delete pMP;
615 SetModel (0);
616 }
617}
618
619G4ModelingParameters* G4VSceneHandler::CreateModelingParameters ()
620{
621 // Create modeling parameters from View Parameters...
622 const G4ViewParameters& vp = fpViewer -> GetViewParameters ();
623
624 // Convert drawing styles...
625 G4ModelingParameters::DrawingStyle modelDrawingStyle =
626 G4ModelingParameters::wf;
627 switch (vp.GetDrawingStyle ()) {
628 default:
629 case G4ViewParameters::wireframe:
630 modelDrawingStyle = G4ModelingParameters::wf;
631 break;
632 case G4ViewParameters::hlr:
633 modelDrawingStyle = G4ModelingParameters::hlr;
634 break;
635 case G4ViewParameters::hsr:
636 modelDrawingStyle = G4ModelingParameters::hsr;
637 break;
638 case G4ViewParameters::hlhsr:
639 modelDrawingStyle = G4ModelingParameters::hlhsr;
640 break;
641 }
642
643 // Decide if covered daughters are really to be culled...
644 G4bool reallyCullCovered =
645 vp.IsCullingCovered() // Culling daughters depends also on...
646 && !vp.IsSection () // Sections (DCUT) not requested.
647 && !vp.IsCutaway () // Cutaways not requested.
648 ;
649
650 G4ModelingParameters* pModelingParams = new G4ModelingParameters
651 (vp.GetDefaultVisAttributes (),
652 modelDrawingStyle,
653 vp.IsCulling (),
654 vp.IsCullingInvisible (),
655 vp.IsDensityCulling (),
656 vp.GetVisibleDensity (),
657 reallyCullCovered,
658 vp.GetNoOfSides ()
659 );
660
661 pModelingParams->SetWarning
662 (G4VisManager::GetInstance()->GetVerbosity() >= G4VisManager::warnings);
663
664 pModelingParams->SetExplodeFactor(vp.GetExplodeFactor());
665 pModelingParams->SetExplodeCentre(vp.GetExplodeCentre());
666
667 pModelingParams->SetSectionPolyhedron(CreateSectionPolyhedron());
668 pModelingParams->SetCutawayPolyhedron(CreateCutawayPolyhedron());
669 // The polyhedron objects are deleted in the modeling parameters destructor.
670
671 return pModelingParams;
672}
673
674const G4Polyhedron* G4VSceneHandler::CreateSectionPolyhedron()
675{
676 /* Disable for now. Boolean processor not up to it.
677 const G4ViewParameters& vp = fpViewer->GetViewParameters();
678 if (vp.IsSection () ) {
679 G4double radius = fpScene->GetExtent().GetExtentRadius();
680 G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
681 G4Box sectionBox("clipper",
682 safe, safe, 1.e-5 * radius); // Thin in z-plane.
683 G4Polyhedron* sectioner = sectionBox.CreatePolyhedron();
684 const G4Plane3D& s = vp.GetSectionPlane ();
685 G4double a = s.a();
686 G4double b = s.b();
687 G4double c = s.c();
688 G4double d = s.d();
689 G4Transform3D transform = G4TranslateZ3D(-d);
690 const G4Normal3D normal(a,b,c);
691 if (normal != G4Normal3D(0,0,1)) {
692 const G4double angle = std::acos(normal.dot(G4Normal3D(0,0,1)));
693 const G4Vector3D axis = G4Normal3D(0,0,1).cross(normal);
694 transform = G4Rotate3D(angle, axis) * transform;
695 }
696 sectioner->Transform(transform);
697 return sectioner;
698 } else {
699 return 0;
700 }
701 */
702 return 0;
703}
704
705const G4Polyhedron* G4VSceneHandler::CreateCutawayPolyhedron()
706{
707 return 0;
708}
709
710const G4Colour& G4VSceneHandler::GetColour (const G4Visible& visible) {
711 // Colour is determined by the applicable vis attributes.
712 const G4Colour& colour = fpViewer ->
713 GetApplicableVisAttributes (visible.GetVisAttributes ()) -> GetColour ();
714 return colour;
715}
716
717const G4Colour& G4VSceneHandler::GetTextColour (const G4Text& text) {
718 const G4VisAttributes* pVA = text.GetVisAttributes ();
719 if (!pVA) {
720 pVA = fpViewer -> GetViewParameters (). GetDefaultTextVisAttributes ();
721 }
722 const G4Colour& colour = pVA -> GetColour ();
723 return colour;
724}
725
726G4double G4VSceneHandler::GetLineWidth(const G4Visible& visible)
727{
728 G4double lineWidth = fpViewer->
729 GetApplicableVisAttributes(visible.GetVisAttributes())->GetLineWidth();
730 if (lineWidth < 1.) lineWidth = 1.;
731 lineWidth *= fpViewer -> GetViewParameters().GetGlobalLineWidthScale();
732 if (lineWidth < 1.) lineWidth = 1.;
733 return lineWidth;
734}
735
736G4ViewParameters::DrawingStyle G4VSceneHandler::GetDrawingStyle
737(const G4VisAttributes* pVisAttribs) {
738 // Drawing style is normally determined by the view parameters, but
739 // it can be overriddden by the ForceDrawingStyle flag in the vis
740 // attributes.
741 G4ViewParameters::DrawingStyle style =
742 fpViewer->GetViewParameters().GetDrawingStyle();
743 if (pVisAttribs -> IsForceDrawingStyle ()) {
744 G4VisAttributes::ForcedDrawingStyle forcedStyle =
745 pVisAttribs -> GetForcedDrawingStyle ();
746 // This is complicated because if hidden line and surface removal
747 // has been requested we wish to preserve this sometimes.
748 switch (forcedStyle) {
749 case (G4VisAttributes::solid):
750 switch (style) {
751 case (G4ViewParameters::hlr):
752 style = G4ViewParameters::hlhsr;
753 break;
754 case (G4ViewParameters::wireframe):
755 style = G4ViewParameters::hsr;
756 break;
757 case (G4ViewParameters::hlhsr):
758 case (G4ViewParameters::hsr):
759 default:
760 break;
761 }
762 break;
763 case (G4VisAttributes::wireframe):
764 default:
765 // But if forced style is wireframe, do it, because one of its
766 // main uses is in displaying the consituent solids of a Boolean
767 // solid and their surfaces overlap with the resulting Booean
768 // solid, making a mess if hlr is specified.
769 style = G4ViewParameters::wireframe;
770 break;
771 }
772 }
773 return style;
774}
775
776G4bool G4VSceneHandler::GetAuxEdgeVisible (const G4VisAttributes* pVisAttribs) {
777 G4bool isAuxEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible ();
778 if (pVisAttribs -> IsForceAuxEdgeVisible()) isAuxEdgeVisible = true;
779 return isAuxEdgeVisible;
780}
781
782G4double G4VSceneHandler::GetMarkerSize
783(const G4VMarker& marker,
784 G4VSceneHandler::MarkerSizeType& markerSizeType)
785{
786 G4bool userSpecified = marker.GetWorldSize() || marker.GetScreenSize();
787 const G4VMarker& defaultMarker =
788 fpViewer -> GetViewParameters().GetDefaultMarker();
789 G4double size = userSpecified ?
790 marker.GetWorldSize() : defaultMarker.GetWorldSize();
791 if (size) {
792 // Draw in world coordinates.
793 markerSizeType = world;
794 }
795 else {
796 size = userSpecified ?
797 marker.GetScreenSize() : defaultMarker.GetScreenSize();
798 // Draw in screen coordinates.
799 markerSizeType = screen;
800 }
801 if (size <= 1.) size = 1.;
802 size *= fpViewer -> GetViewParameters().GetGlobalMarkerScale();
803 if (size <= 1.) size = 1.;
804 return size;
805}
806
807G4int G4VSceneHandler::GetNoOfSides(const G4VisAttributes* pVisAttribs)
808{
809 // No. of sides (lines segments per circle) is normally determined
810 // by the view parameters, but it can be overriddden by the
811 // ForceLineSegmentsPerCircle in the vis attributes.
812 G4int lineSegmentsPerCircle = fpViewer->GetViewParameters().GetNoOfSides();
813 if (pVisAttribs->GetForcedLineSegmentsPerCircle() > 0)
814 lineSegmentsPerCircle = pVisAttribs->GetForcedLineSegmentsPerCircle();
815 const G4int nSegmentsMin = 12;
816 if (lineSegmentsPerCircle < nSegmentsMin) {
817 lineSegmentsPerCircle = nSegmentsMin;
818 G4cout <<
819 "G4VSceneHandler::GetNoOfSides: attempt to set the"
820 "\nnumber of line segements per circle < " << nSegmentsMin
821 << "; forced to " << lineSegmentsPerCircle << G4endl;
822 }
823 return lineSegmentsPerCircle;
824}
825
826std::ostream& operator << (std::ostream& os, const G4VSceneHandler& s) {
827
828 os << "Scene handler " << s.fName << " has "
829 << s.fViewerList.size () << " viewer(s):";
830 for (size_t i = 0; i < s.fViewerList.size (); i++) {
831 os << "\n " << *(s.fViewerList [i]);
832 }
833
834 if (s.fpScene) {
835 os << "\n " << *s.fpScene;
836 }
837 else {
838 os << "\n This scene handler currently has no scene.";
839 }
840
841 return os;
842}
Note: See TracBrowser for help on using the repository browser.