Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Qt3DUtils.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, 18th July 2020
27
28#if defined (G4VIS_BUILD_QT3D_DRIVER) || defined (G4VIS_USE_QT3D)
29
30#include "G4Qt3DUtils.hh"
31
32#include "G4Qt3DQEntity.hh"
34
35Qt3DCore::QTransform* G4Qt3DUtils::CreateQTransformFrom(const G4Transform3D& g)
36{
37 auto* q = new Qt3DCore::QTransform;
38 q->setMatrix
39 (QMatrix4x4
40 (g.xx(),g.xy(),g.xz(),g.dx(),
41 g.yx(),g.yy(),g.yz(),g.dy(),
42 g.zx(),g.zy(),g.zz(),g.dz(),
43 0,0,0,1));
44 q->setObjectName("transform");
45 return q;
46}
47
48QColor G4Qt3DUtils::ConvertToQColor(const G4Colour& c) {
49 QColor qColor;
50 qColor.setRgbF(c.GetRed(),c.GetGreen(),c.GetBlue(),c.GetAlpha());
51 return qColor;
52}
53
54QVector3D G4Qt3DUtils::ConvertToQVector3D(const G4ThreeVector& v) {
55 return QVector3D(v.x(),v.y(),v.z());
56}
57
58// https://stackoverflow.com/questions/45759274/how-can-i-delete-all-nodes-recursively-in-the-root-entity-of-qt3dwindow
59void G4Qt3DUtils::delete_entity_recursively(Qt3DCore::QNode *node){
60#ifdef G4QT3DDEBUG
61 G4Qt3DUtils::LogFile << "node " << node->objectName().toStdString() << std::endl;
62#endif
63 Qt3DCore::QEntity* entity = dynamic_cast<Qt3DCore::QEntity*>(node);
64 if(entity == nullptr){
65 G4String name = node->objectName().toStdString();
66 if (name == "") name = "X";
67#ifdef G4QT3DDEBUG
68 G4Qt3DUtils::LogFile << (void*)node << ": "
69 << "Deleting non-entity node " << name << std::endl;
70#endif
71 delete node;
72 node = nullptr;
73 return;
74 }
75 for (auto component: entity->components()) {
76 G4String name = component->objectName().toStdString();
77 if (name == "") name = "X";
78#ifdef G4QT3DDEBUG
79 G4Qt3DUtils::LogFile << (void*)node << ": " << "Deleting component " << name
80 << " of " << entity->objectName().toStdString() << std::endl;
81#endif
82 entity->removeComponent(component);
83 delete component;
84 component = nullptr;
85 }
86 for (auto child_node: entity->childNodes()) {
87 G4String name = child_node->objectName().toStdString();
88 if (name == "") name = "X";
89#ifdef G4QT3DDEBUG
90 G4Qt3DUtils::LogFile << (void*)child_node << ": " << "Child node " << name
91 << " of " << entity->objectName().toStdString() << std::endl;
92#endif
93 delete_entity_recursively(child_node);
94 }
95 G4String name = entity->objectName().toStdString();
96 if (name == "") name = "X";
97#ifdef G4QT3DDEBUG
98 G4Qt3DUtils::LogFile << (void*)entity << ": " << "Deleting entity " << name << std::endl;
99#endif
100 delete entity;
101 entity = nullptr;
102}
103
104void G4Qt3DUtils::delete_components_and_children_of_entity_recursively(Qt3DCore::QNode *node){
105 Qt3DCore::QEntity* entity = dynamic_cast<Qt3DCore::QEntity*>(node);
106 if(entity == nullptr){
107 G4String name = node->objectName().toStdString();
108 if (name == "") name = "X";
109#ifdef G4QT3DDEBUG
110 G4Qt3DUtils::LogFile << (void*)node << ": " << "Found non-entity node " << name << std::endl;
111#endif
112 return;
113 }
114 for (auto component: entity->components()){
115 G4String name = component->objectName().toStdString();
116 if (name == "") name = "X";
117#ifdef G4QT3DDEBUG
118 G4Qt3DUtils::LogFile << (void*)entity << ": " << "Deleting component " << name
119 << " of " << entity->objectName().toStdString() << std::endl;
120#endif
121 entity->removeComponent(component);
122 delete(component);
123 component = nullptr;
124 }
125 auto child_nodes = entity->childNodes();
126 for (auto child_node: child_nodes) {
127 G4String name = child_node->objectName().toStdString();
128 if (name == "") name = "X";
129#ifdef G4QT3DDEBUG
130 G4Qt3DUtils::LogFile << (void*)child_node << ": " << "Child node " << name
131 << " of " << entity->objectName().toStdString() << std::endl;
132#endif
133 delete_entity_recursively(child_node);
134 }
135 G4String name = entity->objectName().toStdString();
136 if (name == "") name = "X";
137#ifdef G4QT3DDEBUG
138 G4Qt3DUtils::LogFile << (void*)entity << ": " << "Clearing child nodes of " << name << std::endl;
139#endif
140 child_nodes.clear();
141}
142
143#ifdef G4QT3DDEBUG
144std::ofstream G4Qt3DUtils::LogFile("LogFile.txt");
145void G4Qt3DUtils::PrintQObjectTree
146 (const QObject* node,
147 const G4String& where)
148{
149 auto& logFile = G4Qt3DUtils::LogFile;
150 if (where.length()) logFile << "\n===== QObjectTree at " << where << std::endl;
151 static G4int iDep = -1;
152 ++iDep;
153 const auto* g4node = dynamic_cast<const G4Qt3DQEntity*>(node);
154 G4String nodeName = node->objectName().toStdString();
155 if (nodeName == "") nodeName = "X";
156 for (G4int i = 0; i < iDep; ++i) logFile << " ";
157 logFile << (void*)node << ": "
158 << "Node at depth " << iDep << ": " << nodeName << ": ";
159 if (g4node) {
160 logFile << g4node->GetPVNodeID() << std::endl;
161 } else {
162 logFile << typeid(node).name() << std::endl;
163 }
164 if (g4node) {
165 for (const auto& component: g4node->components()) {
166 G4String name = component->objectName().toStdString();
167 if (name == "") name = "X";
168 for (G4int i = 0; i < iDep; ++i) logFile << " ";
169 logFile << (void*)component << ": "<< "Component at depth " << iDep << " "
170 << name << " of " << nodeName << std::endl;
171 }
172 }
173 for (const auto& child: node->children()) {
174 PrintQObjectTree(child);
175 }
176 --iDep;
177 if (where.length()) logFile << "===== End: QObjectTree at " << where << std::endl;
178 return;
179}
180#endif
181
182#endif // #if defined (G4VIS_BUILD_QT3D_DRIVER) || defined (G4VIS_USE_QT3D)
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
G4double GetBlue() const
Definition: G4Colour.hh:152
G4double GetAlpha() const
Definition: G4Colour.hh:153
G4double GetRed() const
Definition: G4Colour.hh:150
G4double GetGreen() const
Definition: G4Colour.hh:151
double dy() const
Definition: Transform3D.h:287
double zz() const
Definition: Transform3D.h:281
double yz() const
Definition: Transform3D.h:272
double dz() const
Definition: Transform3D.h:290
double dx() const
Definition: Transform3D.h:284
double xy() const
Definition: Transform3D.h:260
double zx() const
Definition: Transform3D.h:275
double yx() const
Definition: Transform3D.h:266
double zy() const
Definition: Transform3D.h:278
double xx() const
Definition: Transform3D.h:257
double yy() const
Definition: Transform3D.h:269
double xz() const
Definition: Transform3D.h:263
const char * name(G4int ptype)