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