Struct vtkio::model::Vtk[][src]

pub struct Vtk {
    pub version: Version,
    pub title: String,
    pub byte_order: ByteOrder,
    pub data: DataSet,
    pub file_path: Option<PathBuf>,
}

Model of the VTK file.

This type unifies legacy and XML data representations.

Fields

version: Versiontitle: Stringbyte_order: ByteOrderdata: DataSetfile_path: Option<PathBuf>

The path to the source file of this Vtk file (if any).

This is used to load pieces stored in other files used in “Parallel” XML file types.

Implementations

impl Vtk[src]

pub fn load_all_pieces(&mut self) -> Result<(), Error>[src]

Loads all referenced pieces into the current struct.

This function is useful for “Parallel” XML files like .pvtu, .pvtp, etc. For all other files this is a no-op.

impl Vtk[src]

pub fn try_into_xml_format(
    self,
    compressor: Compressor,
    compression_level: u32
) -> Result<VTKFile, Error>
[src]

Converts the given Vtk model into an XML format represented by VTKFile.

This function allows one to specify the compression level (0-9):

0 -> No compression
1 -> Fastest write
...
5 -> Balanced performance
...
9 -> Slowest but smallest file size.

impl Vtk[src]

pub fn parse_legacy_be(reader: impl Read) -> Result<Vtk, Error>[src]

Parse a legacy VTK file from the given reader.

If the file is in binary format, numeric types will be interpreted in big endian format, which is the most common among VTK files. Note that this function and parse_legacy_le also work equally well for parsing VTK files in ASCII format.

Examples

Parsing an ASCII file:

use vtkio::model::*; // import the model definition of a VTK file
let vtk_ascii: &[u8] = b"
# vtk DataFile Version 2.0
Triangle example
ASCII
DATASET POLYDATA
POINTS 3 float
0.0 0.0 0.0
1.0 0.0 0.0
0.0 0.0 -1.0

POLYGONS 1 4
3 0 1 2
";

let vtk = Vtk::parse_legacy_be(vtk_ascii).expect("Failed to parse vtk file");

assert_eq!(vtk, Vtk {
    version: Version::new((2,0)),
    byte_order: ByteOrder::BigEndian,
    title: String::from("Triangle example"),
    file_path: None,
    data: DataSet::inline(PolyDataPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0].into(),
        polys: Some(VertexNumbers::Legacy {
            num_cells: 1,
            vertices: vec![3, 0, 1, 2]
        }),
        data: Attributes::new(),
        ..Default::default()
    })
});

pub fn parse_legacy_le(reader: impl Read) -> Result<Vtk, Error>[src]

Parse a legacy VTK file from the given reader.

If the file is in binary format, numeric types will be interpreted in little endian format. Note that this function and parse_legacy_be also work equally well for parsing VTK files in ASCII format.

Examples

Parsing an ASCII file:

use vtkio::model::*; // import the model definition of a VTK file
let vtk_ascii: &[u8] = b"
# vtk DataFile Version 2.0
Triangle example
ASCII
DATASET POLYDATA
POINTS 3 float
0.0 0.0 0.0
1.0 0.0 0.0
0.0 0.0 -1.0

POLYGONS 1 4
3 0 1 2
";

let vtk = Vtk::parse_legacy_le(vtk_ascii).expect("Failed to parse vtk file");

assert_eq!(vtk, Vtk {
    version: Version::new((2,0)),
    byte_order: ByteOrder::LittleEndian,
    title: String::from("Triangle example"),
    file_path: None,
    data: DataSet::inline(PolyDataPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0].into(),
        polys: Some(VertexNumbers::Legacy {
            num_cells: 1,
            vertices: vec![3, 0, 1, 2]
        }),
        data: Attributes::new(),
        ..Default::default()
    })
});

pub fn parse_legacy_buf_be(
    reader: impl Read,
    buf: &mut Vec<u8>
) -> Result<Vtk, Error>
[src]

Parse a legacy VTK file in big endian format from the given reader and a buffer.

This is the buffered version of Vtk::parse_legacy_be, which allows one to reuse the same heap allocated space when reading many files.

pub fn parse_legacy_buf_le(
    reader: impl Read,
    buf: &mut Vec<u8>
) -> Result<Vtk, Error>
[src]

Parse a legacy VTK file in little endian format from the given reader and a buffer.

This is the buffered version of parse_legacy_le, which allows one to reuse the same heap allocated space when reading many files.

pub fn parse_xml(reader: impl BufRead) -> Result<Vtk, Error>[src]

Parse a modern XML style VTK file from a given reader.

Examples

Parsing a binary file in big endian format representing a polygon mesh consisting of a single triangle:

use vtkio::model::*; // import the model definition of a VTK file

let input: &[u8] = b"\
<VTKFile type=\"PolyData\" version=\"2.0\" byte_order=\"BigEndian\">\
  <PolyData>\
    <Piece NumberOfPoints=\"3\" NumberOfLines=\"0\" NumberOfStrips=\"0\" NumberOfPolys=\"1\" NumberOfVerts=\"0\">\
      <PointData/>\
      <CellData/>\
      <Points>\
        <DataArray type=\"Float32\" format=\"binary\" NumberOfComponents=\"3\">\
          AAAAAAAAAAQAAAAAAAAAAAAAAAA/gAAAAAAAAAAAAAAAAAAAAAAAAL+AAAA=\
        </DataArray>\
      </Points>\
      <Polys>\
        <DataArray type=\"UInt64\" Name=\"connectivity\" format=\"binary\" NumberOfComponents=\"1\">\
          AAAAAAAAAAgAAAAAAAAAAAAAAAAAAAABAAAAAAAAAAI=\
        </DataArray>\
        <DataArray type=\"UInt64\" Name=\"offsets\" format=\"binary\" NumberOfComponents=\"1\">\
          AAAAAAAAAAgAAAAAAAAAAw==\
        </DataArray>\
      </Polys>\
    </Piece>\
  </PolyData>\
</VTKFile>";

let vtk = Vtk::parse_xml(input).expect("Failed to parse XML VTK file");

assert_eq!(vtk, Vtk {
    version: Version::new((2,0)),
    byte_order: ByteOrder::BigEndian, // This is default
    title: String::new(),
    file_path: None,
    data: DataSet::inline(PolyDataPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0].into(),
        polys: Some(VertexNumbers::XML {
            connectivity: vec![0, 1, 2],
            offsets: vec![3]
        }),
        data: Attributes::new(),
        ..Default::default()
    })
});

pub fn import(file_path: impl AsRef<Path>) -> Result<Vtk, Error>[src]

Import a VTK file at the specified path.

This function determines the vtk file type from the extension as prescribed by the VTK file formats documentation:

  • Legacy (.vtk) – Simple legacy file format (Non-XML)
  • Image data (.vti) – Serial vtkImageData (structured)
  • PolyData (.vtp) – Serial vtkPolyData (unstructured)
  • RectilinearGrid (.vtr) – Serial vtkRectilinearGrid (structured)
  • StructuredGrid (.vts) – Serial vtkStructuredGrid (structured)
  • UnstructuredGrid (.vtu) – Serial vtkUnstructuredGrid (unstructured)
  • PImageData (.pvti) – Parallel vtkImageData (structured)
  • PPolyData (.pvtp) – Parallel vtkPolyData (unstructured)
  • PRectilinearGrid (.pvtr) – Parallel vtkRectilinearGrid (structured)
  • PStructuredGrid (.pvts) – Parallel vtkStructuredGrid (structured)
  • PUnstructuredGrid (.pvtu) – Parallel vtkUnstructuredGrid (unstructured)

Examples

The following example imports a legacy .vtk file called tet.vtk, and panics with an appropriate error message if the file fails to load.

use vtkio::Vtk;
use std::path::PathBuf;

let file_path = PathBuf::from("tet.vtk");

let mut vtk_file = Vtk::import(&file_path)
    .expect(&format!("Failed to load file: {:?}", file_path));

pub fn import_legacy_le(file_path: impl AsRef<Path>) -> Result<Vtk, Error>[src]

Import a legacy VTK file at the specified path.

If the file is in binary format, numeric types will be interpreted in little endian format. For the default byte order used by most .vtk files use import or import_legacy_be.

pub fn import_le(file_path: impl AsRef<Path>) -> Result<Vtk, Error>[src]

👎 Deprecated since 0.6.2:

Please use Vtk::import_legacy_le instead

pub fn import_legacy_be(file_path: impl AsRef<Path>) -> Result<Vtk, Error>[src]

Import a legacy VTK file at the specified path.

If the file is in binary format, numeric types will be interpreted in big endian format. This function behaves the same as import, but expects the given file to be strictly in legacy .vtk format.

pub fn import_be(file_path: impl AsRef<Path>) -> Result<Vtk, Error>[src]

👎 Deprecated since 0.6.2:

Please use Vtk::import_legacy_be instead

pub fn export(self, file_path: impl AsRef<Path>) -> Result<(), Error>[src]

Export given Vtk file to the specified file.

The type of file exported is determined by the extension in file_path.

Files ending in .vtk are exported in binary format. For exporting in ASCII, use export_ascii.

Endianness is determined by the byte_order field of the Vtk type.

Examples

use vtkio::model::*;
use std::path::PathBuf;
let vtk = Vtk {
    version: Version::new((4,1)),
    byte_order: ByteOrder::BigEndian,
    title: String::from("Tetrahedron"),
    file_path: Some(PathBuf::from("./test.vtk")),
    data: DataSet::inline(UnstructuredGridPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 1.0, 0.0].into(),
        cells: Cells {
            cell_verts: VertexNumbers::Legacy {
                num_cells: 1,
                vertices: vec![4, 0, 1, 2, 3]
            },
            types: vec![CellType::Tetra],
        },
        data: Attributes::new(),
    })
};
vtk.export("test.vtk");

pub fn write_legacy(self, writer: impl Write) -> Result<(), Error>[src]

Write the given VTK file in binary legacy format to the specified Writer.

Examples

Writing a binary file in big endian format representing a polygon mesh consisting of a single triangle:

use vtkio::model::*; // import model definition of a VTK file

let mut vtk_bytes = Vec::<u8>::new();

Vtk {
    version: Version::new((2,0)),
    byte_order: ByteOrder::BigEndian,
    title: String::from("Triangle example"),
    file_path: None,
    data: DataSet::inline(PolyDataPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0].into(),
        polys: Some(VertexNumbers::Legacy {
            num_cells: 1,
            vertices: vec![3, 0, 1, 2]
        }),
        data: Attributes::new(),
        ..Default::default()
    })
}.write_legacy(&mut vtk_bytes);

println!("{}", String::from_utf8_lossy(&vtk_bytes));

pub fn write_legacy_ascii(self, writer: impl Write) -> Result<(), Error>[src]

Write the given VTK file in binary legacy format to the specified Writer.

Examples

Writing an ASCII file representing a polygon mesh consisting of a single triangle:

use vtkio::model::*; // import model definition of a VTK file

let mut vtk_string = String::new();

Vtk {
    version: Version::new((2,0)),
    byte_order: ByteOrder::BigEndian, // Ignored
    title: String::from("Triangle example"),
    file_path: None,
    data: DataSet::inline(PolyDataPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0].into(),
        polys: Some(VertexNumbers::Legacy {
            num_cells: 1,
            vertices: vec![3, 0, 1, 2]
        }),
        data: Attributes::new(),
        ..Default::default()
    })
}.write_legacy_ascii(&mut vtk_string);

assert_eq!(vtk_string.as_str(), "\
# vtk DataFile Version 2.0
Triangle example
ASCII

DATASET POLYDATA
POINTS 3 float
0 0 0 1 0 0 0 0 -1

POLYGONS 1 4
3 0 1 2

POINT_DATA 3

CELL_DATA 1

");

pub fn write_xml(self, writer: impl Write) -> Result<(), Error>[src]

Write the given VTK file in modern XML format to the specified Writer.

Examples

Writing a binary file in big endian format representing a polygon mesh consisting of a single triangle:

use vtkio::model::*; // import model definition of a VTK file

let mut vtk_bytes = Vec::<u8>::new();

Vtk {
    version: Version::new((2,0)),
    byte_order: ByteOrder::BigEndian,
    title: String::from("Triangle example"),
    file_path: None,
    data: DataSet::inline(PolyDataPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0].into(),
        polys: Some(VertexNumbers::Legacy {
            num_cells: 1,
            vertices: vec![3, 0, 1, 2]
        }),
        data: Attributes::new(),
        ..Default::default()
    })
}.write_xml(&mut vtk_bytes);

assert_eq!(String::from_utf8_lossy(&vtk_bytes), "\
<VTKFile type=\"PolyData\" version=\"2.0\" byte_order=\"BigEndian\" header_type=\"UInt64\">\
  <PolyData>\
    <Piece NumberOfPoints=\"3\" NumberOfLines=\"0\" NumberOfStrips=\"0\" NumberOfPolys=\"1\" NumberOfVerts=\"0\">\
      <PointData/>\
      <CellData/>\
      <Points>\
        <DataArray type=\"Float32\" format=\"binary\" NumberOfComponents=\"3\">\
          AAAAAAAAACQAAAAAAAAAAAAAAAA/gAAAAAAAAAAAAAAAAAAAAAAAAL+AAAA=\
        </DataArray>\
      </Points>\
      <Polys>\
        <DataArray type=\"UInt64\" Name=\"connectivity\" format=\"binary\" NumberOfComponents=\"1\">\
          AAAAAAAAABgAAAAAAAAAAAAAAAAAAAABAAAAAAAAAAI=\
        </DataArray>\
        <DataArray type=\"UInt64\" Name=\"offsets\" format=\"binary\" NumberOfComponents=\"1\">\
          AAAAAAAAAAgAAAAAAAAAAw==\
        </DataArray>\
      </Polys>\
    </Piece>\
  </PolyData>\
</VTKFile>");

pub fn export_le(self, file_path: impl AsRef<Path>) -> Result<(), Error>[src]

Export the VTK data to the specified path in little endian binary format.

This function is used as export but overrides endiannes.

pub fn export_be(self, file_path: impl AsRef<Path>) -> Result<(), Error>[src]

Export the VTK data to the specified path in big endian binary format.

This function is used as export but overrides endiannes.

pub fn export_ascii(self, file_path: impl AsRef<Path>) -> Result<(), Error>[src]

Export VTK data to the specified file in ASCII format.

Examples

use vtkio::model::*;
use std::path::PathBuf;
let vtk = Vtk {
    version: Version::new((4,1)),
    title: String::from("Tetrahedron"),
    byte_order: ByteOrder::BigEndian,
    file_path: Some(PathBuf::from("./test.vtk")),
    data: DataSet::inline(UnstructuredGridPiece {
        points: vec![0.0f32, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 1.0, 0.0].into(),
        cells: Cells {
            cell_verts: VertexNumbers::Legacy {
                num_cells: 1,
                vertices: vec![4, 0, 1, 2, 3]
            },
            types: vec![CellType::Tetra],
        },
        data: Attributes::new(),
    })
};
vtk.export_ascii("test.vtk");

Trait Implementations

impl Clone for Vtk[src]

impl Debug for Vtk[src]

impl PartialEq<Vtk> for Vtk[src]

impl StructuralPartialEq for Vtk[src]

impl TryFrom<VTKFile> for Vtk[src]

type Error = Error

The type returned in the event of a conversion error.

impl TryFrom<Vtk> for VTKFile[src]

type Error = Error

The type returned in the event of a conversion error.

Auto Trait Implementations

impl RefUnwindSafe for Vtk

impl Send for Vtk

impl Sync for Vtk

impl Unpin for Vtk

impl UnwindSafe for Vtk

Blanket Implementations

impl<T> Any for T where
    T: 'static + ?Sized
[src]

impl<T> Borrow<T> for T where
    T: ?Sized
[src]

impl<T> BorrowMut<T> for T where
    T: ?Sized
[src]

impl<T> From<T> for T[src]

impl<T, U> Into<U> for T where
    U: From<T>, 
[src]

impl<T> ToOwned for T where
    T: Clone
[src]

type Owned = T

The resulting type after obtaining ownership.

impl<T, U> TryFrom<U> for T where
    U: Into<T>, 
[src]

type Error = Infallible

The type returned in the event of a conversion error.

impl<T, U> TryInto<U> for T where
    U: TryFrom<T>, 
[src]

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.