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