Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Qt3DSceneHandler.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// John Allison 17th June 2019
27
28#include "G4Qt3DSceneHandler.hh"
29
32#include "G4VPhysicalVolume.hh"
33#include "G4LogicalVolume.hh"
35#include "G4Box.hh"
36#include "G4Polyline.hh"
37#include "G4Polymarker.hh"
38#include "G4Text.hh"
39#include "G4Circle.hh"
40#include "G4Square.hh"
41#include "G4Polyhedron.hh"
42#include "G4Scene.hh"
43#include "G4Threading.hh"
44#include "G4Mesh.hh"
45#include "G4PseudoScene.hh"
46#include "G4VisManager.hh"
47
48#include "G4Qt3DViewer.hh"
49#include "G4Qt3DUtils.hh"
50#include "G4Qt3DQEntity.hh"
51
52#include <Qt3DCore>
53#include <Qt3DExtras>
54#include <Qt3DRender>
55
56#include <utility>
57
58#define G4warn G4cout
59
60// Qt3D seems to offer a choice of type - float or double. It would be nice
61// to use double since it offers the prospect of higher precision, hopefully
62// avoiding some issues that we see at high zoom. But it currently gives the
63// following warning:
64// Qt5: "findBoundingVolumeComputeData: Position attribute not
65// suited for bounding volume computation",
66// Qt6: "Failed to build graphics pipeline: Geometry doesn't match expected layout
67// An attribute type is not supported "vertexPosition" Qt3DCore::QAttribute::Double"
68// so for now we use float.
69#define PRECISION float
70#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
71#define BASETYPE Qt3DRender::QAttribute::Float
72#else
73#define BASETYPE Qt3DCore::QAttribute::Float
74#endif
75
76// Qt3D types move between namespaces between 5 and 6
77namespace G4Qt3DCompat
78{
79#if (QT_VERSION < QT_VERSION_CHECK(6, 0, 0))
80 using Qt3DRender::QAttribute;
81 using Qt3DRender::QBuffer;
82 using Qt3DRender::QGeometry;
83#else
84 using Qt3DCore::QAttribute;
85 using Qt3DCore::QBuffer;
86 using Qt3DCore::QGeometry;
87#endif
88}
89
91
93 (G4VGraphicsSystem& system, const G4String& name)
94: G4VSceneHandler(system, fSceneIdCount++, name)
95{
96#ifdef G4QT3DDEBUG
97 G4cout << "G4Qt3DSceneHandler::G4Qt3DSceneHandler called" << G4endl;
98#endif
99 fpQt3DScene = new Qt3DCore::QEntity;
100 fpQt3DScene->setObjectName("G4Qt3DSceneRoot");
102}
103
105{
106 // Doesn't like this - it gives BAD_ACCESS in delete_entity_recursively.
107 // Curiously the delete traceback shows three calls to this recursively:
108 /*#1 0x0000000100411906 in (anonymous namespace)::delete_entity_recursively(Qt3DCore::QNode*) at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:60
109 #2 0x0000000100411840 in G4Qt3DSceneHandler::~G4Qt3DSceneHandler() at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:169
110 #3 0x0000000100411fc5 in G4Qt3DSceneHandler::~G4Qt3DSceneHandler() at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:168
111 #4 0x0000000100411fe9 in G4Qt3DSceneHandler::~G4Qt3DSceneHandler() at /Users/johna/Geant4/geant4-dev/source/visualization/Qt3D/src/G4Qt3DSceneHandler.cc:168
112 #5 0x0000000101032510 in G4VisManager::~G4VisManager() at /Users/johna/Geant4/geant4-dev/source/visualization/management/src/G4VisManager.cc:214
113 #6 0x0000000100013885 in G4VisExecutive::~G4VisExecutive() at /Users/johna/Geant4/geant4-dev/source/visualization/management/include/G4VisExecutive.hh:119
114 #7 0x00000001000119a5 in G4VisExecutive::~G4VisExecutive() at /Users/johna/Geant4/geant4-dev/source/visualization/management/include/G4VisExecutive.hh:119
115 #8 0x00000001000119c9 in G4VisExecutive::~G4VisExecutive() at /Users/johna/Geant4/geant4-dev/source/visualization/management/include/G4VisExecutive.hh:119
116 #9 0x00000001000117dd in main at /Users/johna/Geant4/geant4-dev/examples/basic/B1/exampleB1.cc:108
117 */
118 //if (fpQt3DScene) delete_entity_recursively(fpQt3DScene);
119}
120
122{
123 fpTransientObjects = new G4Qt3DQEntity(fpQt3DScene); // Hangs from root
124 fpTransientObjects ->setObjectName("G4Qt3DTORoot");
125 fpPersistentObjects = new G4Qt3DQEntity(fpQt3DScene); // Hangs from root
126 fpPersistentObjects ->setObjectName("G4Qt3DPORoot");
127
128 // Physical volume objects hang from POs
129 if (fpScene) {
130 const auto& sceneModels = fpScene->GetRunDurationModelList();
131 for (const auto& sceneModel : sceneModels) {
132 const auto& pvModel = dynamic_cast<G4PhysicalVolumeModel*>(sceneModel.fpModel);
133 if (pvModel) {
134 auto entity = new G4Qt3DQEntity(fpPersistentObjects);
135 const auto& pv = pvModel->GetTopPhysicalVolume();
136 entity->setObjectName("G4Qt3DPORoot_"+QString(pv->GetName()));
137 entity->SetPVNodeID(G4PhysicalVolumeModel::G4PhysicalVolumeNodeID(pv));
138 fpPhysicalVolumeObjects.push_back(entity);
139 }
140 }
141 }
142}
143
145{
146 // Create a G4Qt3DQEntity node suitable for next solid or primitive
147
148 G4Qt3DQEntity* newNode = nullptr;
149
150 if (fReadyForTransients) { // All transients hang from this node
151 newNode = new G4Qt3DQEntity(fpTransientObjects);
152 G4String name = fpModel? fpModel->GetGlobalTag(): G4String("User");
153 newNode->setObjectName(name.c_str());
154 return newNode;
155 }
156
157 G4PhysicalVolumeModel* pPVModel =
158 dynamic_cast<G4PhysicalVolumeModel*>(fpModel);
159
160 if (!pPVModel) { // Persistent objects (e.g., axes)
161 newNode = new G4Qt3DQEntity(fpPersistentObjects);
162 newNode->setObjectName(fpModel->GetGlobalTag().c_str());
163 return newNode;
164 }
165
166 // So this is a G4PhysicalVolumeModel
167
169 typedef std::vector<PVNodeID> PVPath;
170// const PVPath& drawnPVPath = pPVModel->GetDrawnPVPath();
171 const PVPath& fullPVPath = pPVModel->GetFullPVPath();
172 //G4int currentDepth = pPVModel->GetCurrentDepth();
173 //G4VPhysicalVolume* pCurrentPV = pPVModel->GetCurrentPV();
174 //G4LogicalVolume* pCurrentLV = pPVModel->GetCurrentLV();
175 //G4Material* pCurrentMaterial = pPVModel->GetCurrentMaterial();
176 // Note: pCurrentMaterial may be zero (parallel world).
177
178#ifdef G4QTDEBUG
179 G4cout << "A: " << fullPVPath << G4endl; // DEBUG
180#endif
181
182 // Find appropriate root
183 const std::size_t nWorlds = fpPhysicalVolumeObjects.size();
184 std::size_t iWorld = 0;
185 for (; iWorld < nWorlds; ++iWorld) {
186 if (fullPVPath[0].GetPhysicalVolume() ==
187 fpPhysicalVolumeObjects[iWorld]->GetPVNodeID().GetPhysicalVolume()) break;
188 }
189 if (iWorld == nWorlds) {
190 G4Exception("G4Qt3DSceneHandler::CreateNewNode", "qt3D-0000", FatalException,
191 "World mis-match - not possible(!?)");
192 }
193
194 // (Re-)establish pv path of root entity
196 wrld->SetPVNodeID(fullPVPath[0]);
197
198 // Create nodes as required
199 G4Qt3DQEntity* node = wrld;
200 newNode = node;
201 const std::size_t depth = fullPVPath.size();
202 std::size_t iDepth = 1;
203 while (iDepth < depth) {
204 const auto& children = node->children();
205 const G4int nChildren = (G4int)children.size();
206 G4int iChild = 0;
207 G4Qt3DQEntity* child = nullptr;
208 for (; iChild < nChildren; ++iChild) {
209 child = static_cast<G4Qt3DQEntity*>(children[iChild]);
210 if (child->GetPVNodeID() == fullPVPath[iDepth]) break;
211 }
212 if (iChild != nChildren) { // Existing node found
213 node = child; // Must be the ancestor of new node (subsequent iteration)
214 } else {
215 // Add a new node as child of node
216 newNode = new G4Qt3DQEntity(node);
217 newNode->SetPVNodeID(fullPVPath[iDepth]);
218 std::ostringstream oss;
219 oss << newNode->GetPVNodeID().GetPhysicalVolume()->GetName()
220 << ':' << newNode->GetPVNodeID().GetCopyNo();
221 newNode->setObjectName(oss.str().c_str());
222 node = newNode;
223 }
224 ++iDepth;
225 }
226
227 return node;
228}
229
231 (const G4Transform3D& objectTransformation,
232 const G4VisAttributes& visAttribs)
233{
234 G4VSceneHandler::PreAddSolid(objectTransformation, visAttribs);
235}
236
241
243{
244// The x,y coordinates of the primitives passed to AddPrimitive are
245// intrepreted as screen coordinates, -1 < x,y < 1. The
246// z-coordinate is ignored.
247// IMPORTANT: invoke this from your polymorphic versions, e.g.:
248// void MyXXXSceneHandler::BeginPrimitives2D
249// (const G4Transform3D& objectTransformation) {
250 static G4bool first = true;
251 if (first) {
252 first = false;
253 G4Exception("G4Qt3DSceneHandler::BeginPrimitives2D", "qt3D-0001",
255 "2D drawing not yet implemented");
256 }
257 G4VSceneHandler::BeginPrimitives2D (objectTransformation);
258// ...
259}
260
262{
263// IMPORTANT: invoke this from your polymorphic versions, e.g.:
264// void MyXXXSceneHandler::EndPrimitives2D () {
265// ...
267}
268
270 (const G4Transform3D& objectTransformation)
271{
272 G4VSceneHandler::BeginPrimitives(objectTransformation);
273}
274
279
281{
282#ifdef G4QT3DDEBUG
283 G4cout <<
284 "G4Qt3DSceneHandler::AddPrimitive(const G4Polyline& polyline) called.\n"
285 << polyline
286 << G4endl;
287#endif
288
289 if (polyline.size() == 0) return;
290
291 auto currentNode = CreateNewNode();
292 if (!currentNode) {
293 static G4bool first = true;
294 if (first) {
295 first = false;
296 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Polyline&)",
297 "qt3d-0003", JustWarning,
298 "No available node!");
299 }
300 return;
301 }
302
304
306 transform->setObjectName("transform");
307
308 auto polylineEntity = new Qt3DCore::QEntity(currentNode);
309 polylineEntity->addComponent(transform);
310
311 const auto vertexByteSize = 3*sizeof(PRECISION);
312
313 const std::size_t nLines = polyline.size() - 1;
314 QByteArray polylineByteArray;
315 const auto polylineBufferByteSize = 2*nLines*vertexByteSize;
316 polylineByteArray.resize((G4int)polylineBufferByteSize);
317 auto polylineBufferArray = reinterpret_cast<PRECISION*>(polylineByteArray.data());
318 G4int iLine = 0;
319 for (std::size_t i = 0; i < nLines; ++i) {
320 polylineBufferArray[iLine++] = polyline[i].x();
321 polylineBufferArray[iLine++] = polyline[i].y();
322 polylineBufferArray[iLine++] = polyline[i].z();
323 polylineBufferArray[iLine++] = polyline[i+1].x();
324 polylineBufferArray[iLine++] = polyline[i+1].y();
325 polylineBufferArray[iLine++] = polyline[i+1].z();
326 }
327 auto polylineGeometry = new G4Qt3DCompat::QGeometry();
328 polylineGeometry->setObjectName("polylineGeometry");
329
330 auto polylineBuffer = new G4Qt3DCompat::QBuffer(polylineGeometry);
331 polylineBuffer->setObjectName("Polyline buffer");
332 polylineBuffer->setData(polylineByteArray);
333
334 auto polylineAtt = new G4Qt3DCompat::QAttribute;
335 polylineAtt->setObjectName("Position attribute");
336 polylineAtt->setName(G4Qt3DCompat::QAttribute::defaultPositionAttributeName());
337 polylineAtt->setBuffer(polylineBuffer);
338 polylineAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
339 polylineAtt->setVertexBaseType(BASETYPE);
340 polylineAtt->setVertexSize(3);
341 polylineAtt->setCount((G4int)nLines);
342 polylineAtt->setByteOffset(0);
343 polylineAtt->setByteStride(vertexByteSize);
344 // Normal attribute (a dummy with count==0) (Qt6 seems to require)
345 auto dummyNormalLineAtt = new G4Qt3DCompat::QAttribute;
346 dummyNormalLineAtt->setObjectName("Normal attribute");
347 dummyNormalLineAtt->setName(G4Qt3DCompat::QAttribute::defaultNormalAttributeName());
348 dummyNormalLineAtt->setBuffer(polylineBuffer);
349 dummyNormalLineAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
350 dummyNormalLineAtt->setVertexBaseType(BASETYPE);
351 dummyNormalLineAtt->setVertexSize(3);
352 dummyNormalLineAtt->setCount(0);
353 dummyNormalLineAtt->setByteOffset(0);
354 dummyNormalLineAtt->setByteStride(vertexByteSize);
355
356 const auto& colour = fpVisAttribs->GetColour();
357
358 polylineGeometry->addAttribute(polylineAtt);
359 polylineGeometry->addAttribute(dummyNormalLineAtt);
360
361 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
362 material->setObjectName("materialForPolyline");
363 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
364 material->setShininess(0.);
365 material->setSpecular(0.);
366 polylineEntity->addComponent(material);
367
368 auto renderer = new Qt3DRender::QGeometryRenderer;
369#if QT_VERSION >= QT_VERSION_CHECK(6,0,0)
370 auto geometryView = new Qt3DCore::QGeometryView(polylineGeometry);
371 geometryView->setObjectName("polylineGeometryView");
372 geometryView->setGeometry(polylineGeometry);
373 geometryView->setVertexCount((G4int)(2*nLines));
374 geometryView->setPrimitiveType(Qt3DCore::QGeometryView::Lines);
375 renderer->setView(geometryView);
376#else
377 renderer->setObjectName("polylineRenderer");
378 renderer->setGeometry(polylineGeometry);
379 renderer->setVertexCount(2*(G4int)nLines);
380 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
381#endif
382 polylineEntity->addComponent(renderer);
383}
384
386{
387 if (polymarker.size() == 0) return;
388
389 auto currentNode = CreateNewNode();
390 if (!currentNode) {
391 static G4bool first = true;
392 if (first) {
393 first = false;
394 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Polymarker&)",
395 "qt3d-0003", JustWarning,
396 "No available node!");
397 }
398 return;
399 }
400
402
403 MarkerSizeType markerSizeType;
404 G4double markerSize = GetMarkerSize(polymarker, markerSizeType);
405
406 switch (polymarker.GetMarkerType()) {
407 default:
409 {
410 const std::size_t nDots = polymarker.size();
411
413 transform->setObjectName("transform");
414
415 auto polymarkerEntity = new Qt3DCore::QEntity(currentNode);
416 polymarkerEntity->addComponent(transform);
417
418 const auto vertexByteSize = 3*sizeof(PRECISION);
419
420 QByteArray polymarkerByteArray;
421 const auto polymarkerBufferByteSize = nDots*vertexByteSize;
422 polymarkerByteArray.resize((G4int)polymarkerBufferByteSize);
423 auto polymarkerBufferArray = reinterpret_cast<PRECISION*>(polymarkerByteArray.data());
424 G4int iMarker = 0;
425 for (std::size_t i = 0; i < polymarker.size(); ++i) {
426 polymarkerBufferArray[iMarker++] = polymarker[i].x();
427 polymarkerBufferArray[iMarker++] = polymarker[i].y();
428 polymarkerBufferArray[iMarker++] = polymarker[i].z();
429 }
430 auto polymarkerGeometry = new G4Qt3DCompat::QGeometry();
431 polymarkerGeometry->setObjectName("polymarkerGeometry");
432 auto polymarkerBuffer = new G4Qt3DCompat::QBuffer(polymarkerGeometry);
433 polymarkerBuffer->setObjectName("Polymarker buffer");
434 polymarkerBuffer->setData(polymarkerByteArray);
435
436 auto polymarkerAtt = new G4Qt3DCompat::QAttribute;
437 polymarkerAtt->setObjectName("Position attribute");
438 polymarkerAtt->setName(G4Qt3DCompat::QAttribute::defaultPositionAttributeName());
439 polymarkerAtt->setBuffer(polymarkerBuffer);
440 polymarkerAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
441 polymarkerAtt->setVertexBaseType(BASETYPE);
442 polymarkerAtt->setVertexSize(3);
443 polymarkerAtt->setCount((G4int)nDots);
444 polymarkerAtt->setByteOffset(0);
445 polymarkerAtt->setByteStride(vertexByteSize);
446
447 const auto& colour = fpVisAttribs->GetColour();
448
449 polymarkerGeometry->addAttribute(polymarkerAtt);
450
451 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
452 material->setObjectName("materialForPolymarker");
453 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
454 material->setShininess(0.);
455 material->setSpecular(0.);
456 polymarkerEntity->addComponent(material);
457
458 auto renderer = new Qt3DRender::QGeometryRenderer;
459 renderer->setObjectName("polymarkerWireframeRenderer");
460 renderer->setGeometry(polymarkerGeometry);
461 renderer->setVertexCount((G4int)nDots);
462 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Points);
463 polymarkerEntity->addComponent(renderer);
464 }
465 break;
467 {
468 G4Circle circle (polymarker); // Default circle
469
470 const auto& colour = fpVisAttribs->GetColour();
471 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
472 material->setObjectName("materialForCircle");
473 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
474 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
475
476 auto sphereMesh = new Qt3DExtras::QSphereMesh;
477 sphereMesh->setObjectName("sphereMesh");
478 G4double radius = markerSize/2.;
479 if (markerSizeType == G4VSceneHandler::screen ) {
480 // Not figured out how to do screen-size, so use scene extent
481 const G4double scale = 200.; // Roughly pixels per scene
482 radius *= fpScene->GetExtent().GetExtentRadius()/scale;
483 }
484 sphereMesh->setRadius(radius);
485// sphereMesh->setInstanceCount(polymarker.size()); // Not undertood instancing yet
486
487// auto currentEntity = new Qt3DCore::QEntity(currentNode); // Not undertood instancing yet
488 for (std::size_t iPoint = 0; iPoint < polymarker.size(); ++iPoint) {
489 auto position = fObjectTransformation*G4Translate3D(polymarker[iPoint]);
491 auto currentEntity = new Qt3DCore::QEntity(currentNode); // Not undertood instancing yet
492 currentEntity->addComponent(material);
493 currentEntity->addComponent(transform);
494 currentEntity->addComponent(sphereMesh);
495 }
496 }
497 break;
499 {
500 G4Square square (polymarker); // Default square
501
502 const auto& colour = fpVisAttribs->GetColour();
503 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
504 material->setObjectName("materialForSquare");
505 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
506 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
507
508 auto boxMesh = new Qt3DExtras::QCuboidMesh();
509 boxMesh->setObjectName("boxMesh");
510 G4double side = markerSize;
511 if (markerSizeType == G4VSceneHandler::screen ) {
512 // Not figured out how to do screen-size, so use scene extent
513 const G4double scale = 200.; // Roughly pixles per scene
514 side *= fpScene->GetExtent().GetExtentRadius()/scale;
515 }
516 boxMesh->setXExtent(side);
517 boxMesh->setYExtent(side);
518 boxMesh->setZExtent(side);
519
520 for (std::size_t iPoint = 0; iPoint < polymarker.size(); ++iPoint) {
521 auto position = fObjectTransformation*G4Translate3D(polymarker[iPoint]);
523 auto currentEntity = new Qt3DCore::QEntity(currentNode);
524 currentEntity->addComponent(material);
525 currentEntity->addComponent(transform);
526 currentEntity->addComponent(boxMesh);
527 }
528 }
529 break;
530 }
531}
532
533#ifdef G4QT3DDEBUG
535 G4cout <<
536 "G4Qt3DSceneHandler::AddPrimitive(const G4Text& text) called.\n"
537 << text
538 << G4endl;
539#else
541#endif
542
543 static G4bool first = true;
544 if (first) {
545 first = false;
546 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Text&)",
547 "qt3D-0002", JustWarning,
548 "Text drawing doesn't work yet");
549 } // OK. Not working, but let it execute, which it does without error.
550
551 /* But it crashes after /vis/viewer/rebuild!!!
552 auto currentNode = CreateNewNode();
553 if (!currentNode) {
554 static G4bool first = true;
555 if (first) {
556 first = false;
557 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Text&)",
558 "qt3d-0003", JustWarning,
559 "No available node!");
560 }
561 return;
562 }
563
564 fpVisAttribs = fpViewer->GetApplicableVisAttributes(text.GetVisAttributes());
565
566 auto position = fObjectTransformation*G4Translate3D(text.GetPosition());
567 auto transform = G4Qt3DUtils::CreateQTransformFrom(position);
568// transform->setScale(10);
569 transform->setScale(0.1);
570
571// auto currentEntity = new Qt3DCore::QEntity(currentNode);
572
573 // This simply does not work
574 auto qtext = new Qt3DExtras::QText2DEntity();
575 qtext->setParent(currentNode);
576// qtext->setParent(currentEntity); // ?? Doesn't help
577 qtext->setText(text.GetText().c_str());
578// qtext->setHeight(100);
579// qtext->setWidth(1000);
580 qtext->setHeight(20);
581 qtext->setWidth(100);
582 qtext->setColor(Qt::green);
583 qtext->setFont(QFont("Courier New", 10));
584 qtext->addComponent(transform);
585
586 // This produces text in 3D facing +z - not what we want
587// const auto& colour = GetTextColour(text);
588// auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
589// material->setObjectName("materialForText");
590// material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
591// if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
592//
593// auto textMesh = new Qt3DExtras::QExtrudedTextMesh();
594// textMesh->setText(text.GetText().c_str());
595// textMesh->setFont(QFont("Courier New", 10));
596// textMesh->setDepth(.01f);
597//
598// currentNode->addComponent(material);
599// currentNode->addComponent(transform);
600// currentNode->addComponent(textMesh);
601 */
602}
603
605{
606#ifdef G4QT3DDEBUG
607 G4cout <<
608 "G4Qt3DSceneHandler::AddPrimitive(const G4Circle& circle) called.\n"
609 << circle
610 << G4endl;
611#endif
612
613#ifdef G4QT3DDEBUG
614 MarkerSizeType sizeType;
615 G4double size = GetMarkerSize (circle, sizeType);
616 switch (sizeType) {
617 default:
618 case screen:
619 // Draw in screen coordinates.
620 G4cout << "screen";
621 break;
622 case world:
623 // Draw in world coordinates.
624 G4cout << "world";
625 break;
626 }
627 G4cout << " size: " << size << G4endl;
628#endif
629
630 auto currentNode = CreateNewNode();
631 if (!currentNode) {
632 static G4bool first = true;
633 if (first) {
634 first = false;
635 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Circle&)",
636 "qt3d-0003", JustWarning,
637 "No available node!");
638 }
639 return;
640 }
641
643
646
647 const auto& colour = fpVisAttribs->GetColour();
648 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
649 material->setObjectName("materialForCircle");
650 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
651 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
652
653 auto sphereMesh = new Qt3DExtras::QSphereMesh;
654 sphereMesh->setObjectName("sphereMesh");
655 G4double radius;
656 if (circle.GetSizeType() == G4VMarker::world ) {
657 radius =circle.GetWorldRadius();
658 } else { // Screen-size or none
659 // Not figured out how to do screen-size, so use scene extent
660 const G4double scale = 200.; // Roughly pixles per scene
661 radius = circle.GetScreenRadius()*fpScene->GetExtent().GetExtentRadius()/scale;
662 }
663 sphereMesh->setRadius(radius);
664
665 auto currentEntity = new Qt3DCore::QEntity(currentNode);
666 currentEntity->addComponent(material);
667 currentEntity->addComponent(transform);
668 currentEntity->addComponent(sphereMesh);
669}
670
672{
673#ifdef G4QT3DDEBUG
674 G4cout <<
675 "G4Qt3DSceneHandler::AddPrimitive(const G4Square& square) called.\n"
676 << square
677 << G4endl;
678#endif
679
680#ifdef G4QT3DDEBUG
681 MarkerSizeType sizeType;
682 G4double size = GetMarkerSize (square, sizeType);
683 switch (sizeType) {
684 default:
685 case screen:
686 // Draw in screen coordinates.
687 G4cout << "screen";
688 break;
689 case world:
690 // Draw in world coordinates.
691 G4cout << "world";
692 break;
693 }
694 G4cout << " size: " << size << G4endl;
695#endif
696
697 auto currentNode = CreateNewNode();
698 if (!currentNode) {
699 static G4bool first = true;
700 if (first) {
701 first = false;
702 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Square&)",
703 "qt3d-0003", JustWarning,
704 "No available node!");
705 }
706 return;
707 }
708
710
713
714 const auto& colour = fpVisAttribs->GetColour();
715 auto material = new Qt3DExtras::QDiffuseSpecularMaterial();
716 material->setObjectName("materialForSquare");
717 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
718 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
719
720 auto boxMesh = new Qt3DExtras::QCuboidMesh();
721 boxMesh->setObjectName("boxMesh");
722 G4double side;
723 if (square.GetSizeType() == G4VMarker::world ) {
724 side = square.GetWorldDiameter();
725 } else { // Screen-size or none
726 // Not figured out how to do screen-size, so use scene extent
727 const G4double scale = 200.; // Roughly pixles per scene
728 side = square.GetScreenDiameter()*fpScene->GetExtent().GetExtentRadius()/scale;
729 }
730 boxMesh->setXExtent(side);
731 boxMesh->setYExtent(side);
732 boxMesh->setZExtent(side);
733
734 auto currentEntity = new Qt3DCore::QEntity(currentNode);
735 currentEntity->addComponent(material);
736 currentEntity->addComponent(transform);
737 currentEntity->addComponent(boxMesh);
738}
739
741{
742 auto currentNode = CreateNewNode();
743 if (!currentNode) {
744 static G4bool first = true;
745 if (first) {
746 first = false;
747 G4Exception("G4Qt3DSceneHandler::AddPrimitive(const G4Polyhedron&)",
748 "qt3d-0003", JustWarning,
749 "No available node!");
750 }
751 return;
752 }
753
754 if (polyhedron.GetNoFacets() == 0) return;
755
757
758 // Roll out vertices and normals for the faces. Note that this means vertices
759 // are duplicated. For example a box has 8 vertices, but to define 6 faces
760 // you need 12 triangles and 36 vertices. If it was just a matter of vertices
761 // we could restrict the number to 8 and use the indices to define the
762 // triangles, but we also have to consider the normals. A vertex can be have
763 // more than one normal, depending on which face it is being used to define.
764 // So we roll out all the vertices and normals for each triangle.
765 std::vector<G4Point3D> vertices;
766 std::vector<G4Normal3D> normals;
767
768 // Also roll out edges (as lines) for wireframe. Avoid duplicate lines,
769 // including those that differ only in the order of vertices.
770 typedef std::pair<G4Point3D,G4Point3D> Line;
771 std::vector<Line> lines;
772 auto insertIfNew = [&lines](const Line& newLine) {
773 // For a large polyhedron, eliminating lines like this is prohibitively
774 // expensive. Comment out for now, and maybe unwind altogether in future.
775 // Allow the graphics-reps utilities to optimise things like this.
776// for (const auto& line: lines) {
777// if ((newLine.first==line.first && newLine.second==line.second) ||
778// (newLine.first==line.second && newLine.second==line.first))
779// return;
780// }
781 lines.push_back(newLine);
782 };
783
784 G4bool isAuxilaryEdgeVisible = fpViewer->GetViewParameters().IsAuxEdgeVisible();
785 G4bool notLastFace;
786 do {
787 G4int nEdges;
788 G4Point3D vertex [4];
789 G4int edgeFlag[4];
790 G4Normal3D normal [4];
791 notLastFace = polyhedron.GetNextFacet(nEdges, vertex, edgeFlag, normal);
792 vertices.push_back(vertex[0]);
793 vertices.push_back(vertex[1]);
794 vertices.push_back(vertex[2]);
795 normals.push_back(normal[0]);
796 normals.push_back(normal[1]);
797 normals.push_back(normal[2]);
798 if(isAuxilaryEdgeVisible||edgeFlag[0]>0)insertIfNew(Line(vertex[0],vertex[1]));
799 if(isAuxilaryEdgeVisible||edgeFlag[1]>0)insertIfNew(Line(vertex[1],vertex[2]));
800 if (nEdges == 3) {
801 // Face is a triangle
802 // One more line for wireframe, triangles for surfaces are complete
803 if(isAuxilaryEdgeVisible||edgeFlag[2]>0)insertIfNew(Line(vertex[2],vertex[0]));
804 } else if (nEdges == 4) {
805 // Face is a quadrilateral
806 // Create another triangle for surfaces, add two more lines for wireframe
807 vertices.push_back(vertex[2]);
808 vertices.push_back(vertex[3]);
809 vertices.push_back(vertex[0]);
810 normals.push_back(normal[2]);
811 normals.push_back(normal[3]);
812 normals.push_back(normal[0]);
813 if(isAuxilaryEdgeVisible||edgeFlag[2]>0)insertIfNew(Line(vertex[2],vertex[3]));
814 if(isAuxilaryEdgeVisible||edgeFlag[3]>0)insertIfNew(Line(vertex[3],vertex[0]));
815 } else {
816 G4warn
817 << "ERROR: polyhedron face with unexpected number of edges (" << nEdges << ')'
818 << "\n Tag: " << fpModel->GetCurrentTag()
819 << G4endl;
820 return;
821 }
822 } while (notLastFace);
823 const auto nVerts = vertices.size();
824 const auto nLines = lines.size();
825
826 // Now put stuff into Qt objects
827
829 transform->setObjectName("transform");
830
831 Qt3DCore::QEntity* wireframeEntity = nullptr;
832 Qt3DCore::QEntity* surfaceEntity = nullptr;
833 static G4int errorCount = 0;
835 switch (drawing_style) {
837 wireframeEntity = new Qt3DCore::QEntity(currentNode);
838 wireframeEntity->addComponent(transform);
839 break;
841 wireframeEntity = new Qt3DCore::QEntity(currentNode);
842 wireframeEntity->addComponent(transform);
843 surfaceEntity = new Qt3DCore::QEntity(currentNode);
844 surfaceEntity->addComponent(transform);
845 break;
847 surfaceEntity = new Qt3DCore::QEntity(currentNode);
848 surfaceEntity->addComponent(transform);
849 break;
851 wireframeEntity = new Qt3DCore::QEntity(currentNode);
852 wireframeEntity->addComponent(transform);
853 surfaceEntity = new Qt3DCore::QEntity(currentNode);
854 surfaceEntity->addComponent(transform);
855 break;
857 // Shouldn't happen in this function (it's a polyhedron!)
858 if (errorCount == 0) {
859 ++errorCount;
860 G4warn << "WARNING: Qt3D: cloud drawing not implemented" << G4endl;
861 }
862 return;
863 break;
864 }
865
866 const auto vertexByteSize = 3*sizeof(PRECISION);
867
868 G4Qt3DCompat::QGeometry* vertexGeometry = nullptr;
869 G4Qt3DCompat::QGeometry* lineGeometry = nullptr;
870
871 G4Qt3DCompat::QAttribute* positionAtt = nullptr;
872 G4Qt3DCompat::QAttribute* normalAtt = nullptr;
873 G4Qt3DCompat::QAttribute* lineAtt = nullptr;
874 G4Qt3DCompat::QAttribute* dummyNormalLineAtt = nullptr;
875
876 G4Qt3DCompat::QBuffer* vertexBuffer = nullptr;
877 if (drawing_style == G4ViewParameters::hlr ||
878 drawing_style == G4ViewParameters::hsr ||
879 drawing_style == G4ViewParameters::hlhsr) {
880
881 // Put vertices, normals into QByteArray
882 // Accomodates both vertices and normals - hence 2*
883 QByteArray vertexByteArray;
884 const auto vertexBufferByteSize = 2*nVerts*vertexByteSize;
885 vertexByteArray.resize((G4int)vertexBufferByteSize);
886 auto vertexBufferArray = reinterpret_cast<PRECISION*>(vertexByteArray.data());
887 G4int i1 = 0;
888 for (std::size_t i = 0; i < nVerts; ++i) {
889 vertexBufferArray[i1++] = vertices[i].x();
890 vertexBufferArray[i1++] = vertices[i].y();
891 vertexBufferArray[i1++] = vertices[i].z();
892 vertexBufferArray[i1++] = normals[i].x();
893 vertexBufferArray[i1++] = normals[i].y();
894 vertexBufferArray[i1++] = normals[i].z();
895 }
896 // Vertex buffer (vertices and normals)
897 vertexGeometry = new G4Qt3DCompat::QGeometry();
898 vertexGeometry->setObjectName("vertexGeometry");
899 vertexBuffer = new G4Qt3DCompat::QBuffer(vertexGeometry);
900 vertexBuffer->setObjectName("Vertex buffer");
901 vertexBuffer->setData(vertexByteArray);
902
903 // Position attribute
904 positionAtt = new G4Qt3DCompat::QAttribute;
905 positionAtt->setObjectName("Position attribute");
906 positionAtt->setName(G4Qt3DCompat::QAttribute::defaultPositionAttributeName());
907 positionAtt->setBuffer(vertexBuffer);
908 positionAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
909 positionAtt->setVertexBaseType(BASETYPE);
910 positionAtt->setVertexSize(3);
911 positionAtt->setCount((G4int)nVerts);
912 positionAtt->setByteOffset(0);
913 positionAtt->setByteStride(2*vertexByteSize);
914
915 // Normal attribute
916 normalAtt = new G4Qt3DCompat::QAttribute;
917 normalAtt->setObjectName("Normal attribute");
918 normalAtt->setName(G4Qt3DCompat::QAttribute::defaultNormalAttributeName());
919 normalAtt->setBuffer(vertexBuffer);
920 normalAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
921 normalAtt->setVertexBaseType(BASETYPE);
922 normalAtt->setVertexSize(3);
923 normalAtt->setCount((G4int)nVerts);
924 normalAtt->setByteOffset(vertexByteSize);
925 normalAtt->setByteStride(2*vertexByteSize);
926 }
927
928 G4Qt3DCompat::QBuffer* lineBuffer = nullptr;
929 if (drawing_style == G4ViewParameters::wireframe ||
930 drawing_style == G4ViewParameters::hlr ||
931 drawing_style == G4ViewParameters::hlhsr) {
932
933 // Put lines into a QByteArray
934 QByteArray lineByteArray;
935 const auto lineBufferByteSize = 2*nLines*vertexByteSize;
936 lineByteArray.resize((G4int)lineBufferByteSize);
937 auto lineBufferArray = reinterpret_cast<PRECISION*>(lineByteArray.data());
938 G4int i2 = 0;
939 for (const auto& line: lines) {
940 lineBufferArray[i2++] = line.first.x();
941 lineBufferArray[i2++] = line.first.y();
942 lineBufferArray[i2++] = line.first.z();
943 lineBufferArray[i2++] = line.second.x();
944 lineBufferArray[i2++] = line.second.y();
945 lineBufferArray[i2++] = line.second.z();
946 }
947 // Line loop buffer
948 lineGeometry = new G4Qt3DCompat::QGeometry();
949 lineGeometry->setObjectName("lineGeometry");
950 lineBuffer = new G4Qt3DCompat::QBuffer(lineGeometry);
951 lineBuffer->setObjectName("Line buffer");
952 lineBuffer->setData(lineByteArray);
953
954 // Line attribute
955 lineAtt = new G4Qt3DCompat::QAttribute;
956 lineAtt->setObjectName("Position attribute");
957 lineAtt->setName(G4Qt3DCompat::QAttribute::defaultPositionAttributeName());
958 lineAtt->setBuffer(lineBuffer);
959 lineAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
960 lineAtt->setVertexBaseType(BASETYPE);
961 lineAtt->setVertexSize(3);
962 lineAtt->setCount((G4int)nLines);
963 lineAtt->setByteOffset(0);
964 lineAtt->setByteStride(vertexByteSize);
965 // Normal attribute (a dummy with count==0) (Qt6 seems to require)
966 dummyNormalLineAtt = new G4Qt3DCompat::QAttribute;
967 dummyNormalLineAtt->setObjectName("Normal attribute");
968 dummyNormalLineAtt->setName(G4Qt3DCompat::QAttribute::defaultNormalAttributeName());
969 dummyNormalLineAtt->setBuffer(lineBuffer);
970 dummyNormalLineAtt->setAttributeType(G4Qt3DCompat::QAttribute::VertexAttribute);
971 dummyNormalLineAtt->setVertexBaseType(BASETYPE);
972 dummyNormalLineAtt->setVertexSize(3);
973 dummyNormalLineAtt->setCount(0);
974 dummyNormalLineAtt->setByteOffset(0);
975 dummyNormalLineAtt->setByteStride(vertexByteSize);
976 }
977
978 // Create material and renderer(s)...
979
980 const auto& colour = fpVisAttribs->GetColour();
981 Qt3DExtras::QDiffuseSpecularMaterial* material;
982 Qt3DRender::QGeometryRenderer* renderer;
983 switch (drawing_style) {
984
986
987 lineGeometry->addAttribute(lineAtt);
988 lineGeometry->addAttribute(dummyNormalLineAtt);
989
990 material = new Qt3DExtras::QDiffuseSpecularMaterial();
991 material->setObjectName("materialForWireframe");
992 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
993 material->setShininess(0.);
994 material->setSpecular(0.);
995 wireframeEntity->addComponent(material);
996
997 renderer = new Qt3DRender::QGeometryRenderer;
998 renderer->setObjectName("polyhedronWireframeRenderer");
999 renderer->setGeometry(lineGeometry);
1000 renderer->setVertexCount(2*(G4int)nLines);
1001 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
1002 wireframeEntity->addComponent(renderer);
1003
1004 break;
1005
1007
1008 // Surfaces with background colour to hide the edges
1009
1010 vertexGeometry->addAttribute(positionAtt);
1011 vertexGeometry->addAttribute(normalAtt);
1012
1013 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1014 material->setObjectName("materialForHiddenLines");
1015 material->setAmbient(Qt::white); // White for now (should be from fVP)
1016 material->setShininess(0.);
1017 material->setSpecular(0.);
1018 surfaceEntity->addComponent(material);
1019
1020 renderer = new Qt3DRender::QGeometryRenderer;
1021 renderer->setObjectName("polyhedronSurfaceRenderer");
1022 renderer->setGeometry(vertexGeometry);
1023 renderer->setVertexCount((G4int)nVerts);
1024 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Triangles);
1025 surfaceEntity->addComponent(renderer);
1026
1027 // Edges
1028
1029 lineGeometry->addAttribute(lineAtt);
1030 lineGeometry->addAttribute(dummyNormalLineAtt);
1031
1032 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1033 material->setObjectName("materialForWireFrame");
1034 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
1035 material->setShininess(0.);
1036 material->setSpecular(0.);
1037 wireframeEntity->addComponent(material);
1038
1039 renderer = new Qt3DRender::QGeometryRenderer;
1040 renderer->setObjectName("polyhedronWireframeRenderer");
1041 renderer->setGeometry(lineGeometry);
1042 renderer->setVertexCount(2*(G4int)nLines);
1043 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
1044 wireframeEntity->addComponent(renderer);
1045
1046 break;
1047
1049
1050 vertexGeometry->addAttribute(positionAtt);
1051 vertexGeometry->addAttribute(normalAtt);
1052
1053 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1054 material->setObjectName("materialForSurface");
1055 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
1056 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
1057 surfaceEntity->addComponent(material);
1058
1059 renderer = new Qt3DRender::QGeometryRenderer;
1060 renderer->setObjectName("polyhedronSurfaceRenderer");
1061 renderer->setGeometry(vertexGeometry);
1062 renderer->setVertexCount((G4int)nVerts);
1063 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Triangles);
1064 surfaceEntity->addComponent(renderer);
1065
1066 break;
1067
1069
1070 // Surfaces
1071
1072 vertexGeometry->addAttribute(positionAtt);
1073 vertexGeometry->addAttribute(normalAtt);
1074
1075 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1076 material->setObjectName("materialForSurface");
1077 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
1078 if (colour.GetAlpha() < 1.) material->setAlphaBlendingEnabled(true);
1079 surfaceEntity->addComponent(material);
1080
1081 renderer = new Qt3DRender::QGeometryRenderer;
1082 renderer->setObjectName("polyhedronSurfaceRenderer");
1083 renderer->setGeometry(vertexGeometry);
1084 renderer->setVertexCount((G4int)nVerts);
1085 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Triangles);
1086 surfaceEntity->addComponent(renderer);
1087
1088 // Edges
1089
1090 lineGeometry->addAttribute(lineAtt);
1091 lineGeometry->addAttribute(dummyNormalLineAtt);
1092
1093 material = new Qt3DExtras::QDiffuseSpecularMaterial();
1094 material->setObjectName("materialForWireframe");
1095 material->setAmbient(G4Qt3DUtils::ConvertToQColor(colour));
1096 material->setShininess(0.);
1097 material->setSpecular(0.);
1098 wireframeEntity->addComponent(material);
1099
1100 renderer = new Qt3DRender::QGeometryRenderer;
1101 renderer->setObjectName("polyhedronSurfaceRenderer");
1102 renderer->setGeometry(lineGeometry);
1103 renderer->setVertexCount(2*(G4int)nLines);
1104 renderer->setPrimitiveType(Qt3DRender::QGeometryRenderer::Lines);
1105 wireframeEntity->addComponent(renderer);
1106
1107 break;
1108
1110 // Case trapped at start of function, so no need to implement
1111 break;
1112 }
1113}
1114
1119
1125
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define PRECISION
#define BASETYPE
#define G4warn
Definition G4Scene.cc:41
HepGeom::Translate3D G4Translate3D
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const std::vector< G4PhysicalVolumeNodeID > & GetFullPVPath() const
MarkerType GetMarkerType() const
void SetPVNodeID(const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID &id)
const G4PhysicalVolumeModel::G4PhysicalVolumeNodeID & GetPVNodeID() const
void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
void BeginPrimitives2D(const G4Transform3D &objectTransformation)
Qt3DCore::QEntity * fpTransientObjects
std::vector< G4Qt3DQEntity * > fpPhysicalVolumeObjects
Qt3DCore::QEntity * fpQt3DScene
void AddCompound(const G4Mesh &)
void AddPrimitive(const G4Polyline &)
Qt3DCore::QEntity * fpPersistentObjects
void BeginPrimitives(const G4Transform3D &objectTransformation)
G4Qt3DSceneHandler(G4VGraphicsSystem &system, const G4String &name)
G4Qt3DQEntity * CreateNewNode()
const std::vector< Model > & GetRunDurationModelList() const
const G4VisExtent & GetExtent() const
SizeType GetSizeType() const
Definition G4VMarker.cc:79
G4double GetScreenRadius() const
G4double GetWorldDiameter() const
G4Point3D GetPosition() const
G4double GetScreenDiameter() const
G4double GetWorldRadius() const
virtual G4String GetCurrentTag() const
Definition G4VModel.cc:46
const G4String & GetGlobalTag() const
const G4String & GetName() const
G4Transform3D fObjectTransformation
virtual void EndPrimitives()
G4double GetMarkerSize(const G4VMarker &, MarkerSizeType &)
virtual void PreAddSolid(const G4Transform3D &objectTransformation, const G4VisAttributes &)
virtual void PostAddSolid()
virtual void EndPrimitives2D()
const G4VisAttributes * fpVisAttribs
virtual void BeginPrimitives2D(const G4Transform3D &objectTransformation=G4Transform3D())
G4ViewParameters::DrawingStyle GetDrawingStyle(const G4VisAttributes *)
virtual void BeginPrimitives(const G4Transform3D &objectTransformation=G4Transform3D())
void StandardSpecialMeshRendering(const G4Mesh &)
const G4VisAttributes * GetApplicableVisAttributes(const G4VisAttributes *) const
const G4ViewParameters & GetViewParameters() const
G4bool IsAuxEdgeVisible() const
const G4Colour & GetColour() const
G4double GetExtentRadius() const
const G4VisAttributes * GetVisAttributes() const
G4int GetNoFacets() const
G4bool GetNextFacet(G4int &n, G4Point3D *nodes, G4int *edgeFlags=nullptr, G4Normal3D *normals=nullptr) const
Qt3DCore::QTransform * CreateQTransformFrom(const G4Transform3D &)
QColor ConvertToQColor(const G4Colour &c)
void delete_components_and_children_of_entity_recursively(Qt3DCore::QNode *node)