API reference
Index
WriteVTK.CollectionFile
WriteVTK.DatasetFile
WriteVTK.MultiblockFile
WriteVTK.VTKBlock
WriteVTK.VTKFile
Base.close
Base.isopen
Base.setindex!
Base.setindex!
WriteVTK.add_field_data
WriteVTK.add_loc_attributes
WriteVTK.data_to_xml
WriteVTK.data_to_xml_appended
WriteVTK.data_to_xml_ascii
WriteVTK.data_to_xml_inline
WriteVTK.multiblock_add_block
WriteVTK.paraview_collection
WriteVTK.paraview_collection_load
WriteVTK.pvtk_grid
WriteVTK.pvtk_grid
WriteVTK.save_with_appended_data
WriteVTK.vtk_grid
WriteVTK.vtk_grid
WriteVTK.vtk_grid
WriteVTK.vtk_grid
WriteVTK.vtk_grid
WriteVTK.vtk_grid
WriteVTK.vtk_grid
WriteVTK.vtk_multiblock
WriteVTK.vtk_multiblock
WriteVTK.vtk_surface
WriteVTK.vtk_write_array
WriteVTK
WriteVTK.CollectionFile
— TypeCollectionFile <: VTKFile
Handler for a ParaView collection file (.pvd
).
WriteVTK.DatasetFile
— TypeDatasetFile <: VTKFile
Handler for a data-containing VTK file.
WriteVTK.MultiblockFile
— TypeMultiblockFile <: VTKFile
Handler for a multiblock VTK file (.vtm
).
WriteVTK.VTKBlock
— TypeVTKBlock
Handler for a nested block in a multiblock file.
WriteVTK.VTKFile
— TypeVTKFile
Abstract type describing a VTK file that may be written using close
.
Base.close
— MethodBase.close(vtk::VTKFile) -> Vector{String}
Write and close VTK file.
Returns a list of paths pointing to the written VTK files (typically just one file, but can be more for e.g. MultiblockFile
).
Base.close(vtm::MultiblockFile) -> Vector{String}
Save and close multiblock file (.vtm
). The VTK files included in the multiblock file are also saved.
Base.isopen
— MethodBase.isopen(vtk::VTKFile)
Check if VTK file is still being written.
Base.setindex!
— Methodsetindex!(vtk::DatasetFile, data, name::AbstractString, [field_type]; [component_names])
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 field_type
argument 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.
The optional component_names
keyword argument allows to override the default component names when writing vector or tensor fields. It should be a tuple or a vector of strings (see below for an example).
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(), component_names = ("Ux", "Uy", "Uz")] = vel
vtk["time", VTKFieldData()] = time
# This also works, and will generally give the same result:
vtk["velocity", component_names = ("Ux", "Uy", "Uz")] = vel
vtk["time"] = time
Base.setindex!
— Methodsetindex!(vtk::DatasetFile, attributes, loc::AbstractFieldData)
Add attributes to point, cell, or field dataset type in a VTK file.
Example
Add "HigherOrderDegrees" dataset and cell dataset type attribute to VTK file.
vtk = vtk_grid(...)
vtk["HigherOrderDegrees", VTKCellData()] = [2; 3; 12]
vtk[VTKCellData()] = Dict("HigherOrderDegrees" => "HigherOrderDegrees")
Note that all three are possible and equivalent:
vtk[VTKCellData()] = Dict("HigherOrderDegrees" => "HigherOrderDegrees")
vtk[VTKCellData()] = "HigherOrderDegrees" => "HigherOrderDegrees"
vtk[VTKCellData()] = ("HigherOrderDegrees" => "HigherOrderDegrees",)
WriteVTK.add_field_data
— Methodadd_field_data(
vtk::DatasetFile, data, name::AbstractString, loc::AbstractFieldData;
component_names = nothing,
)
Add either point or cell data to VTK file.
WriteVTK.add_loc_attributes
— Methodadd_loc_attributes(vtk::DatasetFile, attributes, loc::AbstractFieldData)
Add attributes to point, cell, or field dataset type tags in a VTK file.
WriteVTK.data_to_xml
— Functiondata_to_xml(
vtk::DatasetFile, xParent::XMLElement, data,
name::AbstractString, Nc::Union{Int,AbstractFieldData} = 1;
component_names = nothing
)
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 sizes are in bytes. The header itself is not compressed, only the data is. For more details, see:
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.multiblock_add_block
— Functionmultiblock_add_block(
vtm::Union{MultiblockFile, VTKBlock},
vtk::VTKFile,
[name = ""],
) -> nothing
Add a block to a MultiblockFile
or a VTKBlock
.
multiblock_add_block(
vtm::Union{MultiblockFile, VTKBlock},
[name = ""],
) -> VTKBlock
Create a sub-block in a MultiblockFile
or a VTKBlock
.
Returns a new VTKBlock
.
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
end
WriteVTK.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
end
WriteVTK.pvtk_grid
— Methodpvtk_grid(
filename, args...;
part, nparts, extents, ismain = (part == 1), ghost_level = 0,
kwargs...,
)
Returns a handler representing a parallel VTK file, which can be eventually written to file with close
.
Positional and keyword arguments in args
and kwargs
are passed to vtk_grid
verbatim. Note that serial filenames are automatically generated from filename
and from the process id part
.
The following keyword arguments only apply to parallel VTK file formats.
Mandatory ones are:
part
: current (1-based) part id,nparts
: total number of parts (only required for unstructured grids),extents
: array specifying the partitioning of a structured grid across different processes (see below for details).
Optional ones are:
ismain
:true
if the current part idpart
is the main (the only one that will write the.pvtk
file),ghost_level
: ghost level.
Specifying extents for a structured grid
For structured grids, the partitioning of the dataset across different processes must be specified via the extents
argument. This is an array where each element represents the data extent associated to a given process.
For example, for a dataset of global dimensions $15×12×4$ distributed across 4 processes, this array may look like the following:
extents = [
( 1:10, 1:5, 1:4), # process 1
(10:15, 1:5, 1:4), # process 2
( 1:10, 5:12, 1:4), # process 3
(10:15, 5:12, 1:4), # process 4
]
Some important notes:
- the
extents
argument must be the same for all processes; - the extents must overlap, or VTK / ParaView will complain when trying to open the files;
- the length of the
extents
array gives the number of processes. Therefore, thenparts
argument is redundant and does not need to be passed.
WriteVTK.pvtk_grid
— Methodpvtk_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 = pvtk_grid(args...; kwargs...) do vtk
# do stuff with the `vtk` file handler
end
WriteVTK.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(filename,
xs::AbstractVector,
cells::AbstractVector{<:AbstractMeshCell};
kwargs...)
Create an unstructured mesh image data (.vtu
) file.
xs
is a vector of coordinates, such as a vector of SVector{3}
elements.
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:5, 2:1.:3, range(4, 5; length = 3)) # different kinds of ranges
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_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
end
WriteVTK.vtk_grid
— Methodvtk_grid(vtm::Union{MultiblockFile, VTKBlock}, [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 associated to vtm
and the number of existent blocks.
WriteVTK.vtk_grid
— Methodvtk_grid(filename,
X::AbstractMatrix,
cells::AbstractVector{<:AbstractMeshCell};
kwargs...)
Create an unstructured mesh image data (.vtu
) file.
X
is a matrix with each column containing the Cartesian coordinates of a point
WriteVTK.vtk_grid
— Methodvtk_grid(filename,
x::AbstractVector{T}, [y::AbstractVector{T}, [z::AbstractVector{T}]],
cells::AbstractVector{<:AbstractMeshCell};
kwargs...) where {T<:Number}
Create an unstructured mesh image data (.vtu
) file.
x
, y
and z
are vectors of containing the corresponding Cartesian coordinates of each point.
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_multiblock
— Methodvtk_multiblock([f::Function], filename) -> MultiblockFile
Initialise VTK multiblock file, linking multiple VTK dataset files.
Returns a handler for a multiblock file. To recursively save the multiblock file and linked dataset files, call close
on the returned handler.
Note that close
is implicitly called if the optional f
argument is passed. This is in particular what happens when using the do-block syntax.
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
end
WriteVTK.vtk_surface
— Functionvtk_surface([f::Function], filename, xs, ys, zs; kwargs...)
Create unstructured grid file (".vtu") representing a surface plot of a 2D function on coordinates (xs
, ys
) with values zs
.
The coordinates xs
and ys
can be given as:
- 1D arrays of respective dimensions
Nx
andNy
(for regular grids); - 2D arrays of dimensions
(Nx, Ny)
(for irregular grids).
The values zs
should be given in a 2D array of dimensions (Nx, Ny)
.
Including additional data
Optionally, one can write additional data to the generated file via the function f
, which works in the same way as for vtk_grid
.
As an example, it is common to colour surface plots by the height z
. To do this, one should write the zs
matrix as point data:
julia> xs = 0:0.5:10; ys = 0:1.0:20;
julia> zs = @. cos(xs) + sin(ys');
julia> vtk_surface("surf", xs, ys, zs) do vtk
vtk["z_values"] = zs
end
Note that the included data must have dimensions (Nx, Ny)
(for point data) or (Nx - 1, Ny - 1)
(for cell data).
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. For instance,
vtk_write_array(filename, (u, v), ("u", "v"))
In that case, the arrays must have the same dimensions.