Geant4 Cross Reference |
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 #include "G4Qt3DUtils.hh" 29 30 #include "G4Qt3DQEntity.hh" 31 #include "G4PhysicalVolumeModel.hh" 32 33 Qt3DCore::QTransform* G4Qt3DUtils::CreateQTransformFrom(const G4Transform3D& g) 34 { 35 auto* q = new Qt3DCore::QTransform; 36 q->setMatrix 37 (QMatrix4x4 38 (g.xx(),g.xy(),g.xz(),g.dx(), 39 g.yx(),g.yy(),g.yz(),g.dy(), 40 g.zx(),g.zy(),g.zz(),g.dz(), 41 0,0,0,1)); 42 q->setObjectName("transform"); 43 return q; 44 } 45 46 QColor G4Qt3DUtils::ConvertToQColor(const G4Colour& c) { 47 QColor qColor; 48 qColor.setRgbF(c.GetRed(),c.GetGreen(),c.GetBlue(),c.GetAlpha()); 49 return qColor; 50 } 51 52 QVector3D G4Qt3DUtils::ConvertToQVector3D(const G4ThreeVector& v) { 53 return QVector3D(v.x(),v.y(),v.z()); 54 } 55 56 // https://stackoverflow.com/questions/45759274/how-can-i-delete-all-nodes-recursively-in-the-root-entity-of-qt3dwindow 57 void G4Qt3DUtils::delete_entity_recursively(Qt3DCore::QNode *node){ 58 #ifdef G4QT3DDEBUG 59 G4Qt3DUtils::LogFile << "node " << node->objectName().toStdString() << std::endl; 60 #endif 61 Qt3DCore::QEntity* entity = dynamic_cast<Qt3DCore::QEntity*>(node); 62 if(entity == nullptr){ 63 #ifdef G4QT3DDEBUG 64 G4String name = node->objectName().toStdString(); 65 if (name == "") name = "X"; 66 G4Qt3DUtils::LogFile << (void*)node << ": " 67 << "Deleting non-entity node " << name << std::endl; 68 #endif 69 delete node; 70 node = nullptr; 71 return; 72 } 73 for (auto component: entity->components()) { 74 #ifdef G4QT3DDEBUG 75 G4String name = component->objectName().toStdString(); 76 if (name == "") name = "X"; 77 G4Qt3DUtils::LogFile << (void*)node << ": " << "Deleting component " << name 78 << " of " << entity->objectName().toStdString() << std::endl; 79 #endif 80 entity->removeComponent(component); 81 delete component; 82 component = nullptr; 83 } 84 for (auto child_node: entity->childNodes()) { 85 G4String name = child_node->objectName().toStdString(); 86 if (name == "") name = "X"; 87 #ifdef G4QT3DDEBUG 88 G4Qt3DUtils::LogFile << (void*)child_node << ": " << "Child node " << name 89 << " of " << entity->objectName().toStdString() << std::endl; 90 #endif 91 delete_entity_recursively(child_node); 92 } 93 G4String name = entity->objectName().toStdString(); 94 if (name == "") name = "X"; 95 #ifdef G4QT3DDEBUG 96 G4Qt3DUtils::LogFile << (void*)entity << ": " << "Deleting entity " << name << std::endl; 97 #endif 98 delete entity; 99 entity = nullptr; 100 } 101 102 void G4Qt3DUtils::delete_components_and_children_of_entity_recursively(Qt3DCore::QNode *node){ 103 Qt3DCore::QEntity* entity = dynamic_cast<Qt3DCore::QEntity*>(node); 104 if(entity == nullptr){ 105 #ifdef G4QT3DDEBUG 106 G4String name = node->objectName().toStdString(); 107 if (name == "") name = "X"; 108 G4Qt3DUtils::LogFile << (void*)node << ": " << "Found non-entity node " << name << std::endl; 109 #endif 110 return; 111 } 112 for (auto component: entity->components()){ 113 #ifdef G4QT3DDEBUG 114 G4String name = component->objectName().toStdString(); 115 if (name == "") name = "X"; 116 G4Qt3DUtils::LogFile << (void*)entity << ": " << "Deleting component " << name 117 << " of " << entity->objectName().toStdString() << std::endl; 118 #endif 119 entity->removeComponent(component); 120 delete(component); 121 component = nullptr; 122 } 123 auto child_nodes = entity->childNodes(); 124 for (auto child_node: child_nodes) { 125 G4String name = child_node->objectName().toStdString(); 126 if (name == "") name = "X"; 127 #ifdef G4QT3DDEBUG 128 G4Qt3DUtils::LogFile << (void*)child_node << ": " << "Child node " << name 129 << " of " << entity->objectName().toStdString() << std::endl; 130 #endif 131 delete_entity_recursively(child_node); 132 } 133 G4String name = entity->objectName().toStdString(); 134 if (name == "") name = "X"; 135 #ifdef G4QT3DDEBUG 136 G4Qt3DUtils::LogFile << (void*)entity << ": " << "Clearing child nodes of " << name << std::endl; 137 #endif 138 child_nodes.clear(); 139 } 140 141 #ifdef G4QT3DDEBUG 142 std::ofstream G4Qt3DUtils::LogFile("LogFile.txt"); 143 void G4Qt3DUtils::PrintQObjectTree 144 (const QObject* node, 145 const G4String& where) 146 { 147 auto& logFile = G4Qt3DUtils::LogFile; 148 if (where.length()) logFile << "\n===== QObjectTree at " << where << std::endl; 149 static G4int iDep = -1; 150 ++iDep; 151 G4String nodeName = node->objectName().toStdString(); 152 if (nodeName == "") nodeName = "X"; 153 for (G4int i = 0; i < iDep; ++i) logFile << " "; 154 logFile << (void*)node << ": " 155 << "Node at depth " << iDep << ": " << nodeName << ": " 156 << "thread: " << node->thread() << ": " 157 << "parent: " << node->parent() << ": "; 158 const auto* g4node = dynamic_cast<const G4Qt3DQEntity*>(node); 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