API
Roentgen.AbstractBeamLimitingDevice
— TypeAbstract type for Beam Limiting Devices are based on
AbstractBeamLimitingDevice
contains the types for MultiLeafCollimators, Jaws, etc.
Roentgen.AbstractBixel
— TypeAbstractBixel
Roentgen.AbstractDoseAlgorithm
— TypeAbstractDoseAlgorithm
Supertype for Dose Algorithms
Roentgen.AbstractDoseGrid
— TypeAbstractDoseGrid
A type of Dose Positions on a regular grid.
Must have axes
field defined.
Roentgen.AbstractExternalSurface
— TypeAbstractExternalSurface
Supertype for external surfaces.
Roentgen.AbstractFluenceElement
— TypeAbstractFluenceElement
Roentgen.AbstractTreatmentField
— TypeTreatment Field
Abstract treatment field, basis for containing multiple treatment types (e.g. VMAT, IMRT, etc.)
Roentgen.Beamlet
— MethodBeamlet(bixel::Bixel, gantry)
Construct a beamlet from a bixel
Roentgen.Bixel
— TypeBixel{T}
Roentgen.BixelGrid
— MethodBixelGrid(x, y, Δx[, Δy])
Construct a grid of bixels.
Each axis starts at the first element (e.g. x[1]), runs to the last element (x[end]), with uniform spacing (Δx). Same for y. Positions are "snapped" to the spacing value (see snapped_range
for details).
If Δy not specified, assumes Δy = Δx.
Roentgen.BixelGrid
— MethodBixelGrid(jaws::Jaws, Δx[, Δy])
Uses the jaw positions to construct a bixel grid.
If Δy not specified, assumes Δy = Δx
Roentgen.ConstantSurface
— TypeConstantSurface
An external surface with a constant source surface distance, regardless of the position. This surface is largely used for testing purposes as its inherently unphysical.
Roentgen.ControlPoint
— TypeControlPoint
Elements of a TreatmentField
Roentgen.CylinderBounds
— TypeCylinderBounds{T}
Represents a cylinder about the z
axis.
Roentgen.CylindricalSurface
— TypeCylindricalSurface(mesh::SimpleMesh, y::AbstractRange[, nϕ=181])
Construct from mesh
over axial axis y
.
Defaults to a 2° azimuthal spacing (nϕ=181).
Roentgen.CylindricalSurface
— TypeCylindricalSurface(y::TRange, ϕ::TRange, rho::TInterp)
Surface stored on a cylindrical-polar grid.
Roentgen.CylindricalSurface
— MethodCylindricalSurface(ϕ::AbstractVector, y::AbstractVector, rho::AbstractMatrix)
Constructed using vectors for ϕ, y, and rho.
Roentgen.CylindricalSurface
— MethodCylindricalSurface(mesh::SimpleMesh, Δy::Real[, nϕ=181])
Construct from mesh
over with axial spacing y
.
Uses the mesh bounds to compute the axial range. Defaults to a 2° azimuthal spacing (nϕ=181).
Roentgen.DoseGrid
— TypeDoseGrid(Δ, bounds::AbstractBounds)
Construct a DoseGrid
with spacing Δ
within bounds
.
Can specify a coordinate system transformation in transform
to generate dose points in a different coordinate system to the bounds.
Roentgen.DoseGrid
— TypeDoseGrid
Cartesian Dose Grid
Roentgen.DoseGridMasked
— TypeDoseGridMasked(Δ, bounds::AbstractBounds, transform=I)
Construct a DoseGridMasked
with spacing Δ
within bounds
.
The mask applies to all points within bounds
. Can specify a coordinate system transformation in transform
to generate dose points in a different coordinate system to the bounds.
Roentgen.DoseGridMasked
— TypeDoseGridMasked
Cartesian Dose Grid with a mask to reduce the number of dose points used.
Roentgen.DoseVolume
— TypeDoseVolume
Stores dose positions and an external surface
Roentgen.FinitePencilBeamKernel
— TypeFinitePencilBeamKernel(parameters, scalingfactor, depth, tanθ)
Dose calculation for Finite Pencil Beam Kernel algorithm (Jelen 2005).
Takes the following commissioned parameters:
parameters
: Weights and steepness parameters by depthscalingfactor
: Scaling factor matrix by depth and tanθdepth
: Depths of parameters of scaling factor valuestanθ
: Angle from central beam axis of scaling factor value
Roentgen.FinitePencilBeamKernel
— MethodFinitePencilBeamKernel(filename::String)
Load commissioned parameters from a .jld
file.
Roentgen.GantryPosition
— TypeGantryPosition{T}
The position/rotation of the gantry and beam-limiting device.
Stores:
gantry_angle
: As defined by IECcollimator_angle
: As defined by IEC (beam limiting device angle)source_axis_distance
: Distance between source and isocentercentral_beam_axis
: Unit vector pointing from the isocenter to the source in IEC Fixed coordinates.
Roentgen.Jaws
— TypeJaws
Jaws
stores the x and y positions of the jaws.
The x/y positions of the jaws can be accessed through the getx
/gety
methods.
The usual constructor directly takes the jaw x and y position vectors, but a single fieldsize can also be specified.
Roentgen.Jaws
— MethodJaws(fieldsize::T)
Create jaws with a square field of length fieldsize
.
Roentgen.LinearSurface
— TypeLinearSurface(params::AbstractVector)
Constructed with a vector of 6 element parameters corresponding to the plane position and normal at gantry angles.
Roentgen.LinearSurface
— MethodLinearSurface(mesh[, SAD=1000, ΔΘ=deg2rad(1)])
Construct a LinearSurface from a mesh.
Computes a set of planes parallel to the surface of the mesh.
Roentgen.MeshBounds
— TypeMeshBounds{T, TMesh}
Represents a surface from a mesh.
Roentgen.MeshSurface
— TypeMeshSurface
An external surface defined by a 3D mesh.
Roentgen.MockKernel
— TypeMockKernel
A mock dose algorithm used as an example and tests
Roentgen.MultiLeafCollimator
— TypeMultiLeafCollimator
MultiLeafCollimator
stores the leaf positions and edges of an MLC.
The 2 by n
leaf positions are stored in the positions
matrix. The leaf position boundaries are stored in the edges
vector, which is 1 element longer than the number of leaf tracks, n
.
Indexing a MultiLeafCollimator
will either return the x and y bounds of the leaf track, or another MultiLeafCollimator
. It also supports the @view
macro. Iteration goes through each leaf track, returning the lower and upper boundaries.
A MultiLeafCollimator
can also be shifted through addition or subtraction, and scaled with multplication and division.
The locate
method is implemented to return the leaf index containing edge position.
Roentgen.MultiLeafCollimator
— MethodMultiLeafCollimator(n::Int, Δy::Real)
Construct an MLC with leaf edges.
Roentgen.MultiLeafCollimator
— MethodMultiLeafCollimator(n::Int, Δy::Real)
Construct an MLC with n
number of leaves and of leaf width Δy
, centered about zero.
Roentgen.MultiLeafCollimatorSequence
— TypeMultiLeafCollimatorSequence
MultiLeafCollimatorSequence
stores a sequence of MLC apertures.
The 2 by n
leaf positions are stored in the positions
matrix. The leaf position boundaries are stored in the edges
vector, which is 1 element longer than the number of leaf tracks, n
.
Indexing a MultiLeafCollimator
will either return the x and y bounds of the leaf track, or another MultiLeafCollimator
. It also supports the @view
macro. Iteration goes through each leaf track, returning the lower and upper boundaries.
A MultiLeafCollimator
can also be shifted through addition or subtraction, and scaled with multplication and division.
The locate
method is implemented to return the leaf index containing edge position.
Roentgen.PlaneSurface
— TypePlaneSurface
A plane at a constant distance from and normal towards the source
The source-surface distance is stored in surf.source_surface_distance
Roentgen.SurfaceBounds
— TypeSurfaceBounds{TSurface}
Represents a surface from a mesh.
Base.getindex
— MethodBase.getindex(bixel::Bixel, i::Int)
Get the centre position of the bixel.
Base.iterate
— Methoditerate
Iteration of a treatment field, returning ControlPoint every time iterate is called.
Roentgen._assert_size
— Method_assert_size(D, pos, beamlets)
Ensures size of D is correct
Roentgen._get_cells
— Method_get_cells(pos)
Returns a list of cell connectivities
Roentgen._sequential_searchsortedlast
— Method_sequential_searchsortedlast(a, x, jstart)
A searchsortedlast
where the index is near and greater than previous index
This assumes that the next index is larger than the starting index jstart
. It first checks that jstart
already satisfies a[j]<=x<a[j+1]
. If not, it searches for j>jstart
, returning length(a)
if a[end]<x
.
Roentgen._vtk_beamlet_points
— Method_vtk_beamlet_points(beamlet, L)
Computes the vertices of the pyramid covered by a beamlet.
See write_vtk(filename, beamlet::Beamlet, L=2)
for further details
Roentgen._vtk_create_cell
— Method_vtk_create_cell(cell)
Return the VTK cell type for a list of cell indices.
Roentgen.bixels_from_bld
— Functionbixels_from_bld(args::AbstractBeamLimitingDevice...)
Create bixels corresponding to the provided beam limiting devices.
Roentgen.bixels_from_bld
— Methodbixels_from_bld(mlcx, mlc::MultiLeafCollimator, jaws::Jaws)
From MultiLeafCollimator and Jaws.
Roentgen.bld_to_fixed
— Methodbld_to_fixed(ϕg, θb, SAD)
Convert from IEC BLD to IEC fixed coordinates
Roentgen.bld_to_gantry
— Methodbld_to_gantry(θb, SAD)
Convert from IEC BLD to IEC gantry coordinates
Roentgen.calibrate!
— Functioncalibrate!(calc, MU, fieldsize, SAD[, SSD=SAD])
Calibrate a dose algorithm with given MU
, fielsize
and SAD
.
Scales the dose such that the maximum dose is 1 Gy for MU
monitor units, given fieldsize
and source-axis distance (SAD
). Can set source-surface distance SSD
if SSD!=SAD
.
Roentgen.closest_intersection
— Methodclosest_intersection(p1, p2, mesh::Domain)
Return the point on the line p1->p2
and mesh
closest to p1
Finds all points of intersection, then returns the point closest to p1
. If no intersections are found, it returns nothing
.
Roentgen.dose_fluence_matrix!
— Functiondose_fluence_matrix!(D<:AbstractMatrix, pos, beamlets, surf, calc)
Compute a dose-fluence matrix from dose positions, beamlets, external surface and dose calculation algorithm.
Requires the point_dose
method to be defined for the given dose calculation algorithm (calc
). point_dose
computes the dose at a given position from a given beamlet.
It stores result in D<:AbstractMatrix
. Currently, the following matrix types are supported:
D::Matrix
D::SparseMatrixCSC
D::CuArray
Roentgen.dose_fluence_matrix!
— Methoddose_fluence_matrix!(D, vol, beamlets, calc)
Gets positions and surface from vol
Roentgen.dose_fluence_matrix
— Methoddose_fluence_matrix(T, vol, beamlets, calc)
Gets positions and surface from vol
Roentgen.dose_fluence_matrix
— Methoddose_fluence_matrix(T, pos, beamlets, surf, calc)
Compute a dose-fluence matrix from dose positions, beamlets, external surface and dose calculation algorithm.
T
is the matrix type. It currently supports:
Matrix
: Dense CPUSparseMatrixCSC
: Sparse CPUCuArray
: Dense GPU
e.g. dose_fluence_matrix(SparseMatrixCSC, ...)
will create a sparse matrix computed using the CPU.
See dose_fluence_matrix!
for implementation.
Roentgen.dose_kernel!
— Methoddose_kernel!(rowval, nzval, pos::AbstractVector{T}, bixels, surf, calc)
Compute the fluence kernel for a given position.
Designed to be used with a dose-fluence matrix of type SparseMatrixCSC. Stores the row value in rowval
, and dose value in nzval
.
Roentgen.extent
— Methodextent(bounds::CylinderBounds)
For a cylinder
Roentgen.extent
— Methodextent(bounds::MeshBounds)
For a mesh
Roentgen.extent
— Methodextent(bounds::AbstractBounds)
Returns the bounding box of the bounds as two vectors: pmin
, pmax
Roentgen.extent
— Methodextent(bounds::SurfaceBounds)
For a mesh
Roentgen.fixed_to_bld
— Methodfixed_to_bld(ϕg, θb, SAD)
Convert from IEC fixed to IEC BLD coordinates
Roentgen.fixed_to_gantry
— Methodfixed_to_gantry(ϕg)
Convert from IEC fixed to IEC gantry coordinates
Roentgen.fluence!
— Methodfluence!(Ψ::AbstractArray{<:AbstractFloat}, bixels::AbstractArray{<:AbstractBixel}, index::AbstractArray{Int}, args...)
Stores the result in Ψ. See fluence(bixels::AbstractArray{<:AbstractBixel}, index::AbstractArray{Int}, args...)
for details
Roentgen.fluence!
— Methodfluence!(Ψ::AbstractArray{<:AbstractFloat}, bixels::AbstractArray{<:AbstractBixel}, args...)
Compute the fluence on a collection of bixels, storing the result in Ψ.
See fluence(bixels::AbstractArray{<:AbstractBixel}, args...)
for details.
Roentgen.fluence
— Methodfluence(bixels::AbstractArray{<:AbstractBixel}, index::AbstractArray{Int}, args...)
Allows precomputation of the location of the bixel in relation to the beam limiting device.
In the case of an MLC, the index is the leaf index that contains that bixel.
Requires fluence(bixel::AbstractBixel, index::Int, args...)
to be defined for the particular beam limiting device
Roentgen.fluence
— Methodfluence(bixels::AbstractArray{<:AbstractBixel}, bld::AbstractBeamLimitingDevice, args...)
Compute the fluence on a collection of bixels.
Broadcasts over the specific fluence(bixel, ...)
method for the provided beam limiting device. e.g.: fluence(bixels, jaws)
, fluence(bixels, mlcx, mlc)
Roentgen.fluence
— Methodfluence(bixel::Bixel, jaws::Jaws)
From the Jaws.
Roentgen.fluence
— Methodfluence(bixel::Bixel{T}, index::Int, mlcx1, mlcx2)
From an MLC aperture sequence using a given leaf index.
Roentgen.fluence
— Methodfluence(bixel::Bixel, index::Int, mlcx)
From an MLC aperture using a given leaf index.
This method assumes the bixel is entirely within the i
th leaf track, and does not overlap with other leaves. Does not check whether these assumptions are true.
Roentgen.fluence
— Methodfluence(bixel::AbstractBixel, bld::AbstractBeamLimitingDevice, args...)
Compute the fluence of bixel
from beam limiting device (e.g. an MLC or jaws).
Roentgen.fluence
— Methodfluence(bixel::Bixel{T}, mlc1::MultiLeafCollimator, mlc2::MultiLeafCollimator)
From an MLC aperture sequence.
Roentgen.fluence
— Methodfluence(bixel::Bixel, mlc::MultiLeafCollimator)
From an MLC aperture.
Roentgen.fluence_from_moving_aperture
— Methodfluence_from_moving_aperture(bixel::Bixel{T}, mlcx1, mlcx2)
From MLC leaf positions which move from mlcx1
to mlcx2
.
Computes the time-weighted fluence as the MLC moves from position mlcx1
to mlcx2
. Assumes the MLC leaves move in a straight line.
Roentgen.fluence_from_rectangle
— Methodfluence_from_rectangle(bixel::Bixel, xlim, ylim)
Compute the fluence of a rectangle with edges at xlim
and ylim
on a bixel.
Roentgen.fluence_onesided
— Methodfluence_onesided(xs, xf, xL, xU)
Compute fluence for 1 leaf position, assuming the other is infinitely far away
Computes the fluence for a leaf trajectory from xs
to xf
, over a bixel from xL
to xU
. The aperture is considered open to the right of the leaf position (i.e. leaf on the B bank). Assumes the bixel is fully in the leaf track.
Roentgen.gantry_to_bld
— Methodgantry_to_bld(θb, SAD)
Convert from IEC gantry to IEC BLD coordinates
Roentgen.gantry_to_fixed
— Methodgantry_to_fixed(ϕg)
Convert from IEC gantry to IEC fixed coordinates
Roentgen.getSSD
— MethodgetSSD(calc::ConstantSurface, pos, src)
When applied to a ConstantSurface
, it returns a constant source surface distance, regardless of pos
.
Roentgen.getSSD
— MethodgetSSD(surf::PlaneSurface, pos, src)
When applied to a PlaneSurface
, it returns the distance to the plane.
Roentgen.getSSD
— MethodgetSSD(surf::AbstractExternalSurface, pos, src)
Get the Source-Surface Distance (SSD) for position pos
to the radiation source src
.
Roentgen.getarea
— Methodgetarea(bixel::Bixel)
Return the area of the bixel.
Roentgen.getaxes
— Methodgetaxes(pos::DoseGrid[, dim])
Return the axes of the grid.
Optionally, can specify which dimension.
Roentgen.getaxes
— Methodgetaxes(pos::AbstractDoseGrid[, dim])
Return the axes of the grid.
Optionally, can specify which dimension.
Roentgen.getcenter
— Methodgetcenter(bixel::Bixel, i::Int)
Get the ith coordinate of the position.
Roentgen.getcenter
— Methodgetcenter(bixel::Bixel)
Get the position of the bixel.
Roentgen.getdepth
— Methodgetdepth(surf::AbstractExternalSurface, pos, src)
Get the depth of the position pos
below the surface from the radiation source src
.
Computes the depth by subtracting
Roentgen.getedge
— Methodgetedge(bixel::Bixel[, dim::Int])
Return the lower edge of the bixel.
Can specify a dim
for which dimension (x or y).
Roentgen.getpositions
— Methodgetpositions(vol::DoseVolume)
Get dose positions/grid
Roentgen.getsurface
— Methodgetsurface(vol::DoseVolume)
Get surface
Roentgen.getwidth
— Methodgetwidth(bixel::Bixel, i::Int)
Return the width of the bixel along axis i
.
Roentgen.getwidth
— Methodgetwidth(bixel::Bixel)
Return the width of the bixel.
Roentgen.hypotenuse
— Methodhypotenuse(a, b)
Compute the hypotenuse of triangle
Roentgen.interp
— Methodinterp(xg::AbstractVector, yg::AbstractVector, fg::AbstractMatrix, xi, yi)
Bilinear interpolation at position xi,yi
, on grid xg-yg
with values fg
Roentgen.interp
— Methodinterp(f1, f2, α)
Interpolate using given values f1
and f2
and the scaled distance from them α
Roentgen.interp
— Methodinterp(xg::AbstractVector{T}, fg::AbstractVector{T}, xi) where T<:Real
Interpolate xi
within uniform grid positions xg
and grid values fg
Roentgen.intersect_mesh
— Methodintersect_mesh(s, mesh)
Find the intersection points of the line
and the mesh
.
Returns a list of intersection points, and an empty list if none present.
Roentgen.intersect_mesh
— Methodintersect_mesh(line::Segment, mesh::Partition[, boxes])
Intersect a partitioned mesh.
Will check intersection with the bounding box of each partition before checking for intersection between each cell. Can pre-compute bounding boxes for extra performance. Single threaded.
Roentgen.leaf_trajectory
— Methodleaf_trajectory(x, x1, x2)
Compute the height of position x
between (x1, 0)
and (x2, 1)
Used in intersection_area
.
Roentgen.load_beam
— Methodload_beam(beam, total_meterset)
Load a beam from a control point sequence in a DICOM RP file.
Roentgen.load_dicom
— Methodload_dicom(filename)
Load a DICOM RP file into a Vector{TreatmentField}.
Roentgen.load_ref_dose
— Methodload_ref_dose(beam)
Load a reference dose, used for calculating the meterset in a control point sequence
Roentgen.load_structure_from_ply
— Methodload_structure_from_ply(filename)
Load a .ply file into a mesh
Uses MeshIO to load the mesh from the file. Unfortunately, MeshIO loads into GeometryBasics.Mesh, not Meshes.SimpleMesh. The rest of the code converts the mesh into a Meshes.SimpleMesh. From https://github.com/JuliaIO/MeshIO.jl/issues/67#issuecomment-913353708
Roentgen.locate
— Methodlocate(x::T, start::T, step::T) where T<:AbstractFloat
Locate the grid index for position xi
Specify the grid starting position (start
) and grid spacing (step
).
Roentgen.locate
— Methodlocate(mlc::AbstractMultiLeafCollimator, x::Number)
Locate the track index i
which contains the position x
Roentgen.locate
— Methodlocate(xg::AbstractRange{T}, xi::T) where T<:Real
When a range is provided, uniform spacing can be assumed.
Roentgen.locate
— Methodlocate(xg::AbstractVector{T}, xi::T) where T<:Real
Locates the grid index for position x
within a grid xg
. Does not limit out of bounds indices.
Roentgen.overlap
— Methodoverlap(x, Δx, xB, xA)
Compute the overlapping length of a bixel spanning x
-w
/2->x
+w
/2 and the length between xL
and xU
, normalised to the length of the bixel. If the Bixel fully within range (xL<=x-w/2 && x+w/2<=xU), return 1. If the bixel is fully outside the range (xU<=x-w/2 || x+w/2<=xL), return 0.
Roentgen.patient_to_fixed
— Methodpatient_to_fixed(isocenter)
Convert from the patient-based coordinate system to IEC Fixed.
Isocenter is in the patient-based coordinate system, as specified in the DICOM RP file.
Roentgen.point_dose
— Functionpoint_dose(p::SVector{3}, beamlet::Beamlet, surf::AbstractExternalSurface, calc::AbstractDoseAlgorithm)
Compute the dose at position p
from beamlet
using a dose algorithm.
Roentgen.point_dose
— Methodpoint_dose(..., calc::FinitePencilBeamKernel)
Using the Finite Pencil Beam Kernel algorithm
Roentgen.reconstruct_dose
— Methodreconstruct_dose(vol::AbstractDoseVolume, plan::AbstractTreatmentPlan,
calc::AbstractDoseAlgorithm; Δx=5., show_progress=true, maxradius=100.)
Dose reconstruction from a treatment plan.
Requires a dose volume (vol
), a treatment plan (plan
), and a dose calculation algorithm (calc
).
Optional arguments include:
Δx
: Size of each bixel in the fluence grid (defaults to 5.)show_progress
: If true (default), displays the progressmaxradius
: Kernel truncation radius
Roentgen.resample
— Methodresample(field, Δη::T; by=:MU)
Resample at uniform steps Δη
, from start to finish.
See resample(field, ηₛ::AbstractVector{T}; by=:MU)
for details.
Roentgen.resample
— Methodresample(field, ηₛ::AbstractVector{T}; by=:time)
Resample a treatment field onto new times or meterset values.
Can resample either by time (by=:time
) or MU (by=:MU
, default).
Roentgen.scale_to_cell
— Methodscale_to_cell(x1, x2, xi)
Scale the position xi
within the positions x1
and x2
.
Roentgen.scale_to_isoplane
— Methodscale_to_isoplane(pᵢ, z_plane)
Scale the position pᵢ
to the plane at distance z_plane
.
Roentgen.snapped_range
— Methodsnapped_range(x1, x2, Δ)
Create a range from x1 to x2 which is "snapped" to the step Δ.
Positions are "snapped" to the step value (e.g. a starting position of x[1]-0.2Δx snaps to x[1]-Δx). The new range always includes the start and end points of the original range
Examples:
- 0.1:1.:9.4 -> 0:1.:10.
Roentgen.subdivide
— Methodsubdivide(bixel::Bixel, nx::Integer, ny::Integer)
Subdivide a bixel by specifing the number of partitions nx
and ny
.
Returns a grid of bixels.
Roentgen.subdivide
— Methodsubdivide(bixel::Bixel{T}, δx::T, δy::T)
Subdivide by specifing widths δx
and δy
Roentgen.transform!
— Methodtransform!(mesh::SimpleMesh, trans)
Apply general transformation trans
to mesh
, modifying the original mesh.
Roentgen.transform
— Methodtransform(mesh::SimpleMesh, trans)
Apply general transformation trans
to mesh
`.
Roentgen.within
— Methodwithin(bounds::CylinderBounds, p)
Whether p
is within the cylinder
Roentgen.within
— Methodwithin(bounds::MeshBounds, p)
Whether p
is within the mesh
Roentgen.within
— Methodwithin(bounds::AbstractBounds, p)
Returns true
if the point p
is within bounds
Roentgen.within
— Methodwithin(bounds::SurfaceBounds, p)
Whether p
is within the mesh
Roentgen.write_vtk
— Functionwrite_vtk(filename, mesh::Beamlet, L=2)
Save a Beamlet
to a VTK (.vtu) file.
The beamlet is drawn from the source position to a length determined by parameter L
. This is the length of the beamlet in units of source-axis distance.
Roentgen.write_vtk
— Functionwrite_vtk(filename, mesh::Vector{<:AbstractBeamlet}, L=2)
Save a vector of Beamlet
to a VTK (.vtu) file.
See write_vtk(filename, beamlet::Beamlet, L=2)
for further details
Roentgen.write_vtk
— Functionwrite_vtk(filename::String, data, args...)
Write a data structure to the VTK file format.
See the individual methods for information.
Roentgen.write_vtk
— Methodwrite_vtk(filename, mesh::SimpleMesh)
Save a SimpleMesh
to a VTK (.vtu) file.
Roentgen.write_vtk
— Methodwrite_vtk(filename, surf::CylindricalSurface)
Save a CylindricalSurface
to a VTK (.vtu) file.
Roentgen.write_vtk
— Methodwrite_vtk(filename::String, pos::DoseGrid, "data1"=>data1, ...)
Save DoseGrid
to the VTK Image data (vti) format.
Can add point data to visualise on the 3D grid by chaining "name"=>data
pairs: e.g. "dose"=>dose, "depth"=>depth
. Also supports Dict
.
Roentgen.write_vtk
— Methodwrite_vtk(filename::String, pos::DoseGridMasked, "data1"=>data1, ...)
Save DoseGridMasked
to the VTK Unstructured Grid (vtu) format.
Can add point data to visualise on the 3D grid by chaining "name"=>data
pairs: e.g. "dose"=>dose, "depth"=>depth
. Also supports Dict
.
Roentgen.write_vtk
— Methodwrite_vtk(filename, surf::MeshSurface)
Save a MeshSurface
to a VTK (.vtu) file.