Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 25 26 #include <vtkPoints.h> << 26 #include "vtkPoints.h" 27 #include <vtkCellArray.h> << 27 #include "vtkCellArray.h" 28 #include <vtkCharArray.h> << 28 #include "vtkCharArray.h" 29 #include <vtkUnstructuredGrid.h> << 29 #include "vtkUnstructuredGrid.h" 30 #include <vtkDataSetMapper.h> << 30 #include "vtkDataSetMapper.h" 31 #include <vtkUnstructuredGridVolumeMapper.h> << 31 #include "vtkUnstructuredGridVolumeMapper.h" 32 #include <vtkUnstructuredGridVolumeRayCastMapp << 32 #include "vtkUnstructuredGridVolumeRayCastMapper.h" 33 #include <vtkUnstructuredGridVolumeZSweepMappe << 33 #include "vtkUnstructuredGridVolumeZSweepMapper.h" 34 #include <vtkOpenGLProjectedTetrahedraMapper.h << 34 #include "vtkOpenGLProjectedTetrahedraMapper.h" 35 #include <vtkActor.h> << 35 #include "vtkActor.h" 36 #include <vtkColorTransferFunction.h> << 36 #include "vtkColorTransferFunction.h" 37 #include <vtkPiecewiseFunction.h> << 37 #include "vtkPiecewiseFunction.h" 38 #include <vtkLookupTable.h> << 38 #include "vtkLookupTable.h" 39 #include <vtkVolumeProperty.h> << 39 #include "vtkVolumeProperty.h" 40 #include <vtkVolume.h> << 40 #include "vtkVolume.h" 41 #include <vtkNew.h> << 41 #include "vtkNew.h" 42 #include <vtkTetra.h> << 42 #include "vtkTetra.h" 43 #include <vtkVoxel.h> << 43 #include "vtkVoxel.h" 44 // #include <vtkCleanUnstructuredGrid.h> << 44 // #include "vtkCleanUnstructuredGrid.h" 45 #include <vtkDataSetTriangleFilter.h> << 45 #include "vtkDataSetTriangleFilter.h" 46 #include <vtkClipDataSet.h> << 46 #include "vtkClipDataSet.h" 47 47 48 #include "G4VtkUnstructuredGridPipeline.hh" 48 #include "G4VtkUnstructuredGridPipeline.hh" 49 #include "G4VtkVisContext.hh" 49 #include "G4VtkVisContext.hh" 50 #include "G4VtkViewer.hh" 50 #include "G4VtkViewer.hh" 51 51 52 #include "G4Mesh.hh" 52 #include "G4Mesh.hh" 53 #include "G4LogicalVolume.hh" 53 #include "G4LogicalVolume.hh" 54 #include "G4VNestedParameterisation.hh" 54 #include "G4VNestedParameterisation.hh" 55 #include "G4ModelingParameters.hh" 55 #include "G4ModelingParameters.hh" 56 #include "G4PhysicalVolumeModel.hh" 56 #include "G4PhysicalVolumeModel.hh" 57 #include "G4PseudoScene.hh" 57 #include "G4PseudoScene.hh" 58 #include "G4Transform3D.hh" 58 #include "G4Transform3D.hh" 59 #include "G4VSolid.hh" 59 #include "G4VSolid.hh" 60 #include "G4Tet.hh" 60 #include "G4Tet.hh" 61 #include "G4Material.hh" 61 #include "G4Material.hh" 62 62 63 G4VtkUnstructuredGridPipeline::G4VtkUnstructur 63 G4VtkUnstructuredGridPipeline::G4VtkUnstructuredGridPipeline(G4String nameIn, const G4VtkVisContext& vcIn) : 64 G4VVtkPipeline(nameIn, "G4VtkUnstructuredPip 64 G4VVtkPipeline(nameIn, "G4VtkUnstructuredPipeline", vcIn, false, vcIn.fViewer->renderer) { 65 65 66 points = vtkSmartPointer<vtkPoints>::Ne 66 points = vtkSmartPointer<vtkPoints>::New(); 67 67 68 pointColourValues = vtkSmartPointer<vtkDoubl 68 pointColourValues = vtkSmartPointer<vtkDoubleArray>::New(); 69 cellColourValues = vtkSmartPointer<vtkDoubl 69 cellColourValues = vtkSmartPointer<vtkDoubleArray>::New(); 70 pointColourValues->SetNumberOfComponents(4); 70 pointColourValues->SetNumberOfComponents(4); 71 cellColourValues->SetNumberOfComponents(4); 71 cellColourValues->SetNumberOfComponents(4); 72 72 73 pointColourIndices = vtkSmartPointer<vtkDoub 73 pointColourIndices = vtkSmartPointer<vtkDoubleArray>::New(); 74 cellColourIndices = vtkSmartPointer<vtkDoub 74 cellColourIndices = vtkSmartPointer<vtkDoubleArray>::New(); 75 pointColourIndices->SetNumberOfComponents(1) 75 pointColourIndices->SetNumberOfComponents(1); 76 cellColourIndices->SetNumberOfComponents(1); 76 cellColourIndices->SetNumberOfComponents(1); 77 77 78 colourLUT = vtkSmartPointer<vtkDiscretizable 78 colourLUT = vtkSmartPointer<vtkDiscretizableColorTransferFunction>::New(); 79 colourLUT->DiscretizeOn(); 79 colourLUT->DiscretizeOn(); 80 80 81 unstructuredGrid = vtkSmartPointer<vtkUnstru 81 unstructuredGrid = vtkSmartPointer<vtkUnstructuredGrid>::New(); 82 unstructuredGrid->SetPoints(points); 82 unstructuredGrid->SetPoints(points); 83 unstructuredGrid->GetPointData()->SetScalars 83 unstructuredGrid->GetPointData()->SetScalars(pointColourValues); 84 unstructuredGrid->GetCellData()->SetScalars( 84 unstructuredGrid->GetCellData()->SetScalars(cellColourValues); 85 85 86 // clean filter 86 // clean filter 87 #if 0 87 #if 0 88 clean = vtkSmartPointer<vtkStaticCleanUnstru 88 clean = vtkSmartPointer<vtkStaticCleanUnstructuredGrid>::New(); 89 clean->SetInputData(unstructuredGrid); 89 clean->SetInputData(unstructuredGrid); 90 clean->ToleranceIsAbsoluteOff(); 90 clean->ToleranceIsAbsoluteOff(); 91 clean->SetTolerance(1e-6); 91 clean->SetTolerance(1e-6); 92 #endif 92 #endif 93 93 94 // Clip filter 94 // Clip filter 95 clip = vtkSmartPointer<vtkClipDataSet>::New( 95 clip = vtkSmartPointer<vtkClipDataSet>::New(); 96 vtkNew<vtkPlane> plane; 96 vtkNew<vtkPlane> plane; 97 clip->SetClipFunction(plane); 97 clip->SetClipFunction(plane); 98 clip->SetInputData(unstructuredGrid); 98 clip->SetInputData(unstructuredGrid); 99 99 100 // Triangle filter 100 // Triangle filter 101 auto tri = vtkSmartPointer<vtkDataSetTriangl 101 auto tri = vtkSmartPointer<vtkDataSetTriangleFilter>::New(); 102 tri->SetInputData(unstructuredGrid); 102 tri->SetInputData(unstructuredGrid); 103 103 104 // create dataset mapper 104 // create dataset mapper 105 mapper = vtkSmartPointer<vtkDataSetMapper>:: 105 mapper = vtkSmartPointer<vtkDataSetMapper>::New(); 106 mapper->SetScalarModeToUseCellData(); 106 mapper->SetScalarModeToUseCellData(); 107 mapper->SetColorModeToDirectScalars(); 107 mapper->SetColorModeToDirectScalars(); 108 mapper->SetInputData(unstructuredGrid); 108 mapper->SetInputData(unstructuredGrid); 109 //mapper->SetInputConnection(clip->GetOutput 109 //mapper->SetInputConnection(clip->GetOutputPort()); 110 //mapper->SetInputConnection(tri->GetOutputP 110 //mapper->SetInputConnection(tri->GetOutputPort()); 111 111 112 // create volume mapper 112 // create volume mapper 113 volumeMapper = vtkSmartPointer<vtkUnstructur 113 volumeMapper = vtkSmartPointer<vtkUnstructuredGridVolumeRayCastMapper>::New(); 114 volumeMapper->SetScalarModeToUseCellData(); 114 volumeMapper->SetScalarModeToUseCellData(); 115 volumeMapper->SetInputConnection(tri->GetOut 115 volumeMapper->SetInputConnection(tri->GetOutputPort()); 116 116 117 // create volume properties 117 // create volume properties 118 auto volumeProp = vtkSmartPointer<vtkVolumeP 118 auto volumeProp = vtkSmartPointer<vtkVolumeProperty>::New(); 119 volumeProp->SetColor(colourLUT); 119 volumeProp->SetColor(colourLUT); 120 120 121 // create actor 121 // create actor 122 actor = vtkSmartPointer<vtkActor>::New(); 122 actor = vtkSmartPointer<vtkActor>::New(); 123 actor->SetMapper(mapper); 123 actor->SetMapper(mapper); 124 actor->SetVisibility(1); 124 actor->SetVisibility(1); 125 125 126 // create volume 126 // create volume 127 volume = vtkSmartPointer<vtkVolume>::New(); 127 volume = vtkSmartPointer<vtkVolume>::New(); 128 volume->SetMapper(volumeMapper); 128 volume->SetMapper(volumeMapper); 129 //volume->SetProperty(volumeProp); 129 //volume->SetProperty(volumeProp); 130 volume->SetVisibility(1); 130 volume->SetVisibility(1); 131 131 132 // add to renderer 132 // add to renderer 133 vc.fViewer->renderer->AddActor(actor); 133 vc.fViewer->renderer->AddActor(actor); 134 //vc.fViewer->renderer->AddVolume(volume); 134 //vc.fViewer->renderer->AddVolume(volume); 135 } 135 } 136 136 137 void G4VtkUnstructuredGridPipeline::PseudoScen 137 void G4VtkUnstructuredGridPipeline::PseudoSceneForTetCells::AddSolid(const G4VSolid& solid) { 138 if (fpPVModel->GetCurrentDepth() == fDepth) 138 if (fpPVModel->GetCurrentDepth() == fDepth) { // Leaf-level cells only 139 // Need to know it's a tet !!!! or impleme 139 // Need to know it's a tet !!!! or implement G4VSceneHandler::AddSolid (const G4Tet&) !!!! 140 try { 140 try { 141 const G4Tet& tet = dynamic_cast<const G4 141 const G4Tet& tet = dynamic_cast<const G4Tet&>(solid); 142 const auto* pVisAtts = fpPVModel->GetCur 142 const auto* pVisAtts = fpPVModel->GetCurrentLV()->GetVisAttributes(); 143 const auto& colour = pVisAtts->GetColour 143 const auto& colour = pVisAtts->GetColour(); 144 144 145 G4ThreeVector p0, p1, p2, p3; 145 G4ThreeVector p0, p1, p2, p3; 146 tet.GetVertices(p0, p1, p2, p3); 146 tet.GetVertices(p0, p1, p2, p3); 147 147 148 G4int idx[4] = {-1, -1, -1, -1}; 148 G4int idx[4] = {-1, -1, -1, -1}; 149 149 150 if (iCell > 0 && iCell < 100000000) { 150 if (iCell > 0 && iCell < 100000000) { 151 #if 0 151 #if 0 152 G4ThreeVector pts[4] = {p0, p1, p2, p3 152 G4ThreeVector pts[4] = {p0, p1, p2, p3}; 153 for (G4int i = 0; i < 4; i++) { 153 for (G4int i = 0; i < 4; i++) { 154 auto h = MakeHash(pts[i]); 154 auto h = MakeHash(pts[i]); 155 if (fpPointMap.find(h) == fpPointMap 155 if (fpPointMap.find(h) == fpPointMap.end()) { 156 fpPointMap.insert(std::make_pair(h 156 fpPointMap.insert(std::make_pair(h, iPoint)); 157 fpPoints->InsertNextPoint(pts[i].x 157 fpPoints->InsertNextPoint(pts[i].x(), pts[i].y(), pts[i].z()); 158 fpPointVector.push_back(pts[i]); 158 fpPointVector.push_back(pts[i]); 159 idx[i] = iPoint; 159 idx[i] = iPoint; 160 iPoint++; 160 iPoint++; 161 } 161 } 162 else { 162 else { 163 idx[i] = fpPointMap[h]; 163 idx[i] = fpPointMap[h]; 164 } 164 } 165 } 165 } 166 #endif 166 #endif 167 167 168 #if 1 168 #if 1 169 fpPoints->InsertNextPoint(p0.x(), p0.y 169 fpPoints->InsertNextPoint(p0.x(), p0.y(), p0.z()); 170 fpPoints->InsertNextPoint(p1.x(), p1.y 170 fpPoints->InsertNextPoint(p1.x(), p1.y(), p1.z()); 171 fpPoints->InsertNextPoint(p2.x(), p2.y 171 fpPoints->InsertNextPoint(p2.x(), p2.y(), p2.z()); 172 fpPoints->InsertNextPoint(p3.x(), p3.y 172 fpPoints->InsertNextPoint(p3.x(), p3.y(), p3.z()); 173 iPoint += 4; 173 iPoint += 4; 174 174 175 idx[0] = 4 * iCellAdd; 175 idx[0] = 4 * iCellAdd; 176 idx[1] = 4 * iCellAdd + 1; 176 idx[1] = 4 * iCellAdd + 1; 177 idx[2] = 4 * iCellAdd + 2; 177 idx[2] = 4 * iCellAdd + 2; 178 idx[3] = 4 * iCellAdd + 3; 178 idx[3] = 4 * iCellAdd + 3; 179 #endif 179 #endif 180 180 181 vtkNew<vtkTetra> tetra; 181 vtkNew<vtkTetra> tetra; 182 tetra->GetPointIds()->SetId(0, idx[0]) 182 tetra->GetPointIds()->SetId(0, idx[0]); 183 tetra->GetPointIds()->SetId(1, idx[1]) 183 tetra->GetPointIds()->SetId(1, idx[1]); 184 tetra->GetPointIds()->SetId(2, idx[2]) 184 tetra->GetPointIds()->SetId(2, idx[2]); 185 tetra->GetPointIds()->SetId(3, idx[3]) 185 tetra->GetPointIds()->SetId(3, idx[3]); 186 186 187 fpGrid->InsertNextCell(tetra->GetCellT 187 fpGrid->InsertNextCell(tetra->GetCellType(), tetra->GetPointIds()); 188 188 189 auto hash = MakeHash(colour); 189 auto hash = MakeHash(colour); 190 unsigned short int iColour = 0; 190 unsigned short int iColour = 0; 191 if (fpColourMap.find(hash) == fpColour 191 if (fpColourMap.find(hash) == fpColourMap.end()) { 192 iColour = fpColourMap.size(); 192 iColour = fpColourMap.size(); 193 fpColourMap.insert(std::make_pair(ha 193 fpColourMap.insert(std::make_pair(hash,iColour)); 194 double c[4] = {colour.GetRed(), colo 194 double c[4] = {colour.GetRed(), colour.GetGreen(), colour.GetBlue(), colour.GetAlpha()}; 195 fpColourLUT->SetIndexedColorRGBA(iCo 195 fpColourLUT->SetIndexedColorRGBA(iColour,c); 196 } 196 } 197 else { 197 else { 198 iColour = fpColourMap[hash]; 198 iColour = fpColourMap[hash]; 199 } 199 } 200 200 201 fpCellColourValues->InsertNextTuple4(c 201 fpCellColourValues->InsertNextTuple4(colour.GetRed(), colour.GetGreen(), colour.GetBlue(), colour.GetAlpha()); 202 fpCellColourIndices->InsertNextTuple1( 202 fpCellColourIndices->InsertNextTuple1(iColour); 203 203 204 iCellAdd++; 204 iCellAdd++; 205 } 205 } 206 iCell++; 206 iCell++; 207 } 207 } 208 catch (const std::bad_cast&) { 208 catch (const std::bad_cast&) { 209 G4ExceptionDescription ed; 209 G4ExceptionDescription ed; 210 ed << "Called for a mesh that is not a t 210 ed << "Called for a mesh that is not a tetrahedron mesh: " << solid.GetName(); 211 G4Exception("PseudoSceneForTetVertices", 211 G4Exception("PseudoSceneForTetVertices","visman0108",JustWarning,ed); 212 } 212 } 213 } 213 } 214 } 214 } 215 215 216 void G4VtkUnstructuredGridPipeline::PseudoScen 216 void G4VtkUnstructuredGridPipeline::PseudoSceneForCubicalCells::AddSolid(const G4Box& box) { 217 if (fpPVModel->GetCurrentDepth() == fDepth) 217 if (fpPVModel->GetCurrentDepth() == fDepth) { // Leaf-level cells only 218 try { 218 try { 219 const auto* pVisAtts = fpPVModel->GetCur 219 const auto* pVisAtts = fpPVModel->GetCurrentLV()->GetVisAttributes(); 220 const auto& colour = pVisAtts->GetColour 220 const auto& colour = pVisAtts->GetColour(); 221 const G4ThreeVector& position = fpCurren 221 const G4ThreeVector& position = fpCurrentObjectTransformation->getTranslation(); 222 222 223 fpPoints->InsertNextPoint(position.x()-b 223 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()-box.GetZHalfLength()); 224 fpPoints->InsertNextPoint(position.x()-b 224 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()-box.GetZHalfLength()); 225 fpPoints->InsertNextPoint(position.x()+b 225 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()-box.GetZHalfLength()); 226 fpPoints->InsertNextPoint(position.x()+b 226 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()-box.GetZHalfLength()); 227 227 228 fpPoints->InsertNextPoint(position.x()-b 228 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()+box.GetZHalfLength()); 229 fpPoints->InsertNextPoint(position.x()-b 229 fpPoints->InsertNextPoint(position.x()-box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()+box.GetZHalfLength()); 230 fpPoints->InsertNextPoint(position.x()+b 230 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()-box.GetYHalfLength(),position.z()+box.GetZHalfLength()); 231 fpPoints->InsertNextPoint(position.x()+b 231 fpPoints->InsertNextPoint(position.x()+box.GetXHalfLength(),position.y()+box.GetYHalfLength(),position.z()+box.GetZHalfLength()); 232 232 233 vtkNew<vtkVoxel> voxel; 233 vtkNew<vtkVoxel> voxel; 234 voxel->GetPointIds()->SetId(0, 8*iCell); 234 voxel->GetPointIds()->SetId(0, 8*iCell); 235 voxel->GetPointIds()->SetId(1, 8*iCell+1 235 voxel->GetPointIds()->SetId(1, 8*iCell+1); 236 voxel->GetPointIds()->SetId(2, 8*iCell+2 236 voxel->GetPointIds()->SetId(2, 8*iCell+2); 237 voxel->GetPointIds()->SetId(3, 8*iCell+3 237 voxel->GetPointIds()->SetId(3, 8*iCell+3); 238 voxel->GetPointIds()->SetId(4, 8*iCell+4 238 voxel->GetPointIds()->SetId(4, 8*iCell+4); 239 voxel->GetPointIds()->SetId(5, 8*iCell+5 239 voxel->GetPointIds()->SetId(5, 8*iCell+5); 240 voxel->GetPointIds()->SetId(6, 8*iCell+6 240 voxel->GetPointIds()->SetId(6, 8*iCell+6); 241 voxel->GetPointIds()->SetId(7, 8*iCell+7 241 voxel->GetPointIds()->SetId(7, 8*iCell+7); 242 242 243 fpGrid->InsertNextCell(voxel->GetCellTyp 243 fpGrid->InsertNextCell(voxel->GetCellType(), voxel->GetPointIds()); 244 244 245 double cols[] = {colour.GetRed(), 245 double cols[] = {colour.GetRed(), 246 colour.GetGreen(), 246 colour.GetGreen(), 247 colour.GetBlue(), 247 colour.GetBlue(), 248 colour.GetAlpha()}; 248 colour.GetAlpha()}; 249 fpCellColourValues->InsertNextTuple(cols 249 fpCellColourValues->InsertNextTuple(cols); 250 //fpCellColourIndices->InsertNextTuple1( 250 //fpCellColourIndices->InsertNextTuple1(iColour); 251 251 252 252 253 iCell++; 253 iCell++; 254 254 255 } 255 } 256 catch (const std::bad_cast&) { 256 catch (const std::bad_cast&) { 257 G4ExceptionDescription ed; 257 G4ExceptionDescription ed; 258 ed << "Called for a mesh that is not a c 258 ed << "Called for a mesh that is not a cubical mesh: " << box.GetName(); 259 G4Exception("PseudoSceneForCubicalVertic 259 G4Exception("PseudoSceneForCubicalVertices", "visman0108", JustWarning, ed); 260 } 260 } 261 } 261 } 262 } 262 } 263 263 264 void G4VtkUnstructuredGridPipeline::SetUnstruc 264 void G4VtkUnstructuredGridPipeline::SetUnstructuredGridData(const G4Mesh &mesh) { 265 265 266 // Modelling parameters 266 // Modelling parameters 267 G4ModelingParameters tmpMP; 267 G4ModelingParameters tmpMP; 268 tmpMP.SetCulling(true); // This 268 tmpMP.SetCulling(true); // This avoids drawing transparent... 269 tmpMP.SetCullingInvisible(true); // ... o 269 tmpMP.SetCullingInvisible(true); // ... or invisble volumes. 270 270 271 // Physical volume model 271 // Physical volume model 272 const auto& container = mesh.GetContainerVol 272 const auto& container = mesh.GetContainerVolume(); 273 const G4bool useFullExtent = true; // To av 273 const G4bool useFullExtent = true; // To avoid calculating the extent 274 G4PhysicalVolumeModel tmpPVModel(container, 274 G4PhysicalVolumeModel tmpPVModel(container, 275 G4PhysicalV 275 G4PhysicalVolumeModel::UNLIMITED, 276 G4Transform 276 G4Transform3D(), // so that positions are in local coordinates 277 &tmpMP, 277 &tmpMP, 278 useFullExte 278 useFullExtent); 279 279 280 // Graphics scene 280 // Graphics scene 281 if(mesh.GetMeshType() == G4Mesh::tetrahedron 281 if(mesh.GetMeshType() == G4Mesh::tetrahedron) { 282 PseudoSceneForTetCells pseudoScene(&tmpPVM 282 PseudoSceneForTetCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points, 283 pointCo 283 pointColourValues, cellColourValues, 284 pointCo 284 pointColourIndices, cellColourIndices, 285 colourL 285 colourLUT, unstructuredGrid); 286 tmpPVModel.DescribeYourselfTo(pseudoScene) 286 tmpPVModel.DescribeYourselfTo(pseudoScene); 287 //G4cout << pseudoScene.GetNumberOfPoints( 287 //G4cout << pseudoScene.GetNumberOfPoints() << " " << pseudoScene.GetNumberOfCells() << " " << pseudoScene.GetNumberOfAddedCells() << G4endl; 288 //pseudoScene.DumpColourMap(); 288 //pseudoScene.DumpColourMap(); 289 } 289 } 290 else if(mesh.GetMeshType() == G4Mesh::nested 290 else if(mesh.GetMeshType() == G4Mesh::nested3DRectangular) { 291 PseudoSceneForCubicalCells pseudoScene(&tm 291 PseudoSceneForCubicalCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points, 292 poi 292 pointColourValues, cellColourValues, 293 poi 293 pointColourIndices, cellColourIndices, 294 col 294 colourLUT, unstructuredGrid); 295 tmpPVModel.DescribeYourselfTo(pseudoScene) 295 tmpPVModel.DescribeYourselfTo(pseudoScene); 296 } 296 } 297 else if(mesh.GetMeshType() == G4Mesh::rectan 297 else if(mesh.GetMeshType() == G4Mesh::rectangle) { 298 PseudoSceneForCubicalCells pseudoScene(&tm 298 PseudoSceneForCubicalCells pseudoScene(&tmpPVModel, mesh.GetMeshDepth(), points, 299 poi 299 pointColourValues, cellColourValues, 300 poi 300 pointColourIndices, cellColourIndices, 301 col 301 colourLUT, unstructuredGrid); 302 tmpPVModel.DescribeYourselfTo(pseudoScene) 302 tmpPVModel.DescribeYourselfTo(pseudoScene); 303 } 303 } 304 304 305 } 305 }