|
iMSTK
Interactive Medical Simulation Toolkit
|
Contains a set of free functions for processing geometry also contains a set of conversion and coupling functions for VTK. More...
Enumerations | |
| enum | MeshNodeRenumberingStrategy { ReverseCuthillMckee } |
| Enumeration for reordering method. | |
Functions | |
| std::shared_ptr< PointSet > | copyToPointSet (vtkSmartPointer< vtkPointSet > vtkMesh) |
| Converts vtk polydata into a imstk point set. | |
| std::shared_ptr< SurfaceMesh > | copyToSurfaceMesh (vtkSmartPointer< vtkPolyData > vtkMesh) |
| Converts vtk polydata into a imstk surface mesh. | |
| std::shared_ptr< LineMesh > | copyToLineMesh (vtkSmartPointer< vtkPolyData > vtkMesh) |
| Converts vtk polydata into a imstk surface mesh. | |
| std::shared_ptr< AbstractCellMesh > | copyToCellMesh (vtkSmartPointer< vtkUnstructuredGrid > vtkMesh) |
| Get imstk cell mesh given vtkUnstructuredGrid as input iMSTK only supports homogenous cells, not unstructured. Uses the last cell type in array. Drops others. | |
| vtkSmartPointer< vtkPointSet > | copyToVtkPointSet (std::shared_ptr< PointSet > imstkMesh) |
| Converts imstk point set into a vtk polydata. | |
| vtkSmartPointer< vtkPolyData > | copyToVtkPolyData (std::shared_ptr< LineMesh > imstkMesh) |
| Converts imstk line mesh into a vtk polydata. | |
| vtkSmartPointer< vtkPolyData > | copyToVtkPolyData (std::shared_ptr< SurfaceMesh > imstkMesh) |
| Converts imstk surface mesh into a vtk polydata. | |
| vtkSmartPointer< vtkUnstructuredGrid > | copyToVtkUnstructuredGrid (std::shared_ptr< TetrahedralMesh > imstkMesh) |
| Converts imstk tetrahedral mesh into a vtk unstructured grid. | |
| vtkSmartPointer< vtkUnstructuredGrid > | copyToVtkUnstructuredGrid (std::shared_ptr< HexahedralMesh > imstkMesh) |
| Converts imstk hexahedral mesh into a vtk unstructured grid. | |
| std::shared_ptr< VecDataArray< double, 3 > > | copyToVecDataArray (vtkSmartPointer< vtkPoints > points) |
| Copy from vtk points to a imstk vertices array. | |
| vtkSmartPointer< vtkPoints > | copyToVtkPoints (std::shared_ptr< VecDataArray< double, 3 >> vertices) |
| Copies vertices from imstk structure to VTK one. | |
| template<size_t dim> | |
| vtkSmartPointer< vtkCellArray > | copyToVtkCellArray (std::shared_ptr< VecDataArray< int, dim >> cells) |
| Copies cells of the given dimension from imstk structure to VTK one. | |
| template<size_t dim> | |
| std::shared_ptr< VecDataArray< int, dim > > | copyToVecDataArray (vtkCellArray *vtkCells) |
| Copy from vtk cells to an imstk index array. | |
| void | copyToDataMap (vtkDataSetAttributes *pointData, std::unordered_map< std::string, std::shared_ptr< AbstractDataArray >> &dataMap) |
| Copy vtkPointData to a data map. | |
| void | copyToVtkDataAttributes (vtkDataSetAttributes *pointData, const std::unordered_map< std::string, std::shared_ptr< AbstractDataArray >> &dataMap) |
| Copy from a data map to vtkDataAttributes (used for vtkCellData and vtkPointData) warning: Components are lost and DataArray's presented as single component. | |
| std::shared_ptr< SurfaceMesh > | toUVSphereSurfaceMesh (std::shared_ptr< Sphere > sphere, const unsigned int phiDivisions, const unsigned int thetaDivisions) |
| UV sphere from imstkSphere. | |
| std::shared_ptr< SurfaceMesh > | toSurfaceMesh (std::shared_ptr< AnalyticalGeometry > geom) |
| Produces SurfaceMesh from an analytical geometry. | |
| std::shared_ptr< TetrahedralMesh > | toTetGrid (const Vec3d ¢er, const Vec3d &size, const Vec3i &divisions, const Quatd orientation=Quatd::Identity()) |
| Produces a tetrahedral grid given the OrientedBox with the given divisions. More... | |
| std::shared_ptr< SurfaceMesh > | toTriangleGrid (const Vec3d ¢er, const Vec2d &size, const Vec2i &dim, const Quatd orientation=Quatd::Identity(), const double uvScale=1.0) |
| Produces a triangle grid on a plane given the imstkPlane. More... | |
| std::shared_ptr< LineMesh > | toLineGrid (const Vec3d &start, const Vec3d &dir, const double length, const int dim) |
| Creates a set of connected lines. More... | |
| int | getOpenEdgeCount (std::shared_ptr< SurfaceMesh > surfMesh) |
| Returns the number of open edges, use to tell if manifold (==0) | |
| bool | isClosed (std::shared_ptr< SurfaceMesh > surfMesh) |
| Returns if the surface is closed or not. | |
| double | getVolume (std::shared_ptr< SurfaceMesh > surfMesh) |
| Returns volume estimate of closed SurfaceMesh. | |
| std::shared_ptr< TetrahedralMesh > | createUniformMesh (const Vec3d &aabbMin, const Vec3d &aabbMax, const int nx, const int ny, const int nz) |
| Create a tetrahedral mesh based on a uniform Cartesian mesh. More... | |
| std::shared_ptr< TetrahedralMesh > | createTetrahedralMeshCover (std::shared_ptr< SurfaceMesh > surfMesh, const int nx, const int ny, int nz) |
| Create a tetrahedral mesh cover. More... | |
| template<typename NeighborContainer > | |
| std::vector< size_t > | reorderConnectivity (const std::vector< NeighborContainer > &neighbors, const MeshNodeRenumberingStrategy &method=MeshNodeRenumberingStrategy::ReverseCuthillMckee) |
| Reorder indices in a connectivity to reduce bandwidth. More... | |
| template<typename ElemConn > | |
| std::vector< size_t > | reorderConnectivity (const std::vector< ElemConn > &conn, const size_t numVerts, const MeshNodeRenumberingStrategy &method=MeshNodeRenumberingStrategy::ReverseCuthillMckee) |
| Reorder using Reverse Cuthill-Mckee. More... | |
| vtkSmartPointer< vtkDataArray > | coupleVtkDataArray (std::shared_ptr< AbstractDataArray > imstkArray) |
| Coupling functions, these create vtk data objects that point to our data objects thus no copying is done here. | |
| vtkSmartPointer< vtkImageData > | coupleVtkImageData (std::shared_ptr< ImageData > imstkImageData) |
| vtkSmartPointer< vtkDataArray > | copyToVtkDataArray (std::shared_ptr< AbstractDataArray > imstkArray) |
| Copy functions, these copy to/from vtk data objects. | |
| std::shared_ptr< AbstractDataArray > | copyToDataArray (vtkSmartPointer< vtkDataArray > vtkArray) |
| std::shared_ptr< ImageData > | copyToImageData (vtkSmartPointer< vtkImageData > imageDataVtk) |
| vtkSmartPointer< vtkImageData > | copyToVtkImageData (std::shared_ptr< ImageData > imageData) |
Contains a set of free functions for processing geometry also contains a set of conversion and coupling functions for VTK.
| std::shared_ptr< TetrahedralMesh > imstk::GeometryUtils::createTetrahedralMeshCover | ( | std::shared_ptr< SurfaceMesh > | surfMesh, |
| const int | nx, | ||
| const int | ny, | ||
| int | nz | ||
| ) |
Create a tetrahedral mesh cover.
Definition at line 1776 of file imstkGeometryUtilities.cpp.


| std::shared_ptr< TetrahedralMesh > imstk::GeometryUtils::createUniformMesh | ( | const Vec3d & | aabbMin, |
| const Vec3d & | aabbMax, | ||
| const int | nx, | ||
| const int | ny, | ||
| const int | nz | ||
| ) |
Create a tetrahedral mesh based on a uniform Cartesian mesh.
| aabbMin | the small corner of a box |
| aabbMax | the large corner of a box |
| nx | number of elements in the x-direction |
| ny | number of elements in the y-direction |
| nz | number of elements in the z-direction |
Definition at line 1701 of file imstkGeometryUtilities.cpp.


| std::vector< size_t > imstk::GeometryUtils::reorderConnectivity | ( | const std::vector< NeighborContainer > & | neighbors, |
| const MeshNodeRenumberingStrategy & | method = MeshNodeRenumberingStrategy::ReverseCuthillMckee |
||
| ) |
Reorder indices in a connectivity to reduce bandwidth.
| [in] | neighbors | array of neighbors of each vertex; eg, neighbors[i] is an object containing all neighbors of vertex-i |
Definition at line 1960 of file imstkGeometryUtilities.cpp.
| std::vector< size_t > imstk::GeometryUtils::reorderConnectivity | ( | const std::vector< ElemConn > & | conn, |
| const size_t | numVerts, | ||
| const MeshNodeRenumberingStrategy & | method = MeshNodeRenumberingStrategy::ReverseCuthillMckee |
||
| ) |
Reorder using Reverse Cuthill-Mckee.
| [in] | conn | element-to-vertex connectivity |
| [in] | numVerts | number of vertices |
| [in] | method | reordering method; see ReorderMethod |
Definition at line 1974 of file imstkGeometryUtilities.cpp.
| std::shared_ptr< LineMesh > imstk::GeometryUtils::toLineGrid | ( | const Vec3d & | start, |
| const Vec3d & | dir, | ||
| const double | length, | ||
| const int | dim | ||
| ) |
Creates a set of connected lines.
| start | of the line mesh |
| direction | to build the lines |
| Total | length of the line mesh |
| divisions |
Definition at line 493 of file imstkGeometryUtilities.cpp.

| std::shared_ptr< TetrahedralMesh > imstk::GeometryUtils::toTetGrid | ( | const Vec3d & | center, |
| const Vec3d & | size, | ||
| const Vec3i & | divisions, | ||
| const Quatd | orientation = Quatd::Identity() |
||
| ) |
Produces a tetrahedral grid given the OrientedBox with the given divisions.
| Center | of the grid |
| Size | of the grid |
| x,y,z | divisions of the grid |
| orientation | of the grid |
Definition at line 348 of file imstkGeometryUtilities.cpp.

| std::shared_ptr< SurfaceMesh > imstk::GeometryUtils::toTriangleGrid | ( | const Vec3d & | center, |
| const Vec2d & | size, | ||
| const Vec2i & | dim, | ||
| const Quatd | orientation = Quatd::Identity(), |
||
| const double | uvScale = 1.0 |
||
| ) |
Produces a triangle grid on a plane given the imstkPlane.
| Center | of the grid plane |
| Size | of the grid |
| x,y | divisions of the grid |
| orientation | of the grid |
Definition at line 425 of file imstkGeometryUtilities.cpp.

1.8.13