source: trunk/source/visualization/management/src/G4VSceneHandler.cc@ 1301

Last change on this file since 1301 was 1288, checked in by garnier, 15 years ago

update from CVS

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