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

Last change on this file since 1176 was 1168, checked in by garnier, 16 years ago

Fix a delete problem

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