Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VSceneHandler.cc
Go to the documentation of this file.
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//
28//
29// John Allison 19th July 1996
30// Abstract interface class for graphics scenes.
31
32#include "G4VSceneHandler.hh"
33
34#include "G4ios.hh"
35#include <sstream>
36
37#include "G4VisManager.hh"
38#include "G4VGraphicsSystem.hh"
39#include "G4VViewer.hh"
40#include "G4VSolid.hh"
41#include "G4RotationMatrix.hh"
42#include "G4ThreeVector.hh"
43#include "G4VPhysicalVolume.hh"
44#include "G4Material.hh"
45#include "G4Polyline.hh"
46#include "G4Text.hh"
47#include "G4Circle.hh"
48#include "G4Square.hh"
49#include "G4Polymarker.hh"
50#include "G4Polyhedron.hh"
51#include "G4Visible.hh"
52#include "G4VisAttributes.hh"
53#include "G4VModel.hh"
55#include "G4Box.hh"
56#include "G4Cons.hh"
57#include "G4Orb.hh"
58#include "G4Para.hh"
59#include "G4Sphere.hh"
60#include "G4Torus.hh"
61#include "G4Trap.hh"
62#include "G4Trd.hh"
63#include "G4Tubs.hh"
64#include "G4Ellipsoid.hh"
65#include "G4Polycone.hh"
66#include "G4Polyhedra.hh"
67#include "G4DisplacedSolid.hh"
68#include "G4LogicalVolume.hh"
71#include "G4VTrajectory.hh"
72#include "G4VTrajectoryPoint.hh"
73#include "G4HitsModel.hh"
74#include "G4VHit.hh"
75#include "G4VDigi.hh"
76#include "G4ScoringManager.hh"
77#include "G4VScoringMesh.hh"
78#include "G4Mesh.hh"
80#include "Randomize.hh"
81#include "G4StateManager.hh"
82#include "G4RunManager.hh"
84#include "G4Run.hh"
85#include "G4Transform3D.hh"
86#include "G4AttHolder.hh"
87#include "G4AttDef.hh"
88#include "G4VVisCommand.hh"
90#include "G4SystemOfUnits.hh"
91
92#include <set>
93
95 fSystem (system),
96 fSceneHandlerId (id),
97 fViewCount (0),
98 fpViewer (0),
99 fpScene (0),
100 fMarkForClearingTransientStore (true), // Ready for first
101 // ClearTransientStoreIfMarked(),
102 // e.g., at end of run (see
103 // G4VisManager.cc).
104 fReadyForTransients (true), // Only false while processing scene.
105 fProcessingSolid (false),
106 fProcessing2D (false),
107 fpModel (0),
108 fNestingDepth (0),
109 fpVisAttribs (0)
112 fpScene = pVMan -> GetCurrentScene ();
113 if (name == "") {
114 std::ostringstream ost;
115 ost << fSystem.GetName () << '-' << fSceneHandlerId;
116 fName = ost.str();
117 }
118 else {
119 fName = name;
120 }
121 fTransientsDrawnThisEvent = pVMan->GetTransientsDrawnThisEvent();
122 fTransientsDrawnThisRun = pVMan->GetTransientsDrawnThisRun();
123}
124
126 G4VViewer* last;
127 while( ! fViewerList.empty() ) {
128 last = fViewerList.back();
129 fViewerList.pop_back();
130 delete last;
133
135{
136 if (fpScene) {
137 return fpScene->GetExtent();
138 } else {
139 static const G4VisExtent defaultExtent = G4VisExtent();
140 return defaultExtent;
141 }
142}
143
144void G4VSceneHandler::PreAddSolid (const G4Transform3D& objectTransformation,
145 const G4VisAttributes& visAttribs) {
146 fObjectTransformation = objectTransformation;
147 fpVisAttribs = &visAttribs;
148 fProcessingSolid = true;
149}
150
152 fpVisAttribs = 0;
153 fProcessingSolid = false;
157 }
158}
159
161(const G4Transform3D& objectTransformation) {
162 //static G4int count = 0;
163 //G4cout << "G4VSceneHandler::BeginPrimitives: " << count++ << G4endl;
165 if (fNestingDepth > 1)
167 ("G4VSceneHandler::BeginPrimitives",
168 "visman0101", FatalException,
169 "Nesting detected. It is illegal to nest Begin/EndPrimitives.");
170 fObjectTransformation = objectTransformation;
171}
172
174 if (fNestingDepth <= 0)
175 G4Exception("G4VSceneHandler::EndPrimitives",
176 "visman0102", FatalException, "Nesting error.");
181 }
182}
183
185(const G4Transform3D& objectTransformation) {
187 if (fNestingDepth > 1)
189 ("G4VSceneHandler::BeginPrimitives2D",
190 "visman0103", FatalException,
191 "Nesting detected. It is illegal to nest Begin/EndPrimitives.");
192 fObjectTransformation = objectTransformation;
193 fProcessing2D = true;
194}
195
197 if (fNestingDepth <= 0)
198 G4Exception("G4VSceneHandler::EndPrimitives2D",
199 "visman0104", FatalException, "Nesting error.");
204 }
205 fProcessing2D = false;
206}
207
209}
210
212{
213 fpModel = 0;
214}
215
217
219
220template <class T> void G4VSceneHandler::AddSolidT
221(const T& solid)
222{
223 // Get and check applicable vis attributes.
225 RequestPrimitives (solid);
226}
227
229(const T& solid)
230{
231 // Get and check applicable vis attributes.
233 // Draw with auxiliary edges unless otherwise specified.
235 // Create a vis atts object for the modified vis atts.
236 // It is static so that we may return a reliable pointer to it.
237 static G4VisAttributes visAttsWithAuxEdges;
238 // Initialise it with the current vis atts and reset the pointer.
239 visAttsWithAuxEdges = *fpVisAttribs;
240 // Force auxiliary edges visible.
241 visAttsWithAuxEdges.SetForceAuxEdgeVisible();
242 fpVisAttribs = &visAttsWithAuxEdges;
243 }
244 RequestPrimitives (solid);
245}
246
248 AddSolidT (box);
249 // If your graphics system is sophisticated enough to handle a
250 // particular solid shape as a primitive, in your derived class write a
251 // function to override this.
252 // Your function might look like this...
253 // void G4MySceneHandler::AddSolid (const G4Box& box) {
254 // Get and check applicable vis attributes.
255 // fpVisAttribs = fpViewer->GetApplicableVisAttributes(fpVisAttribs);
256 // Do not draw if not visible.
257 // if (fpVisAttribs->IsVisible()) {
258 // Get parameters of appropriate object, e.g.:
259 // G4double dx = box.GetXHalfLength ();
260 // G4double dy = box.GetYHalfLength ();
261 // G4double dz = box.GetZHalfLength ();
262 // ...
263 // and Draw or Store in your display List.
264}
265
267 AddSolidT (cons);
268}
269
272}
273
275 AddSolidT (para);
276}
277
280}
281
284}
285
287 AddSolidT (trap);
288}
289
291 AddSolidT (trd);
292}
293
295 AddSolidT (tubs);
296}
297
298void G4VSceneHandler::AddSolid (const G4Ellipsoid& ellipsoid) {
299 AddSolidWithAuxiliaryEdges (ellipsoid);
300}
301
302void G4VSceneHandler::AddSolid (const G4Polycone& polycone) {
303 AddSolidT (polycone);
304}
305
306void G4VSceneHandler::AddSolid (const G4Polyhedra& polyhedra) {
307 AddSolidT (polyhedra);
308}
309
311 AddSolidT (tess);
312}
313
315 AddSolidT (solid);
316}
317
319 G4TrajectoriesModel* trajectoriesModel =
320 dynamic_cast<G4TrajectoriesModel*>(fpModel);
321 if (trajectoriesModel)
322 traj.DrawTrajectory();
323 else {
325 ("G4VSceneHandler::AddCompound(const G4VTrajectory&)",
326 "visman0105", FatalException, "Not a G4TrajectoriesModel.");
327 }
328}
329
331 // Cast away const because Draw is non-const!!!!
332 const_cast<G4VHit&>(hit).Draw();
333}
334
336 // Cast away const because Draw is non-const!!!!
337 const_cast<G4VDigi&>(digi).Draw();
338}
339
341 using MeshScoreMap = G4VScoringMesh::MeshScoreMap;
342 //G4cout << "AddCompound: hits: " << &hits << G4endl;
343 G4bool scoreMapHits = false;
345 if (scoringManager) {
346 size_t nMeshes = scoringManager->GetNumberOfMesh();
347 for (size_t iMesh = 0; iMesh < nMeshes; ++iMesh) {
348 G4VScoringMesh* mesh = scoringManager->GetMesh(iMesh);
349 if (mesh && mesh->IsActive()) {
350 MeshScoreMap scoreMap = mesh->GetScoreMap();
351 const G4String& mapNam = const_cast<G4THitsMap<G4double>&>(hits).GetName();
352 for(MeshScoreMap::const_iterator i = scoreMap.begin();
353 i != scoreMap.end(); ++i) {
354 const G4String& scoreMapName = i->first;
355 if (scoreMapName == mapNam) {
356 G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
357 scoreMapHits = true;
358 mesh->DrawMesh(scoreMapName, &colorMap);
359 }
360 }
361 }
362 }
363 }
364 if (scoreMapHits) {
365 static G4bool first = true;
366 if (first) {
367 first = false;
368 G4cout <<
369 "Scoring map drawn with default parameters."
370 "\n To get gMocren file for gMocren browser:"
371 "\n /vis/open gMocrenFile"
372 "\n /vis/viewer/flush"
373 "\n Many other options available with /score/draw... commands."
374 "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
375 << G4endl;
376 }
377 } else { // Not score map hits. Just call DrawAllHits.
378 // Cast away const because DrawAllHits is non-const!!!!
379 const_cast<G4THitsMap<G4double>&>(hits).DrawAllHits();
380 }
381}
382
384 using MeshScoreMap = G4VScoringMesh::MeshScoreMap;
385 //G4cout << "AddCompound: hits: " << &hits << G4endl;
386 G4bool scoreMapHits = false;
388 if (scoringManager) {
389 size_t nMeshes = scoringManager->GetNumberOfMesh();
390 for (size_t iMesh = 0; iMesh < nMeshes; ++iMesh) {
391 G4VScoringMesh* mesh = scoringManager->GetMesh(iMesh);
392 if (mesh && mesh->IsActive()) {
393 MeshScoreMap scoreMap = mesh->GetScoreMap();
394 for(MeshScoreMap::const_iterator i = scoreMap.begin();
395 i != scoreMap.end(); ++i) {
396 const G4String& scoreMapName = i->first;
397 const G4THitsMap<G4StatDouble>* foundHits = i->second;
398 if (foundHits == &hits) {
399 G4DefaultLinearColorMap colorMap("G4VSceneHandlerColorMap");
400 scoreMapHits = true;
401 mesh->DrawMesh(scoreMapName, &colorMap);
402 }
403 }
404 }
405 }
406 }
407 if (scoreMapHits) {
408 static G4bool first = true;
409 if (first) {
410 first = false;
411 G4cout <<
412 "Scoring map drawn with default parameters."
413 "\n To get gMocren file for gMocren browser:"
414 "\n /vis/open gMocrenFile"
415 "\n /vis/viewer/flush"
416 "\n Many other options available with /score/draw... commands."
417 "\n You might want to \"/vis/viewer/set/autoRefresh false\"."
418 << G4endl;
419 }
420 } else { // Not score map hits. Just call DrawAllHits.
421 // Cast away const because DrawAllHits is non-const!!!!
422 const_cast<G4THitsMap<G4StatDouble>&>(hits).DrawAllHits();
423 }
424}
425
427{
429 ed << "There has been an attempt to draw a mesh (a nested parameterisation),"
430 "\nbut it is not implemented by the current graphics driver. Here we simply"
431 "\ndraw the container, \"" << mesh.GetContainerVolume()->GetName() << "\".";
432 G4Exception("G4VSceneHandler::AddCompound(const G4Mesh&)",
433 "visman0107", JustWarning, ed);
434
435 const auto& pv = mesh.GetContainerVolume();
436 const auto& lv = pv->GetLogicalVolume();
437 const auto& solid = lv->GetSolid();
438 const auto& transform = mesh.GetTransform();
439 // Make sure container is visible
440 const auto& saveVisAtts = lv->GetVisAttributes();
441 auto tmpVisAtts = *saveVisAtts;
442 tmpVisAtts.SetVisibility(true);
443 auto colour = saveVisAtts->GetColour();
444 colour.SetAlpha(1.);
445 tmpVisAtts.SetColour(colour);
446 // Draw container
447 PreAddSolid(transform,tmpVisAtts);
448 solid->DescribeYourselfTo(*this);
449 PostAddSolid();
450 // Restore vis attributes
451 lv->SetVisAttributes(saveVisAtts);
452}
453
455 fViewerList.push_back (pViewer);
456}
457
459 switch (polymarker.GetMarkerType()) {
460 default:
462 {
463 G4Circle dot (polymarker);
464 dot.SetWorldSize (0.);
465 dot.SetScreenSize (0.1); // Very small circle.
466 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
467 dot.SetPosition (polymarker[iPoint]);
468 AddPrimitive (dot);
469 }
470 }
471 break;
473 {
474 G4Circle circle (polymarker); // Default circle
475 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
476 circle.SetPosition (polymarker[iPoint]);
477 AddPrimitive (circle);
478 }
479 }
480 break;
482 {
483 G4Square square (polymarker); // Default square
484 for (size_t iPoint = 0; iPoint < polymarker.size (); iPoint++) {
485 square.SetPosition (polymarker[iPoint]);
486 AddPrimitive (square);
487 }
488 }
489 break;
490 }
491}
492
494 fViewerList.remove(pViewer); // Does nothing if already removed
495 // And reset current viewer
496 auto visManager = G4VisManager::GetInstance();
497 visManager->SetCurrentViewer(nullptr);
498}
499
500
502 G4cerr << "WARNING: Plotter not implemented for " << fSystem.GetName() << G4endl;
503 G4cerr << " Open a plotter-aware graphics system or remove plotter with" << G4endl;
504 G4cerr << " /vis/scene/removeModel Plotter" << G4endl;
505}
506
508 fpScene = pScene;
509 // Notify all viewers that a kernel visit is required.
511 for (i = fViewerList.begin(); i != fViewerList.end(); i++) {
512 (*i) -> SetNeedKernelVisit (true);
513 }
514}
515
517{
520
521 switch (style) {
522 default:
527 {
528 // Use polyhedral representation
530 G4Polyhedron* pPolyhedron = solid.GetPolyhedron ();
532 if (pPolyhedron) {
533 pPolyhedron -> SetVisAttributes (fpVisAttribs);
535 AddPrimitive (*pPolyhedron);
536 EndPrimitives ();
537 break;
538 } else { // Print warnings and drop through to cloud
540 static std::set<const G4VSolid*> problematicSolids;
541 if (verbosity >= G4VisManager::errors &&
542 problematicSolids.find(&solid) == problematicSolids.end()) {
543 problematicSolids.insert(&solid);
544 G4cerr <<
545 "ERROR: G4VSceneHandler::RequestPrimitives"
546 "\n Polyhedron not available for " << solid.GetName ();
547 G4PhysicalVolumeModel* pPVModel = dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
548 if (pPVModel) {
549 G4cerr << "\n Touchable path: " << pPVModel->GetFullPVPath();
550 }
551 static G4bool explanation = false;
552 if (!explanation) {
553 explanation = true;
554 G4cerr <<
555 "\n This means it cannot be visualized in the usual way on most systems."
556 "\n 1) The solid may not have implemented the CreatePolyhedron method."
557 "\n 2) For Boolean solids, the BooleanProcessor, which attempts to create"
558 "\n the resultant polyhedron, may have failed."
559 "\n Try RayTracer. It uses Geant4's tracking algorithms instead.";
560 }
561 G4cerr << "\n Drawing solid with cloud of points.";
562 G4cerr << G4endl;
563 }
564 }
565 }
566 [[fallthrough]];
567
569 {
570 // Form solid out of cloud of dots on surface of solid
571 G4Polymarker dots;
572 // Note: OpenGL has a fast implementation of polymarker so it's better
573 // to build a polymarker rather than add a succession of circles.
574 // And anyway, in Qt, in the latter case each circle would be a scene-tree
575 // entry, something we would want to avoid.
578 dots.SetSize(G4VMarker::screen,1.);
579 G4int numberOfCloudPoints = GetNumberOfCloudPoints(fpVisAttribs);
580 if (numberOfCloudPoints <= 0) numberOfCloudPoints = vp.GetNumberOfCloudPoints();
581 for (G4int i = 0; i < numberOfCloudPoints; ++i) {
583 dots.push_back(p);
584 }
586 AddPrimitive(dots);
587 EndPrimitives ();
588 break;
589 }
590 }
591}
592
593//namespace {
594// void DrawExtent(const G4VModel* pModel)
595// {
596// // Show extent boxes - debug only, OGLSX only (OGLSQt problem?)
597// if (pModel->GetExtent() != G4VisExtent::GetNullExtent()) {
598// const auto& extent = pModel->GetExtent();
599// const auto& centre = extent.GetExtentCenter();
600// const auto& position = G4Translate3D(centre);
601// const auto& dx = (extent.GetXmax()-extent.GetXmin())/2.;
602// const auto& dy = (extent.GetYmax()-extent.GetYmin())/2.;
603// const auto& dz = (extent.GetZmax()-extent.GetZmin())/2.;
604// auto visAtts = G4VisAttributes();
605// visAtts.SetForceWireframe();
606// G4Box extentBox("Extent",dx,dy,dz);
607// G4VisManager::GetInstance()->Draw(extentBox,visAtts,position);
608// }
609// }
610//}
611
613{
614 // Assumes graphics database store has already been cleared if
615 // relevant for the particular scene handler.
616
617 if(!fpScene)
618 return;
619
621 {
622 G4Exception("G4VSceneHandler::ProcessScene", "visman0106", JustWarning,
623 "The scene has no extent.");
624 }
625
627
628 if(!visManager->GetConcreteInstance())
629 return;
630
631 G4VisManager::Verbosity verbosity = visManager->GetVerbosity();
632
633 fReadyForTransients = false;
634
635 // Reset fMarkForClearingTransientStore. (Leaving
636 // fMarkForClearingTransientStore true causes problems with
637 // recomputing transients below.) Restore it again at end...
638 G4bool tmpMarkForClearingTransientStore = fMarkForClearingTransientStore;
640
641 // Traverse geometry tree and send drawing primitives to window(s).
642
643 const std::vector<G4Scene::Model>& runDurationModelList =
645
646 if(runDurationModelList.size())
647 {
648 if(verbosity >= G4VisManager::confirmations)
649 {
650 G4cout << "Traversing scene data..." << G4endl;
651 }
652
654
655 // Create modeling parameters from view parameters...
657
658 for(size_t i = 0; i < runDurationModelList.size(); i++)
659 {
660 if(runDurationModelList[i].fActive)
661 {
662 fpModel = runDurationModelList[i].fpModel;
665 // To see the extents of each model represented as wireframe boxes,
666 // uncomment the next line and DrawExtent in namespace above
667 // DrawExtent(fpModel);
669 }
670 }
671
672 fpModel = 0;
673 delete pMP;
674
675 EndModeling();
676 }
677
678 fReadyForTransients = true;
679
680 // Refresh event from end-of-event model list.
681 // Allow only in Idle or GeomClosed state...
683 G4ApplicationState state = stateManager->GetCurrentState();
684 if(state == G4State_Idle || state == G4State_GeomClosed)
685 {
686 visManager->SetEventRefreshing(true);
687
688 if(visManager->GetRequestedEvent())
689 {
690 DrawEvent(visManager->GetRequestedEvent());
691 }
692 else
693 {
695 if(runManager)
696 {
697 const G4Run* run = runManager->GetCurrentRun();
698 const std::vector<const G4Event*>* events =
699 run ? run->GetEventVector() : 0;
700 size_t nKeptEvents = 0;
701 if(events)
702 nKeptEvents = events->size();
703 if(nKeptEvents)
704 {
706 {
707 if(verbosity >= G4VisManager::confirmations)
708 {
709 G4cout << "Refreshing event..." << G4endl;
710 }
711 const G4Event* event = 0;
712 if(events && events->size())
713 event = events->back();
714 if(event)
715 DrawEvent(event);
716 }
717 else
718 { // Accumulating events.
719
720 if(verbosity >= G4VisManager::confirmations)
721 {
722 G4cout << "Refreshing events in run..." << G4endl;
723 }
724 for(const auto& event : *events)
725 {
726 if(event)
727 DrawEvent(event);
728 }
729
731 {
732 if(verbosity >= G4VisManager::warnings)
733 {
734 G4cout << "WARNING: Cannot refresh events accumulated over more"
735 "\n than one runs. Refreshed just the last run."
736 << G4endl;
737 }
738 }
739 }
740 }
741 }
742 }
743 visManager->SetEventRefreshing(false);
744 }
745
746 // Refresh end-of-run model list.
747 // Allow only in Idle or GeomClosed state...
748 if(state == G4State_Idle || state == G4State_GeomClosed)
749 {
751 }
752
753 fMarkForClearingTransientStore = tmpMarkForClearingTransientStore;
754}
755
757{
758 const std::vector<G4Scene::Model>& EOEModelList =
759 fpScene -> GetEndOfEventModelList ();
760 size_t nModels = EOEModelList.size();
761 if (nModels) {
763 pMP->SetEvent(event);
764 for (size_t i = 0; i < nModels; i++) {
765 if (EOEModelList[i].fActive) {
766 fpModel = EOEModelList[i].fpModel;
767 fpModel -> SetModelingParameters(pMP);
768 fpModel -> DescribeYourselfTo (*this);
769 fpModel -> SetModelingParameters(0);
770 }
771 }
772 fpModel = 0;
773 delete pMP;
774 }
775}
776
778{
779 const std::vector<G4Scene::Model>& EORModelList =
780 fpScene -> GetEndOfRunModelList ();
781 size_t nModels = EORModelList.size();
782 if (nModels) {
784 pMP->SetEvent(0);
785 for (size_t i = 0; i < nModels; i++) {
786 if (EORModelList[i].fActive) {
787 fpModel = EORModelList[i].fpModel;
788 fpModel -> SetModelingParameters(pMP);
789 fpModel -> DescribeYourselfTo (*this);
790 fpModel -> SetModelingParameters(0);
791 }
792 }
793 fpModel = 0;
794 delete pMP;
795 }
796}
797
799{
800 // Create modeling parameters from View Parameters...
801 if (!fpViewer) return NULL;
802
803 const G4ViewParameters& vp = fpViewer -> GetViewParameters ();
804
805 // Convert drawing styles...
806 G4ModelingParameters::DrawingStyle modelDrawingStyle =
808 switch (vp.GetDrawingStyle ()) {
809 default:
811 modelDrawingStyle = G4ModelingParameters::wf;
812 break;
814 modelDrawingStyle = G4ModelingParameters::hlr;
815 break;
817 modelDrawingStyle = G4ModelingParameters::hsr;
818 break;
820 modelDrawingStyle = G4ModelingParameters::hlhsr;
821 break;
823 modelDrawingStyle = G4ModelingParameters::cloud;
824 break;
825 }
826
827 // Decide if covered daughters are really to be culled...
828 G4bool reallyCullCovered =
829 vp.IsCullingCovered() // Culling daughters depends also on...
830 && !vp.IsSection () // Sections (DCUT) not requested.
831 && !vp.IsCutaway () // Cutaways not requested.
832 ;
833
834 G4ModelingParameters* pModelingParams = new G4ModelingParameters
836 modelDrawingStyle,
837 vp.IsCulling (),
838 vp.IsCullingInvisible (),
839 vp.IsDensityCulling (),
840 vp.GetVisibleDensity (),
841 reallyCullCovered,
842 vp.GetNoOfSides ()
843 );
844
845 pModelingParams->SetNumberOfCloudPoints(vp.GetNumberOfCloudPoints());
846 pModelingParams->SetWarning
848
849 pModelingParams->SetCBDAlgorithmNumber(vp.GetCBDAlgorithmNumber());
850 pModelingParams->SetCBDParameters(vp.GetCBDParameters());
851
852 pModelingParams->SetExplodeFactor(vp.GetExplodeFactor());
853 pModelingParams->SetExplodeCentre(vp.GetExplodeCentre());
854
855 pModelingParams->SetSectionSolid(CreateSectionSolid());
856 pModelingParams->SetCutawaySolid(CreateCutawaySolid());
857 // The polyhedron objects are deleted in the modeling parameters destructor.
858
860
861 pModelingParams->SetSpecialMeshRendering(vp.IsSpecialMeshRendering());
862 pModelingParams->SetSpecialMeshVolumes(vp.GetSpecialMeshVolumes());
863
864 return pModelingParams;
865}
866
868{
869 G4DisplacedSolid* sectioner = 0;
870
872 if (vp.IsSection () ) {
873
875 G4double safe = radius + fpScene->GetExtent().GetExtentCentre().mag();
876 G4VSolid* sectionBox =
877 new G4Box("_sectioner", safe, safe, 1.e-5 * radius); // Thin in z-plane...
878 const G4Normal3D originalNormal(0,0,1); // ...so this is original normal.
879
880 const G4Plane3D& sp = vp.GetSectionPlane ();
881 const G4double& a = sp.a();
882 const G4double& b = sp.b();
883 const G4double& c = sp.c();
884 const G4double& d = sp.d();
885 const G4Normal3D newNormal(a,b,c);
886
887 G4Transform3D requiredTransform;
888 // Rotate
889 if (newNormal != originalNormal) {
890 const G4double& angle = std::acos(newNormal.dot(originalNormal));
891 const G4Vector3D& axis = originalNormal.cross(newNormal);
892 requiredTransform = G4Rotate3D(angle, axis);
893 }
894 // Translate
895 requiredTransform = requiredTransform * G4TranslateZ3D(-d);
896
897 sectioner = new G4DisplacedSolid
898 ("_displaced_sectioning_box", sectionBox, requiredTransform);
899 }
900
901 return sectioner;
902}
903
905{
906 // To be reviewed.
907 return 0;
908 /*** An alternative way of getting a cutaway is to use
909 Command /vis/scene/add/volume
910 Guidance :
911 Adds a physical volume to current scene, with optional clipping volume.
912 If physical-volume-name is "world" (the default), the top of the
913 main geometry tree (material world) is added. If "worlds", the
914 top of all worlds - material world and parallel worlds, if any - are
915 added. Otherwise a search of all worlds is made, taking the first
916 matching occurrence only. To see a representation of the geometry
917 hierarchy of the worlds, try "/vis/drawTree [worlds]" or one of the
918 driver/browser combinations that have the required functionality, e.g., HepRep.
919 If clip-volume-type is specified, the subsequent parameters are used to
920 to define a clipping volume. For example,
921 "/vis/scene/add/volume ! ! ! -box km 0 1 0 1 0 1" will draw the world
922 with the positive octant cut away. (If the Boolean Processor issues
923 warnings try replacing 0 by 0.000000001 or something.)
924 If clip-volume-type is prepended with '-', the clip-volume is subtracted
925 (cutaway). (This is the default if there is no prepended character.)
926 If '*' is prepended, the intersection of the physical-volume and the
927 clip-volume is made. (You can make a section/DCUT with a thin box, for
928 example).
929 For "box", the parameters are xmin,xmax,ymin,ymax,zmin,zmax.
930 Only "box" is programmed at present.
931 ***/
932}
933
935{
936 // Load G4Atts from G4VisAttributes, if any...
937 const G4VisAttributes* va = visible.GetVisAttributes();
938 if (va) {
939 const std::map<G4String,G4AttDef>* vaDefs =
940 va->GetAttDefs();
941 if (vaDefs) {
942 holder->AddAtts(visible.GetVisAttributes()->CreateAttValues(), vaDefs);
943 }
944 }
945
946 G4PhysicalVolumeModel* pPVModel =
947 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
948 if (pPVModel) {
949 // Load G4Atts from G4PhysicalVolumeModel...
950 const std::map<G4String,G4AttDef>* pvDefs = pPVModel->GetAttDefs();
951 if (pvDefs) {
952 holder->AddAtts(pPVModel->CreateCurrentAttValues(), pvDefs);
953 }
954 }
955
956 G4TrajectoriesModel* trajModel = dynamic_cast<G4TrajectoriesModel*>(fpModel);
957 if (trajModel) {
958 // Load G4Atts from trajectory model...
959 const std::map<G4String,G4AttDef>* trajModelDefs = trajModel->GetAttDefs();
960 if (trajModelDefs) {
961 holder->AddAtts(trajModel->CreateCurrentAttValues(), trajModelDefs);
962 }
963 // Load G4Atts from trajectory...
964 const G4VTrajectory* traj = trajModel->GetCurrentTrajectory();
965 if (traj) {
966 const std::map<G4String,G4AttDef>* trajDefs = traj->GetAttDefs();
967 if (trajDefs) {
968 holder->AddAtts(traj->CreateAttValues(), trajDefs);
969 }
970 G4int nPoints = traj->GetPointEntries();
971 for (G4int i = 0; i < nPoints; ++i) {
972 G4VTrajectoryPoint* trajPoint = traj->GetPoint(i);
973 if (trajPoint) {
974 const std::map<G4String,G4AttDef>* pointDefs = trajPoint->GetAttDefs();
975 if (pointDefs) {
976 holder->AddAtts(trajPoint->CreateAttValues(), pointDefs);
977 }
978 }
979 }
980 }
981 }
982
983 G4HitsModel* hitsModel = dynamic_cast<G4HitsModel*>(fpModel);
984 if (hitsModel) {
985 // Load G4Atts from hit...
986 const G4VHit* hit = hitsModel->GetCurrentHit();
987 const std::map<G4String,G4AttDef>* hitsDefs = hit->GetAttDefs();
988 if (hitsDefs) {
989 holder->AddAtts(hit->CreateAttValues(), hitsDefs);
990 }
991 }
992}
993
996 const G4Colour& colour = fpVisAttribs -> GetColour ();
997 return colour;
998}
999
1001 auto pVA = visible.GetVisAttributes();
1003 return pVA->GetColour();
1004}
1005
1007 auto pVA = text.GetVisAttributes();
1009 return pVA->GetColour();
1010}
1011
1013{
1014 G4double lineWidth = pVisAttribs->GetLineWidth();
1015 if (lineWidth < 1.) lineWidth = 1.;
1016 lineWidth *= fpViewer -> GetViewParameters().GetGlobalLineWidthScale();
1017 if (lineWidth < 1.) lineWidth = 1.;
1018 return lineWidth;
1019}
1020
1022(const G4VisAttributes* pVisAttribs) {
1023 // Drawing style is normally determined by the view parameters, but
1024 // it can be overriddden by the ForceDrawingStyle flag in the vis
1025 // attributes.
1027 const G4ViewParameters::DrawingStyle viewerStyle = vp.GetDrawingStyle();
1028 G4ViewParameters::DrawingStyle resultantStyle = viewerStyle;
1029 if (pVisAttribs -> IsForceDrawingStyle ()) {
1031 pVisAttribs -> GetForcedDrawingStyle ();
1032 // This is complicated because if hidden line and surface removal
1033 // has been requested we wish to preserve this sometimes.
1034 switch (forcedStyle) {
1036 switch (viewerStyle) {
1037 case (G4ViewParameters::hlr):
1038 resultantStyle = G4ViewParameters::hlhsr;
1039 break;
1041 resultantStyle = G4ViewParameters::hsr;
1042 break;
1044 resultantStyle = G4ViewParameters::hsr;
1045 break;
1047 case (G4ViewParameters::hsr):
1048 break;
1049 }
1050 break;
1052 resultantStyle = G4ViewParameters::cloud;
1053 break;
1055 default:
1056 // But if forced style is wireframe, do it, because one of its
1057 // main uses is in displaying the consituent solids of a Boolean
1058 // solid and their surfaces overlap with the resulting Booean
1059 // solid, making a mess if hlr is specified.
1060 resultantStyle = G4ViewParameters::wireframe;
1061 break;
1062 }
1063 }
1064 return resultantStyle;
1065}
1066
1068(const G4VisAttributes* pVisAttribs) const {
1069 // Returns no of cloud points from current view parameters, unless the user
1070 // has forced through the vis attributes, thereby over-riding the
1071 // current view parameter.
1072 G4int numberOfCloudPoints = fpViewer->GetViewParameters().GetNumberOfCloudPoints();
1073 if (pVisAttribs -> IsForceDrawingStyle() &&
1074 pVisAttribs -> GetForcedDrawingStyle() == G4VisAttributes::cloud &&
1075 pVisAttribs -> GetForcedNumberOfCloudPoints() > 0) {
1076 numberOfCloudPoints = pVisAttribs -> GetForcedNumberOfCloudPoints();
1077 }
1078 return numberOfCloudPoints;
1079}
1080
1082 G4bool isAuxEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible ();
1083 if (pVisAttribs -> IsForceAuxEdgeVisible()) {
1084 isAuxEdgeVisible = pVisAttribs->IsForcedAuxEdgeVisible();
1085 }
1086 return isAuxEdgeVisible;
1087}
1088
1090(const G4VMarker& marker,
1091 G4VSceneHandler::MarkerSizeType& markerSizeType)
1092{
1093 G4bool userSpecified = marker.GetWorldSize() || marker.GetScreenSize();
1094 const G4VMarker& defaultMarker =
1095 fpViewer -> GetViewParameters().GetDefaultMarker();
1096 G4double size = userSpecified ?
1097 marker.GetWorldSize() : defaultMarker.GetWorldSize();
1098 if (size) {
1099 // Draw in world coordinates.
1100 markerSizeType = world;
1101 }
1102 else {
1103 size = userSpecified ?
1104 marker.GetScreenSize() : defaultMarker.GetScreenSize();
1105 // Draw in screen coordinates.
1106 markerSizeType = screen;
1107 }
1108 size *= fpViewer -> GetViewParameters().GetGlobalMarkerScale();
1109 if (markerSizeType == screen && size < 1.) size = 1.;
1110 return size;
1111}
1112
1114{
1115 // No. of sides (lines segments per circle) is normally determined
1116 // by the view parameters, but it can be overriddden by the
1117 // ForceLineSegmentsPerCircle in the vis attributes.
1118 G4int lineSegmentsPerCircle = fpViewer->GetViewParameters().GetNoOfSides();
1119 if (pVisAttribs) {
1120 if (pVisAttribs->IsForceLineSegmentsPerCircle())
1121 lineSegmentsPerCircle = pVisAttribs->GetForcedLineSegmentsPerCircle();
1122 if (lineSegmentsPerCircle < pVisAttribs->GetMinLineSegmentsPerCircle()) {
1123 lineSegmentsPerCircle = pVisAttribs->GetMinLineSegmentsPerCircle();
1124 G4cout <<
1125 "G4VSceneHandler::GetNoOfSides: attempt to set the"
1126 "\nnumber of line segments per circle < " << lineSegmentsPerCircle
1127 << "; forced to " << pVisAttribs->GetMinLineSegmentsPerCircle() << G4endl;
1128 }
1129 }
1130 return lineSegmentsPerCircle;
1131}
1132
1133std::ostream& operator << (std::ostream& os, const G4VSceneHandler& sh) {
1134
1135 os << "Scene handler " << sh.fName << " has "
1136 << sh.fViewerList.size () << " viewer(s):";
1137 for (size_t i = 0; i < sh.fViewerList.size (); i++) {
1138 os << "\n " << *(sh.fViewerList [i]);
1139 }
1140
1141 if (sh.fpScene) {
1142 os << "\n " << *sh.fpScene;
1143 }
1144 else {
1145 os << "\n This scene handler currently has no scene.";
1146 }
1147
1148 return os;
1149}
G4ApplicationState
@ G4State_Idle
@ G4State_GeomClosed
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static const G4double angle[DIMMOTT]
HepGeom::Rotate3D G4Rotate3D
HepGeom::TranslateZ3D G4TranslateZ3D
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
std::vector< G4VViewer * >::iterator G4ViewerListIterator
Definition: G4ViewerList.hh:42
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void AddAtts(const std::vector< G4AttValue > *values, const std::map< G4String, G4AttDef > *defs)
Definition: G4AttHolder.hh:64
Definition: G4Box.hh:56
Definition: G4Cons.hh:78
const G4VHit * GetCurrentHit() const
Definition: G4HitsModel.hh:57
G4VSolid * GetSolid() const
Definition: G4Mesh.hh:47
G4VPhysicalVolume * GetContainerVolume() const
Definition: G4Mesh.hh:61
const G4Transform3D & GetTransform() const
Definition: G4Mesh.hh:64
void SetCBDParameters(const std::vector< G4double > &)
void SetWarning(G4bool)
void SetNumberOfCloudPoints(G4int)
void SetCBDAlgorithmNumber(G4int)
void SetExplodeFactor(G4double explodeFactor)
void SetVisAttributesModifiers(const std::vector< VisAttributesModifier > &)
void SetExplodeCentre(const G4Point3D &explodeCentre)
void SetCutawaySolid(G4DisplacedSolid *pCutawaySolid)
void SetSectionSolid(G4DisplacedSolid *pSectionSolid)
void SetEvent(const G4Event *pEvent)
void SetSpecialMeshVolumes(const std::vector< PVNameCopyNo > &)
void SetSpecialMeshRendering(G4bool)
Definition: G4Orb.hh:56
Definition: G4Para.hh:79
std::vector< G4AttValue > * CreateCurrentAttValues() const
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
void SetMarkerType(MarkerType)
MarkerType GetMarkerType() const
static G4RunManager * GetMasterRunManager()
const G4Run * GetCurrentRun() const
Definition: G4Run.hh:49
const std::vector< Model > & GetRunDurationModelList() const
G4bool GetRefreshAtEndOfEvent() const
const G4VisExtent & GetExtent() const
G4bool GetRefreshAtEndOfRun() const
G4VScoringMesh * GetMesh(G4int i) const
size_t GetNumberOfMesh() const
static G4ScoringManager * GetScoringManagerIfExist()
const G4ApplicationState & GetCurrentState() const
static G4StateManager * GetStateManager()
Definition: G4Text.hh:72
const G4VTrajectory * GetCurrentTrajectory() const
std::vector< G4AttValue > * CreateCurrentAttValues() const
const std::map< G4String, G4AttDef > * GetAttDefs() const
Definition: G4Trd.hh:63
Definition: G4Tubs.hh:75
const G4String & GetName() const
Definition: G4VHit.hh:48
virtual std::vector< G4AttValue > * CreateAttValues() const
Definition: G4VHit.hh:64
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
Definition: G4VHit.hh:58
G4double GetScreenSize() const
void SetSize(SizeType, G4double)
Definition: G4VMarker.cc:94
void SetScreenSize(G4double)
void SetWorldSize(G4double)
void SetPosition(const G4Point3D &)
G4double GetWorldSize() const
void SetModelingParameters(const G4ModelingParameters *)
virtual void DescribeYourselfTo(G4VGraphicsScene &)=0
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
virtual void BeginModeling()
G4int GetNumberOfCloudPoints(const G4VisAttributes *) const
G4int GetNoOfSides(const G4VisAttributes *)
virtual void ClearTransientStore()
void LoadAtts(const G4Visible &, G4AttHolder *)
void DrawEvent(const G4Event *)
G4ModelingParameters * CreateModelingParameters()
const G4Colour & GetTextColour(const G4Text &)
const G4Colour & GetColour()
void AddSolidT(const T &solid)
G4Transform3D fObjectTransformation
virtual void EndPrimitives()
G4bool fTransientsDrawnThisEvent
void AddSolidWithAuxiliaryEdges(const T &solid)
virtual G4DisplacedSolid * CreateSectionSolid()
virtual void EndModeling()
virtual const G4VisExtent & GetExtent() const
G4ViewerList fViewerList
virtual void ProcessScene()
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
G4VSceneHandler(G4VGraphicsSystem &system, G4int id, const G4String &name="")
G4bool fTransientsDrawnThisRun
G4VViewer * fpViewer
virtual void PostAddSolid()
const G4String & GetName() const
void AddViewerToList(G4VViewer *pView)
virtual void EndPrimitives2D()
virtual void SetScene(G4Scene *)
G4bool fMarkForClearingTransientStore
const G4VisAttributes * fpVisAttribs
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation=G4Transform3D())
virtual void RequestPrimitives(const G4VSolid &solid)
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
void RemoveViewerFromList(G4VViewer *pView)
virtual G4DisplacedSolid * CreateCutawaySolid()
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())
G4double GetLineWidth(const G4VisAttributes *)
G4VGraphicsSystem & fSystem
virtual void AddSolid(const G4Box &)
virtual void ClearStore()
virtual void AddCompound(const G4VTrajectory &)
virtual ~G4VSceneHandler()
virtual void AddPrimitive(const G4Polyline &)=0
G4bool GetAuxEdgeVisible(const G4VisAttributes *)
G4bool IsActive() const
std::map< G4String, RunScore * > MeshScoreMap
void DrawMesh(const G4String &psName, G4VScoreColorMap *colorMap, G4int axflg=111)
MeshScoreMap GetScoreMap() const
G4String GetName() const
virtual G4ThreeVector GetPointOnSurface() const
Definition: G4VSolid.cc:152
virtual G4Polyhedron * GetPolyhedron() const
Definition: G4VSolid.cc:705
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual void DrawTrajectory() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4ViewParameters & GetViewParameters() const
static G4VVisManager * GetConcreteInstance()
const std::vector< G4ModelingParameters::VisAttributesModifier > & GetVisAttributesModifiers() const
G4int GetNoOfSides() const
G4bool IsSpecialMeshRendering() const
G4double GetExplodeFactor() const
G4int GetNumberOfCloudPoints() const
G4bool IsCutaway() const
G4bool IsSection() const
G4bool IsCulling() const
const G4VisAttributes * GetDefaultTextVisAttributes() const
const std::vector< G4double > & GetCBDParameters() const
G4int GetCBDAlgorithmNumber() const
const std::vector< G4ModelingParameters::PVNameCopyNo > & GetSpecialMeshVolumes() const
G4bool IsCullingInvisible() const
const G4VisAttributes * GetDefaultVisAttributes() const
G4bool IsDensityCulling() const
G4double GetVisibleDensity() const
const G4Point3D & GetExplodeCentre() const
G4bool IsCullingCovered() const
const G4Plane3D & GetSectionPlane() const
DrawingStyle GetDrawingStyle() const
G4bool IsAuxEdgeVisible() const
void remove(G4VViewer *)
Definition: G4ViewerList.cc:30
const std::map< G4String, G4AttDef > * GetAttDefs() const
G4bool IsForceLineSegmentsPerCircle() const
G4double GetLineWidth() const
void SetForceAuxEdgeVisible(G4bool=true)
G4int GetForcedLineSegmentsPerCircle() const
const std::vector< G4AttValue > * CreateAttValues() const
const G4Colour & GetColour() const
G4bool IsForceAuxEdgeVisible() const
G4bool IsForcedAuxEdgeVisible() const
static G4int GetMinLineSegmentsPerCircle()
static const G4VisExtent & GetNullExtent()
Definition: G4VisExtent.cc:60
G4double GetExtentRadius() const
Definition: G4VisExtent.cc:75
const G4Point3D & GetExtentCentre() const
Definition: G4VisExtent.cc:65
void SetEventRefreshing(G4bool)
G4bool GetTransientsDrawnThisEvent() const
G4bool GetTransientsDrawnThisRun() const
static Verbosity GetVerbosity()
const G4Event * GetRequestedEvent() const
static G4VisManager * GetInstance()
void SetVisAttributes(const G4VisAttributes *)
Definition: G4Visible.cc:96
const G4VisAttributes * GetVisAttributes() const
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
T dot(const BasicVector3D< T > &v) const
static void SetNumberOfRotationSteps(G4int n)
static void ResetNumberOfRotationSteps()
const char * name(G4int ptype)
G4bool transform(G4String &input, const G4String &type)
Definition: run.py:1