iMSTK
Interactive Medical Simulation Toolkit
GeometryProcessingExample.cpp
1 /*
2 ** This file is part of the Interactive Medical Simulation Toolkit (iMSTK)
3 ** iMSTK is distributed under the Apache License, Version 2.0.
4 ** See accompanying NOTICE for details.
5 */
6 
7 #include "imstkCamera.h"
8 #include "imstkDataArray.h"
9 #include "imstkDirectionalLight.h"
10 #include "imstkImageData.h"
11 #include "imstkMeshIO.h"
12 #include "imstkNew.h"
13 #include "imstkQuadricDecimate.h"
14 #include "imstkRenderMaterial.h"
15 #include "imstkScene.h"
16 #include "imstkSceneManager.h"
17 #include "imstkSceneObject.h"
18 #include "imstkSimulationManager.h"
19 #include "imstkSurfaceMesh.h"
20 #include "imstkSurfaceMeshDistanceTransform.h"
21 #include "imstkSurfaceMeshFlyingEdges.h"
22 #include "imstkTetrahedralMesh.h"
23 #include "imstkVisualModel.h"
24 #include "imstkVTKViewer.h"
25 
26 using namespace imstk;
27 
31 int
32 main()
33 {
35 
36  // simManager and Scene
37  imstkNew<Scene> scene("GeometryProcessing");
38  scene->getActiveCamera()->setPosition(Vec3d(0.0, 12.0, 12.0));
39 
40  auto coarseTetMesh = MeshIO::read<TetrahedralMesh>(iMSTK_DATA_ROOT "/asianDragon/asianDragon.veg");
41  std::shared_ptr<SurfaceMesh> coarseSurfMesh = coarseTetMesh->extractSurfaceMesh();
42 
43  // Compute DT
45  createSdf->setInputMesh(coarseSurfMesh);
46  createSdf->setDimensions(50, 50, 50);
47  createSdf->update();
48 
49  // Erode
50  const double erosionDist = 0.2;
51  DataArray<double>& scalars = *std::dynamic_pointer_cast<DataArray<double>>(createSdf->getOutputImage()->getScalars());
52  for (int i = 0; i < scalars.size(); i++)
53  {
54  scalars[i] += erosionDist;
55  }
56 
57  // Extract surface
59  isoExtract->setInputImage(createSdf->getOutputImage());
60  isoExtract->update();
61 
62  // Reduce surface
64  reduce->setInputMesh(isoExtract->getOutputMesh());
65  reduce->setTargetReduction(0.5);
66  reduce->update();
67 
68  // Create the scene object
69  imstkNew<SceneObject> sceneObj("Mesh");
70  // Create the eroded visual model
71  {
72  imstkNew<VisualModel> surfMeshModel;
73  surfMeshModel->setGeometry(reduce->getOutput());
74  imstkNew<RenderMaterial> material;
75  material->setDisplayMode(RenderMaterial::DisplayMode::Surface);
76  material->setLineWidth(4.0);
77  material->setEdgeColor(Color::Color::Orange);
78  surfMeshModel->setRenderMaterial(material);
79  sceneObj->addVisualModel(surfMeshModel);
80  }
81  // Create the original mesh visual model
82  {
83  imstkNew<VisualModel> surfMeshModel;
84  surfMeshModel->setGeometry(coarseSurfMesh);
85  imstkNew<RenderMaterial> material;
86  material->setColor(Color::Red);
87  material->setDisplayMode(RenderMaterial::DisplayMode::Surface);
88  material->setLineWidth(1.0);
89  material->setOpacity(0.2f);
90  surfMeshModel->setRenderMaterial(material);
91  sceneObj->addVisualModel(surfMeshModel);
92  }
93  scene->addSceneObject(sceneObj);
94 
95  // Light
97  light->setFocalPoint(Vec3d(5.0, -8.0, -5.0));
98  light->setIntensity(1);
99  scene->addLight("light", light);
100 
101  // Run the simulation
102  {
103  // Setup a viewer to render in its own thread
104  imstkNew<VTKViewer> viewer;
105  viewer->setActiveScene(scene);
106 
107  // Setup a scene manager to advance the scene in its own thread
108  imstkNew<SceneManager> sceneManager;
109  sceneManager->setActiveScene(scene);
110 
112  driver->addModule(viewer);
113  driver->addModule(sceneManager);
114  driver->start();
115  }
116 
117  return 0;
118 }
void setActiveScene(std::shared_ptr< Scene > scene) override
Set scene to be rendered.
int size() const
Get number of values/tuples.
void setInputMesh(std::shared_ptr< SurfaceMesh > inputMesh)
Required input, port 0.
Compound Geometry.
void addModule(std::shared_ptr< Module > module) override
Add a module to run.
void setIntensity(const double intensity)
Set the light intensity. This value is unbounded.
Definition: imstkLight.h:71
std::shared_ptr<T> obj = std::make_shared<T>(); equivalent, convenience class for STL shared allocati...
Definition: imstkNew.h:29
std::shared_ptr< Geometry > getOutput(size_t port=0) const
Returns output geometry given port, returns nullptr if doesn&#39;t exist.
void update()
Do the actual algorithm.
void setInputImage(std::shared_ptr< ImageData > inputImage)
Required input, port 0.
void setRenderMaterial(std::shared_ptr< RenderMaterial > renderMaterial)
Set/Get render material.
void setFocalPoint(const Vec3d &p)
Get/Set the light focal point.
Definition: imstkLight.h:33
void setActiveScene(std::string newSceneName)
Sets the currently updating scene.
static LoggerG3 & startLogger()
Starts logger with default sinks, use getInstance to create a logger with no sinks.
void setInputMesh(std::shared_ptr< SurfaceMesh > surfMesh)
Required input, port 0.