WriteVTK.jl
Documentation for WriteVTK.jl
WriteVTK.VTKDataType — TypeVTKDataTypeUnion of data types allowed by VTK.
WriteVTK.AbstractFieldData — TypeAbstractFieldDataAbstract type representing any kind of dataset.
WriteVTK.AbstractVTKDataset — TypeAbstractVTKDatasetAbstract type representing any structured or unstructured VTK dataset.
The dataset classification is described in the VTK file format specification, page 12.
WriteVTK.MeshCell — TypeMeshCellSingle cell element in unstructured or polygonal grid.
It is characterised by a cell type (for instance, VTKCellType.TRIANGLE or PolyData.Strips) and by a connectivity vector determining the points on the grid defining this cell.
MeshCell(cell_type, connectivity::AbstractVector)Define a single cell element of an unstructured grid.
The cell_type argument characterises the type of cell (e.g. vertex, triangle, hexaedron, ...):
- cell types for unstructured datasets are defined in the
VTKCellTypes
module;
- cell types for polygonal datasets are defined in the
PolyDatamodule.
The connectivity argument is a vector containing the indices of the points passed to vtk_grid which define this cell.
Example
Define a triangular cell passing by points with indices [3, 5, 42].
julia> cell = MeshCell(VTKCellTypes.VTK_TRIANGLE, [3, 5, 42])
MeshCell{VTKCellType, Vector{Int64}}(VTKCellType("VTK_TRIANGLE", 0x05, 3), [3, 5, 42])WriteVTK.StructuredVTKDataset — TypeStructuredVTKDataset <: AbstractVTKDatasetAbstract type representing a structured VTK dataset.
Subtypes are VTKImageData, VTKRectilinearGrid and VTKStructuredGrid.
WriteVTK.UnstructuredVTKDataset — TypeUnstructuredVTKDataset <: AbstractVTKDatasetAbstract type representing an unstructured VTK dataset.
Subtypes are VTKPolyData and VTKUnstructuredGrid.
Base.close — Methodclose(vtk::VTKFile)Write and close VTK file.
Base.isopen — Methodisopen(vtk::VTKFile)Check if VTK file is still being written.
Base.setindex! — Methodsetindex!(vtk::DatasetFile, data, name::AbstractString, [field_type])Add a new dataset to VTK file.
The number of components of the dataset (e.g. for scalar or vector fields) is determined automatically from the input data dimensions.
The optional argument field_type should be an instance of VTKPointData, VTKCellData or VTKFieldData. It determines whether the data should be associated to grid points, cells or none. If not given, this is guessed from the input data size and the grid dimensions.
Example
Add "velocity" dataset and time scalar to VTK file.
vel = rand(3, 12, 14, 42) # vector field
time = 42.0
vtk = vtk_grid(...)
vtk["velocity", VTKPointData()] = vel
vtk["time", VTKFieldData()] = time
# This also works, and will generally give the same result:
vtk["velocity"] = vel
vtk["time"] = timeWriteVTK.add_field_data — Methodadd_field_data(vtk::DatasetFile, data,
name::AbstractString, loc::AbstractFieldData)Add either point or cell data to VTK file.
WriteVTK.data_to_xml — Functiondata_to_xml(
vtk::DatasetFile, xParent::XMLElement, data,
name::AbstractString, Nc::Union{Int,AbstractFieldData} = 1,
)Add numerical data to VTK XML file.
Data is written under the xParent XML node.
Nc may be either the number of components, or the type of field data. In the latter case, the number of components will be deduced from the data dimensions and the type of field data.
WriteVTK.data_to_xml_appended — Methoddata_to_xml_appended(vtk::DatasetFile, xDA::XMLElement, data)Add appended raw binary data to VTK XML file.
Data is written to the vtk.buf buffer.
When compression is enabled:
the data array is written in compressed form (obviously);
the header, written before the actual numerical data, is an array of
HeaderType (UInt32 / UInt64) values: [num_blocks, blocksize, last_blocksize, compressed_blocksizes] All the sizes are in bytes. The header itself is not compressed, only the data is. For more details, see: http://public.kitware.com/pipermail/paraview/2005-April/001391.html http://mathema.tician.de/what-they-dont-tell-you-about-vtk-xml-binary-formats (This is not really documented in the VTK specification...)
Otherwise, if compression is disabled, the header is just a single HeaderType value containing the size of the data array in bytes.
WriteVTK.data_to_xml_ascii — Methoddata_to_xml_ascii(vtk::DatasetFile, xDA::XMLElement, data)Add inline data to VTK XML file in ASCII format.
WriteVTK.data_to_xml_inline — Methoddata_to_xml_inline(vtk::DatasetFile, xDA::XMLElement, data)Add inline, base64-encoded data to VTK XML file.
WriteVTK.paraview_collection — Methodparaview_collection(f::Function, args...; kwargs...)Create VTK file and apply f to it. The file is automatically closed by the end of the call.
This allows to use the do-block syntax for creating VTK files:
saved_files = paraview_collection(args...; kwargs...) do vtk
# do stuff with the `vtk` file handler
endWriteVTK.paraview_collection_load — Methodparaview_collection_load(f::Function, args...; kwargs...)Create VTK file and apply f to it. The file is automatically closed by the end of the call.
This allows to use the do-block syntax for creating VTK files:
saved_files = paraview_collection_load(args...; kwargs...) do vtk
# do stuff with the `vtk` file handler
endWriteVTK.save_with_appended_data — MethodWrite VTK XML file containing appended binary data to disk.
In this case, the XML file is written manually instead of using the save_file function of LightXML, which doesn't allow to write raw binary data.
WriteVTK.vtk_grid — Methodvtk_grid(f::Function, args...; kwargs...)Create VTK file and apply f to it. The file is automatically closed by the end of the call.
This allows to use the do-block syntax for creating VTK files:
saved_files = vtk_grid(args...; kwargs...) do vtk
# do stuff with the `vtk` file handler
endWriteVTK.vtk_grid — Methodvtk_grid(vtm::MultiblockFile, [filename], griddata...; kwargs...)Create new dataset file that is added to an existent multiblock file. The VTK grid is specified by the elements of griddata.
If the filename is not given, it is determined automatically from the filename of the vtm file and the number of existent blocks.
WriteVTK.vtk_grid — Methodvtk_grid(vtb::vtkBlock, [filename], griddata...; kwargs...)Create new dataset file that is added to an existent VTKBlock. The VTK grid is specified by the elements of griddata.
If the filename is not given, it is determined automatically from the filename of the vtb file and the number of existent blocks.
WriteVTK.vtk_grid — Methodvtk_grid(filename::AbstractString,
x::AbstractVector{T}, y::AbstractVector{T}, [z::AbstractVector{T}];
kwargs...)Create 2D or 3D rectilinear grid (.vtr) file.
Coordinates are specified by separate vectors x, y, z.
Examples
julia> vtk = vtk_grid("abc", [0., 0.2, 0.5], collect(-2.:0.2:3), [1., 2.1, 2.3])
VTK file 'abc.vtr' (RectilinearGrid file, open)WriteVTK.vtk_grid — Methodvtk_grid(filename, x::AbstractRange{T}, y::AbstractRange{T}, [z::AbstractRange{T}];
kwargs...)Create image data (.vti) file.
Along each direction, the grid is specified in terms of an AbstractRange object.
Examples
julia> vtk = vtk_grid("abc", 1:0.2:5, 2:1.:3, 4:1.:5) # 3D dataset
VTK file 'abc.vti' (ImageData file, open)
julia> vtk = vtk_grid("abc", 1:0.2:5, 2:1.:3) # 2D dataset
VTK file 'abc.vti' (ImageData file, open)
julia> vtk = vtk_grid("def",
LinRange(0., 5., 10),
LinRange(0., 2π, 16),
LinRange(1., 10., 12))
VTK file 'def.vti' (ImageData file, open)
WriteVTK.vtk_multiblock — Methodvtk_multiblock(f::Function, args...; kwargs...)Create VTK file and apply f to it. The file is automatically closed by the end of the call.
This allows to use the do-block syntax for creating VTK files:
saved_files = vtk_multiblock(args...; kwargs...) do vtk
# do stuff with the `vtk` file handler
endWriteVTK.vtk_save — Methodvtk_save(vtm::MultiblockFile)Save and close multiblock file (.vtm). The VTK files included in the multiblock file are also saved.
WriteVTK.vtk_write_array — Methodvtk_write_array(filename, arrays, labels)
vtk_write_array(filename, array, label="array")Write Julia arrays to a VTK image data file (.vti).
Useful for general visualisation of arrays. The input can be a 2D or 3D array.
Multiple arrays can be given as a tuple, e.g.
vtk_write_array(filename, (x, y), ("x", "y"))In that case, the arrays must have the same dimensions.
WriteVTK.VTKCellTypes — ModuleVTKCellTypesModule defining cell types for unstructured datasets.
Definitions are adapted from the VTK source code.
WriteVTK.VTKCellTypes.nodes — Methodnodes(c::VTKCellTypes)Returns the number of nodes (or grid points) required by the cell type.
For instance, this returns 3 for VTK_TRIANGLE.
For cell types that can take any number of nodes, such as VTK_POLY_LINE, this returns -1.
WriteVTK.PolyData — ModulePolyDataDefines cell types for polygonal datasets.
The following singleton types are defined:
PolyData.Vertsfor vertices,PolyData.Linesfor lines,PolyData.Stripsfor triangular strips,PolyData.Polysfor polygons.