ASE-VTK

For ASE, the vtk interface consists of Python modules for automatic visualization of positions, bonds, forces and volume data (e.g. wave functions) from an Atoms object, provided such data is made available by the calculator.

Note

The Python modules in ASE are intended to wrap lower-level functionality of the VTK object models in small and easy-to-comprehend classes. To be able to distinguish between build-in VTK objects and their wrappers, and because VTK uses the CamelCase naming convention whereas ASE uses lower-case cf. our coding conventions, all variables referring to VTK built-in types are prefixed by vtk_. However, both VTK and wrapper classes are named according to the standard vtkFooBar.

Representing atoms

class ase.visualize.vtk.atoms.vtkAtoms(atoms, scale=1)[source]

Bases: ase.visualize.vtk.module.vtkModuleAnchor, ase.visualize.vtk.grid.vtkAtomicPositions

Provides fundamental representation for Atoms-specific data in VTK.

The vtkAtoms class plots atoms during simulations, extracting the relevant information from the list of atoms. It is created using the list of atoms as an argument to the constructor. Then one or more visualization modules can be attached using add_module(name, module).

Example:

>>> va = vtkAtoms(atoms)
>>> va.add_forces()
>>> va.add_axes()
>>> XXX va.add_to_renderer(vtk_ren)

Construct a fundamental VTK-representation of atoms.

atoms: Atoms object or list of Atoms objects
The atoms to be plotted.
scale = 1: float or int
Relative scaling of all Atoms-specific visualization.
add_axes()[source]

Add an orientation indicator for the cartesian axes. An appropriate vtkAxesModule is added to the module anchor under axes.

add_cell()[source]

Add a box outline of the cell using atoms.get_cell(). The existing vtkUnitCellModule is added to the module anchor under cell.

add_forces()[source]

Add force vectors for the atoms using atoms.get_forces(). A vtkGlyphModule is added to the module anchor under force.

add_velocities()[source]

Add velocity vectors for the atoms using atoms.get_velocities(). A vtkGlyphModule is added to the module anchor under velocity.

Atom-centered data

The superclass vtkAtomicPositions implements the basic concepts for representing atomic-centered data in VTK.

class ase.visualize.vtk.grid.vtkAtomicPositions(pos, cell)[source]

Provides an interface for adding Atoms-centered data to VTK modules. Atomic positions, e.g. obtained using atoms.get_positions(), constitute an unstructured grid in VTK, to which scalar and vector can be added as point data sets.

Just like Atoms, instances of vtkAtomicPositions can be divided into subsets, which makes it easy to select atoms and add properties.

Example:

>>> cell = vtkUnitCellModule(atoms)
>>> apos = vtkAtomicPositions(atoms.get_positions(), cell)
>>> apos.add_scalar_property(atoms.get_charges(), 'charges')
>>> apos.add_vector_property(atoms.get_forces(), 'forces')

Construct basic VTK-representation of a set of atomic positions.

pos: NumPy array of dtype float and shape (n,3)
Cartesian positions of the atoms.
cell: Instance of vtkUnitCellModule of subclass thereof
Holds information equivalent to that of atoms.get_cell().
add_scalar_property(data, name=None, active=True)[source]

Add VTK-representation of scalar data at the atomic positions.

data: NumPy array of dtype float and shape (n,)
Scalar values corresponding to the atomic positions.
name=None: str
Unique identifier for the scalar data.
active=True: bool
Flag indicating whether to use as active scalar data.
add_vector_property(data, name=None, active=True)[source]

Add VTK-representation of vector data at the atomic positions.

data: NumPy array of dtype float and shape (n,3)
Vector components corresponding to the atomic positions.
name=None: str
Unique identifier for the vector data.
active=True: bool
Flag indicating whether to use as active vector data.
get_points(subset=None)[source]

Return (subset of) vtkPoints containing atomic positions.

subset=None: list of int
A list of indices into the atomic positions; ignored if None.
get_unstructured_grid(subset=None)[source]

Return (subset of) an unstructured grid of the atomic positions.

subset=None: list of int
A list of indices into the atomic positions; ignored if None.

Predefined shapes

The class vtkGlyphModule implements the lower-level objects for representing predefined shapes (glyphs) in VTK.

class ase.visualize.vtk.module.vtkGlyphModule(vtk_pointset, glyph_source, scalemode=None, colormode=None)[source]

Modules represent a unified collection of VTK objects needed for introducing basic visualization concepts such as surfaces or shapes.

A common trait of all modules is the need for an actor representation and corresponding generic properties such as lighting, color and transparency.

Poly data modules are based on polygonal data sets, which can be mapped into graphics primitives suitable for rendering within the VTK framework.

Glyph modules construct these polygonal data sets by replicating a glyph source across a specific set of points, using available scalar or vector point data to scale and orientate the glyph source if desired.

Example:

>>> atoms = molecule('C60')
>>> cell = vtkUnitCellModule(atoms)
>>> apos = vtkAtomicPositions(atoms.get_positions(), cell)
>>> vtk_ugd = apos.get_unstructured_grid()
>>> glyph_source = vtkAtomSource('C')
>>> glyph_module = vtkGlyphModule(vtk_ugd, glyph_source)

Construct VTK-representation of a module containing glyphs. These glyphs share a common source, defining their geometrical shape, which is cloned and oriented according to the input data.

vtk_pointset: Instance of vtkPointSet or subclass thereof
A vtkPointSet defines a set of positions, which may then be assigned specific scalar of vector data across the entire set.
glyph_source: Instance of ~vtk.vtkCustomGlyphSource or subclass thereof
Provides the basic shape to distribute over the point set.
get_actor()

Return the actor representing this module in a rendering scene.

set_actor(vtk_act)

Set the actor representing this module in a rendering scene.

set_property(vtk_property)

Set the property of the actor representing this module.