svFSIplus
Classes | Enumerations | Functions | Variables
vtk_xml_parser Namespace Reference

Classes

class  VtkFileExtentions
 

Enumerations

enum class  VtkFileFormat {
  VTP ,
  VTU
}
 

Functions

const std::string NODE_IDS_NAME ("GlobalNodeID")
 Names of data arrays store in VTK mesh files. More...
 
const std::string ELEMENT_IDS_NAME ("GlobalElementID")
 
void store_element_conn (vtkSmartPointer< vtkPolyData > vtk_polydata, faceType &face)
 Store element connectivity into the Face object. More...
 
void store_element_conn (vtkSmartPointer< vtkPolyData > vtk_polydata, mshType &mesh)
 
void store_element_conn (vtkSmartPointer< vtkUnstructuredGrid > vtk_ugrid, mshType &mesh)
 Store VTK vtkUnstructuredGrid cell connectivity into a mesh element connectivity. More...
 
void store_element_ids (vtkSmartPointer< vtkUnstructuredGrid > vtk_ugrid, mshType &mesh)
 Store element IDs into a mshType object. More...
 
void store_element_ids (vtkSmartPointer< vtkPolyData > vtk_polydata, faceType &face)
 Store element IDs into a faceType object. More...
 
void store_nodal_coords (vtkPoints *points, faceType &face)
 Store nodal coordinates from the VTK points into a face. More...
 
void store_nodal_coords (vtkPoints *points, mshType &mesh)
 Store VTK points into a mesh nodal coordinates. More...
 
void store_nodal_ids (vtkSmartPointer< vtkUnstructuredGrid > vtk_ugrid, mshType &mesh)
 Store VTK node IDs data array into a mesh nodal IDs. More...
 
void store_nodal_ids (vtkSmartPointer< vtkPolyData > vtk_polydata, faceType &face)
 Store VTK node IDs data array into a face nodal IDs. More...
 
void store_nodal_ids (vtkSmartPointer< vtkPolyData > vtk_polydata, mshType &mesh)
 
void load_fiber_direction_vtu (const std::string &file_name, const std::string &data_name, const int idx, const int nsd, mshType &mesh)
 Read fiber direction data from a VTK VTU file and copy it into a mesh.. More...
 
void load_vtp (const std::string &file_name, faceType &face)
 Store a surface mesh read from a VTK .vtp file into a Face object. More...
 
void load_vtp (const std::string &file_name, mshType &mesh)
 Store a surface mesh read from a VTK .vtp file into a Mesh object. More...
 
void load_vtu (const std::string &file_name, mshType &mesh)
 Read a mesh from a .vtu file. More...
 

Variables

std::map< unsigned char, int > vtk_cell_to_elem
 Map used to convert VTK cell types to number of nodes per element. More...
 
std::map< unsigned char, std::vector< std::vector< int > > > vtk_cell_ordering
 

Detailed Description

Copyright (c) Stanford University, The Regents of the University of California, and others.

All Rights Reserved.

See Copyright-SimVascular.txt for additional details.

Permission is hereby granted, free of charge, to any person obtaining a copy of this software and associated documentation files (the "Software"), to deal in the Software without restriction, including without limitation the rights to use, copy, modify, merge, publish, distribute, sublicense, and/or sell copies of the Software, and to permit persons to whom the Software is furnished to do so, subject to the following conditions:

The above copyright notice and this permission notice shall be included in all copies or substantial portions of the Software.

THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.

Function Documentation

◆ load_fiber_direction_vtu()

void vtk_xml_parser::load_fiber_direction_vtu ( const std::string &  file_name,
const std::string &  data_name,
const int  idx,
const int  nsd,
mshType mesh 
)

Read fiber direction data from a VTK VTU file and copy it into a mesh..

Data is stored in the mesh for all fiber direction files.

Mesh variables set mesh.fN - Fiber orientations stored at the element level.
mesh.fN shape is Array<double>(num_fiber_files*nsd, num_elems)

Arguments: file_name - The name of the VTK VTU file storing fiber direction data data_name - The name of the VTK Cell Data Array storing fiber direction data idx - The id used to identify a fiber file (0, 1, ..., number of fiber files - 1) mesh - The mesh to store the data into.

◆ load_vtp() [1/2]

void vtk_xml_parser::load_vtp ( const std::string &  file_name,
faceType face 
)

Store a surface mesh read from a VTK .vtp file into a Face object.

◆ load_vtp() [2/2]

void vtk_xml_parser::load_vtp ( const std::string &  file_name,
mshType mesh 
)

Store a surface mesh read from a VTK .vtp file into a Mesh object.

◆ load_vtu()

void vtk_xml_parser::load_vtu ( const std::string &  file_name,
mshType mesh 
)

Read a mesh from a .vtu file.

Mesh variables set mesh.gnNo - number of nodes mesh.x - node coordinates mesh.gN - node IDs mesh.gnEl - number of elements mesh.eNoN - number of noders per element mesh.gIEN - element connectivity (num_nodes_per_elem, num_elems)

Replicates 'subroutine loadVTK(vtk,fName,istat)' defined in vtkXMLParser.f90.

◆ NODE_IDS_NAME()

const std::string vtk_xml_parser::NODE_IDS_NAME ( "GlobalNodeID"  )

Names of data arrays store in VTK mesh files.

◆ store_element_conn() [1/2]

void vtk_xml_parser::store_element_conn ( vtkSmartPointer< vtkPolyData >  vtk_polydata,
faceType face 
)

Store element connectivity into the Face object.

Face variables set face.nEl - number of elements face.eNoN - number of noders per element face.IEN - element connectivity (num_nodes_per_elem, num_elems)

Note: The face.IEN array is allocated as in the Fortran code as (face.eNoN, face.nEl).

◆ store_element_conn() [2/2]

void vtk_xml_parser::store_element_conn ( vtkSmartPointer< vtkUnstructuredGrid >  vtk_ugrid,
mshType mesh 
)

Store VTK vtkUnstructuredGrid cell connectivity into a mesh element connectivity.

Mesh variables set mesh.gnEl - number of elements mesh.eNoN - number of noders per element mesh.gIEN - element connectivity (num_nodes_per_elem, num_elems)

Note: The mesh.gIEN array is allocated as in the Fortran code as (np_elem, num_elems).

◆ store_element_ids() [1/2]

void vtk_xml_parser::store_element_ids ( vtkSmartPointer< vtkPolyData >  vtk_polydata,
faceType face 
)

Store element IDs into a faceType object.

Face data set face.gE

◆ store_element_ids() [2/2]

void vtk_xml_parser::store_element_ids ( vtkSmartPointer< vtkUnstructuredGrid >  vtk_ugrid,
mshType mesh 
)

Store element IDs into a mshType object.

Todo:
[NOTE] Are element IDs used?

◆ store_nodal_coords() [1/2]

void vtk_xml_parser::store_nodal_coords ( vtkPoints *  points,
faceType face 
)

Store nodal coordinates from the VTK points into a face.

Face data set face.nNo - number of nodes face.x - node coordinates

◆ store_nodal_coords() [2/2]

void vtk_xml_parser::store_nodal_coords ( vtkPoints *  points,
mshType mesh 
)

Store VTK points into a mesh nodal coordinates.

Mesh variables set mesh.gnNo - number of nodes mesh.x - node coordinates

◆ store_nodal_ids() [1/2]

void vtk_xml_parser::store_nodal_ids ( vtkSmartPointer< vtkPolyData >  vtk_polydata,
faceType face 
)

Store VTK node IDs data array into a face nodal IDs.

Note: This array appeats to be optional.

Face variables set face.gN - nodal IDs

◆ store_nodal_ids() [2/2]

void vtk_xml_parser::store_nodal_ids ( vtkSmartPointer< vtkUnstructuredGrid >  vtk_ugrid,
mshType mesh 
)

Store VTK node IDs data array into a mesh nodal IDs.

Note: This array appeats to be optional.

Mesh variables set mesh.gN - nodal IDs

Variable Documentation

◆ vtk_cell_ordering

std::map<unsigned char, std::vector<std::vector<int> > > vtk_xml_parser::vtk_cell_ordering
Initial value:
{
{VTK_HEXAHEDRON, {{0,3,2,1},
{4,5,6,7},
{0,1,5,4},
{1,2,6,5},
{2,3,7,6},
{3,0,4,7}}},
{VTK_LINE, {{0},
{1}}},
{VTK_QUAD, {{0,1},
{1,2},
{2,3},
{3,0}}},
{VTK_TETRA, {{0,1,2},
{0,1,3},
{1,2,3},
{2,0,3}}},
{VTK_TRIANGLE, {{0,1},
{1,2},
{2,0}}},
{VTK_WEDGE, {{0,1,2},
{3,4,5},
{0,1,4,3},
{1,2,5,4},
{2,0,3,5}}}
}

◆ vtk_cell_to_elem

std::map<unsigned char,int> vtk_xml_parser::vtk_cell_to_elem
Initial value:
{
{VTK_HEXAHEDRON, 8},
{VTK_LINE, 2},
{VTK_QUAD, 4},
{VTK_TETRA, 4},
{VTK_TRIANGLE, 3},
{VTK_WEDGE, 6}
}

Map used to convert VTK cell types to number of nodes per element.