Vtk

Struct Vtk 

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

Model of the VTK file.

This type unifies legacy and XML data representations.

Fields§

§version: Version§title: String§byte_order: ByteOrder§data: DataSet§file_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§

Source§

impl Vtk

Source

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

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.

Source§

impl Vtk

Source

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

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.
Source§

impl Vtk

Source

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

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()
    })
});
Source

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

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()
    })
});
Source

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

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.

Source

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

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.

Source

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

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()
    })
});
Source

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

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));
Source

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

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.

Source

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

👎Deprecated since 0.6.2: Please use Vtk::import_legacy_le instead
Source

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

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.

Source

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

👎Deprecated since 0.6.2: Please use Vtk::import_legacy_be instead
Source

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

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");
Source

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

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));
Source

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

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

");
Source

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

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>");
Source

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

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

This function is used as export but overrides endiannes.

Source

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

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

This function is used as export but overrides endiannes.

Source

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

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§

Source§

impl Clone for Vtk

Source§

fn clone(&self) -> Vtk

Returns a duplicate of the value. Read more
1.0.0 · Source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
Source§

impl Debug for Vtk

Source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result<(), Error>

Formats the value using the given formatter. Read more
Source§

impl<T: Real> MeshExtractor<T> for Vtk

Source§

fn extract_mesh(&self) -> Result<Mesh<T>, Error>

Constructs an unstructured Mesh from this VTK model.

This function will clone the given model as necessary.

Source§

fn extract_polymesh(&self) -> Result<PolyMesh<T>, Error>

Constructs a PolyMesh from this VTK model.

This function will clone the given model as necessary.

Source§

fn extract_tetmesh(&self) -> Result<TetMesh<T>, Error>

Constructs a TetMesh from this VTK model.

This function will clone the given model as necessary.

Source§

fn extract_pointcloud(&self) -> Result<PointCloud<T>, Error>

Constructs a PointCloud from this VTK model.

This function will clone the given model as necessary.

Source§

impl PartialEq for Vtk

Source§

fn eq(&self, other: &Vtk) -> bool

Tests for self and other values to be equal, and is used by ==.
1.0.0 · Source§

fn ne(&self, other: &Rhs) -> bool

Tests for !=. The default implementation is almost always sufficient, and should not be overridden without very good reason.
Source§

impl TryFrom<VTKFile> for Vtk

Source§

type Error = Error

The type returned in the event of a conversion error.
Source§

fn try_from(xml: VTKFile) -> Result<Vtk, <Vtk as TryFrom<VTKFile>>::Error>

Performs the conversion.
Source§

impl TryFrom<Vtk> for VTKFile

Source§

type Error = Error

The type returned in the event of a conversion error.
Source§

fn try_from(vtk: Vtk) -> Result<VTKFile, Error>

Performs the conversion.
Source§

impl StructuralPartialEq for Vtk

Auto Trait Implementations§

§

impl Freeze for Vtk

§

impl RefUnwindSafe for Vtk

§

impl Send for Vtk

§

impl Sync for Vtk

§

impl Unpin for Vtk

§

impl UnwindSafe for Vtk

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> Bytes for T

Source§

fn as_bytes(&self) -> &[u8]

Get a slice of bytes representing Self.
Source§

fn interpret_bytes(bytes: &[u8]) -> &Self

Panics if the size of the given bytes slice is not equal to the size of Self.
Source§

impl<T> CloneBytes for T
where T: Clone + 'static,

Source§

unsafe fn clone_bytes(src: &[MaybeUninit<u8>]) -> Box<[MaybeUninit<u8>]>

Source§

unsafe fn clone_from_bytes(dst: &mut [MaybeUninit<u8>], src: &[MaybeUninit<u8>])

Source§

unsafe fn clone_into_raw_bytes( src: &[MaybeUninit<u8>], dst: &mut [MaybeUninit<u8>], )

Source§

impl<T> CloneToUninit for T
where T: Clone,

Source§

unsafe fn clone_to_uninit(&self, dest: *mut u8)

🔬This is a nightly-only experimental API. (clone_to_uninit)
Performs copy-assignment from self to dest. Read more
Source§

impl<T> DebugBytes for T
where T: Debug + 'static,

Source§

unsafe fn fmt_bytes( bytes: &[MaybeUninit<u8>], f: &mut Formatter<'_>, ) -> Result<(), Error>

Source§

impl<T> Downcast for T
where T: Any,

Source§

fn into_any(self: Box<T>) -> Box<dyn Any>

Converts Box<dyn Trait> (where Trait: Downcast) to Box<dyn Any>, which can then be downcast into Box<dyn ConcreteType> where ConcreteType implements Trait.
Source§

fn into_any_rc(self: Rc<T>) -> Rc<dyn Any>

Converts Rc<Trait> (where Trait: Downcast) to Rc<Any>, which can then be further downcast into Rc<ConcreteType> where ConcreteType implements Trait.
Source§

fn as_any(&self) -> &(dyn Any + 'static)

Converts &Trait (where Trait: Downcast) to &Any. This is needed since Rust cannot generate &Any’s vtable from &Trait’s.
Source§

fn as_any_mut(&mut self) -> &mut (dyn Any + 'static)

Converts &mut Trait (where Trait: Downcast) to &Any. This is needed since Rust cannot generate &mut Any’s vtable from &mut Trait’s.
Source§

impl<T> DowncastSend for T
where T: Any + Send,

Source§

fn into_any_send(self: Box<T>) -> Box<dyn Any + Send>

Converts Box<Trait> (where Trait: DowncastSend) to Box<dyn Any + Send>, which can then be downcast into Box<ConcreteType> where ConcreteType implements Trait.
Source§

impl<T> DowncastSync for T
where T: Any + Send + Sync,

Source§

fn into_any_sync(self: Box<T>) -> Box<dyn Any + Sync + Send>

Converts Box<Trait> (where Trait: DowncastSync) to Box<dyn Any + Send + Sync>, which can then be downcast into Box<ConcreteType> where ConcreteType implements Trait.
Source§

fn into_any_arc(self: Arc<T>) -> Arc<dyn Any + Sync + Send>

Converts Arc<Trait> (where Trait: DowncastSync) to Arc<Any>, which can then be downcast into Arc<ConcreteType> where ConcreteType implements Trait.
Source§

impl<T> DropBytes for T
where T: 'static,

Source§

unsafe fn drop_bytes(bytes: &mut [MaybeUninit<u8>])

Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<'a, S, I> Get<'a, I> for S
where I: GetIndex<'a, S>,

Source§

type Output = <I as GetIndex<'a, S>>::Output

Source§

fn get(&self, idx: I) -> Option<<I as GetIndex<'a, S>>::Output>

Source§

fn at(&self, idx: I) -> Self::Output

Return a value at the given index. This is provided as the checked version of get that will panic if the equivalent get call is None, which typically means that the given index is out of bounds. Read more
Source§

unsafe fn at_unchecked(&self, idx: I) -> Self::Output

Return a value at the given index. Read more
Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> IntoEither for T

Source§

fn into_either(self, into_left: bool) -> Either<Self, Self>

Converts self into a Left variant of Either<Self, Self> if into_left is true. Converts self into a Right variant of Either<Self, Self> otherwise. Read more
Source§

fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
where F: FnOnce(&Self) -> bool,

Converts self into a Left variant of Either<Self, Self> if into_left(&self) returns true. Converts self into a Right variant of Either<Self, Self> otherwise. Read more
Source§

impl<S, I> Isolate<I> for S
where I: IsolateIndex<S>,

Source§

type Output = <I as IsolateIndex<S>>::Output

Source§

unsafe fn isolate_unchecked(self, idx: I) -> <S as Isolate<I>>::Output

Unchecked version of isolate. Read more
Source§

fn try_isolate(self, idx: I) -> Option<<S as Isolate<I>>::Output>

Source§

fn isolate(self, idx: I) -> Self::Output
where Self: Sized,

Return a value at the given index. This is provided as the checked version of try_isolate that will panic if the equivalent try_isolate call is None, which typically means that the given index is out of bounds. Read more
Source§

impl<T> PartialEqBytes for T
where T: PartialEq + 'static,

Source§

unsafe fn eq_bytes(a: &[MaybeUninit<u8>], b: &[MaybeUninit<u8>]) -> bool

Source§

impl<T> Pointable for T

Source§

const ALIGN: usize

The alignment of pointer.
Source§

type Init = T

The type for initializers.
Source§

unsafe fn init(init: <T as Pointable>::Init) -> usize

Initializes a with the given initializer. Read more
Source§

unsafe fn deref<'a>(ptr: usize) -> &'a T

Dereferences the given pointer. Read more
Source§

unsafe fn deref_mut<'a>(ptr: usize) -> &'a mut T

Mutably dereferences the given pointer. Read more
Source§

unsafe fn drop(ptr: usize)

Drops the object pointed to by the given pointer. Read more
Source§

impl<T, N> PushArrayToVec<N> for T
where T: Clone, N: Array<T>,

Source§

fn push_to_vec(element: <N as Array<T>>::Array, set: &mut Vec<T>)

This method tells this type how it can be pushed to a Vec as an array.
Source§

impl<T> Same for T

Source§

type Output = T

Should always be Self
Source§

impl<SS, SP> SupersetOf<SS> for SP
where SS: SubsetOf<SP>,

Source§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
Source§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
Source§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
Source§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
Source§

impl<T> ToOwned for T
where T: Clone,

Source§

type Owned = T

The resulting type after obtaining ownership.
Source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
Source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

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

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.
Source§

impl<T> AttributeValue for T
where T: Clone + PartialEq + Debug + Send + Sync + 'static,

Source§

impl<T> Elem for T
where T: Any + DropBytes,

Source§

impl<T> Scalar for T
where T: 'static + Clone + PartialEq + Debug,