ovito.data
¶
This Python module defines various data object types, which are produced and processed within OVITO’s data pipeline system.
It also provides the DataCollection
class as a container for such data objects as well as several utility classes for
computing neighbor lists or iterating over the bonds of particles.
Data containers:
DataObject
(base class for all data object types)
DataCollection
(an entire dataset made of several data objects)
PropertyContainer
(base container class storing a set ofProperty
arrays)
Particles
(specializedPropertyContainer
for particles)
Bonds
(specializedPropertyContainer
for bonds)
VoxelGrid
(specializedPropertyContainer
for 2d and 3d data grids)
DataTable
(specializedPropertyContainer
for tabulated data)
Data objects:
Property
(array of per-data-element property values)
SimulationCell
(simulation box geometry and boundary conditions)
SurfaceMesh
(polyhedral mesh representation of the boundary surface of a spatial volume)
TrajectoryLines
(set of particle trajectory lines)
TriangleMesh
(general mesh data structure made of vertices and triangular faces)
DislocationNetwork
(set of discrete dislocation lines with Burgers vector information)
Auxiliary data objects:
ElementType
(base class for element type definitions)
ParticleType
(describes a single particle/atom type)
BondType
(describes a single bond type)
DislocationSegment
(a dislocation line in aDislocationNetwork
)
Utility classes:
CutoffNeighborFinder
(finds all neighboring particles within a cutoff distance)
NearestNeighborFinder
(finds N nearest neighbor particles)
BondsEnumerator
(used for efficiently iterating over the bonds of individual particles)
-
class
ovito.data.
BondType
¶ Base:
ovito.data.ElementType
Represents a bond type. This class inherits all its fields from the
ElementType
base class.You can enumerate the list of defined bond types by accessing the
bond_types
bond property object:bond_type_list = data.particles.bonds.bond_types.types for bond_type in bond_type_list: print(bond_type.id, bond_type.name, bond_type.color)
-
class
ovito.data.
Bonds
¶ Base:
ovito.data.PropertyContainer
Stores the list of bonds and their properties. A
Bonds
object is always part of a parentParticles
object. You can access it as follows:data = pipeline.compute() print("Number of bonds:", data.particles.bonds.count)
The
Bonds
class inherits thecount
attribute from itsPropertyContainer
base class. This attribute returns the number of bonds.Bond properties
Bonds can be associated with arbitrary bond properties, which are managed in the
Bonds
container as a set ofProperty
data arrays. Each bond property has a unique name by which it can be looked up:print("Bond property names:") print(data.particles.bonds.keys()) if 'Length' in data.particles.bonds: length_prop = data.particles.bonds['Length'] assert(len(length_prop) == data.particles.bonds.count)
New bond properties can be added using the
PropertyContainer.create_property()
method.Bond topology
The
Topology
bond property, which is always present, defines the connectivity between particles in the form of a N x 2 array of indices into theParticles
array. In other words, each bond is defined by a pair of particle indices.for a,b in data.particles.bonds.topology: print("Bond from particle %i to particle %i" % (a,b))
Note that the bonds of a system are not stored in any particular order. If you need to enumerate all bonds connected to a certain particle, you can use the
BondsEnumerator
utility class for that.Bonds visualization
The
Bonds
data object has aBondsVis
element attached to it, which controls the visual appearance of the bonds in rendered images. It can be accessed through thevis
attribute:data.particles.bonds.vis.enabled = True data.particles.bonds.vis.shading = BondsVis.Shading.Flat data.particles.bonds.vis.width = 0.3
Computing bond vectors
Since each bond is defined by two indices into the particles array, we can use these indices to determine the corresponding spatial bond vectors connecting the particles. They can be computed from the positions of the particles:
topology = data.particles.bonds.topology positions = data.particles.positions bond_vectors = positions[topology[:,1]] - positions[topology[:,0]]
Here, the first and the second column of the bonds topology array are used to index into the particle positions array. The subtraction of the two indexed arrays yields the list of bond vectors. Each vector in this list points from the first particle to the second particle of the corresponding bond.
Finally, we may have to correct for the effect of periodic boundary conditions when a bond connects two particles on opposite sides of the box. OVITO keeps track of such cases by means of the the special
Periodic Image
bond property. It stores a shift vector for each bond, specifying the directions in which the bond crosses periodic boundaries. We make use of this information to correct the bond vectors computed above. This is done by adding the product of the cell matrix and the shift vectors from thePeriodic Image
bond property:bond_vectors += numpy.dot(data.cell[:3,:3], data.particles.bonds.pbc_vectors.T).T
The shift vectors array is transposed here to facilitate the transformation of the entire array of vectors with a single 3x3 cell matrix. To summarize: In the two code snippets above we have performed the following calculation for every bond (a, b) in parallel:
v = x(b) - x(a) + dot(H, pbc)
where H is the cell matrix and pbc is the bond’s PBC shift vector of the form (n\ x, n\ y, n\ z).
Standard bond properties
The following standard properties are defined for bonds:
Property name
Python attribute
Data type
Component names
Bond Type
int
Color
float
R, G, B
Length
float
Periodic Image
int
X, Y, Z
Selection
int
Topology
int64
A, B
Transparency
float
-
property
bond_types
¶ The
Property
data array for theBond Type
standard bond property; orNone
if that property is undefined.
-
property
colors
¶ The
Property
data array for theColor
standard bond property; orNone
if that property is undefined.
-
create_bond
(a, b, type=None, pbcvec=None)¶ Creates a new bond between two particles a and b, both parameters being indices into the
Particles
list to which thisBonds
container belongs.- Parameters
a (int) – Index of first particle connected by the new bond. Particle indices start at 0.
b (int) – Index of second particle connected by the new bond.
type (int) – Optional type ID to be assigned to the new bond. This value will be stored to the
bond_types
array.pbcvec (tuple) – Three integers specifying the bond’s crossings of periodic cell boundaries. The information will be stored in the
pbc_vectors
array.
- Returns
The index of the newly created bond, i.e.
(Bonds.count-1)
.
The method does not check if there already is an existing bond connecting the same pair of particles.
The method does not check if the particle indices a and b do exist. Thus, it is your responsibility to ensure that both indices are in the range 0 to
(Particles.count-1)
.In case the
SimulationCell
has periodic boundary conditions enabled, and the two particles connected by the bond are located in different periodic images, make sure you provide the pbcvec argument. It is required so that OVITO does not draw the bond as a direct line from particle a to particle b but as a line passing through the periodic cell faces. You can use theParticles.delta_vector()
function to compute pbcvec or use thepbc_shift
vector returned by theCutoffNeighborFinder
utility.
-
property
pbc_vectors
¶ The
Property
data array for thePeriodic Image
standard bond property; orNone
if that property is undefined.
-
property
-
class
ovito.data.
BondsEnumerator
¶ Utility class that permits efficient iteration over the bonds connected to specific particles.
The constructor takes a
Bonds
object as input. From the generally unordered list of bonds, theBondsEnumerator
will build a lookup table for quick enumeration of bonds of particular particles.All bonds connected to a specific particle can be subsequently visited using the
bonds_of_particle()
method.Warning: Do not modify the underlying
Bonds
object while theBondsEnumerator
is in use. Adding or deleting bonds would render the internal lookup table of theBondsEnumerator
invalid.Usage example
from ovito.io import import_file from ovito.data import BondsEnumerator from ovito.modifiers import ComputePropertyModifier # Load a dataset containing atoms and bonds. pipeline = import_file('input/bonds.data.gz', atom_style='bond') # For demonstration purposes, lets here define a compute modifier that calculates the length # of each bond, storing the results in a new bond property named 'Length'. pipeline.modifiers.append(ComputePropertyModifier(operate_on='bonds', output_property='Length', expressions=['BondLength'])) # Obtain pipeline results. data = pipeline.compute() positions = data.particles.positions # array with atomic positions bond_topology = data.particles.bonds.topology # array with bond topology bond_lengths = data.particles.bonds['Length'] # array with bond lengths # Create bonds enumerator object. bonds_enum = BondsEnumerator(data.particles.bonds) # Loop over atoms. for particle_index in range(data.particles.count): # Loop over bonds of current atom. for bond_index in bonds_enum.bonds_of_particle(particle_index): # Obtain the indices of the two particles connected by the bond: a = bond_topology[bond_index, 0] b = bond_topology[bond_index, 1] # Bond directions can be arbitrary (a->b or b->a): assert(a == particle_index or b == particle_index) # Obtain the length of the bond from the 'Length' bond property: length = bond_lengths[bond_index] print("Bond from atom %i to atom %i has length %f" % (a, b, length))
-
class
ovito.data.
CutoffNeighborFinder
(cutoff, data_collection)¶ A utility class that computes particle neighbor lists.
This class lets you iterate over all neighbors of a particle that are located within a specified spherical cutoff. You can use it to build neighbors lists or perform computations that require neighbor vector information.
The constructor takes a positive cutoff radius and a
DataCollection
providing the input particles and theSimulationCell
(needed for periodic systems).Once the
CutoffNeighborFinder
has been constructed, you can call itsfind()
method to iterate over the neighbors of a particle, for example:from ovito.io import import_file from ovito.data import CutoffNeighborFinder # Load input simulation file. pipeline = import_file("input/simulation.dump") data = pipeline.compute() # Initialize neighbor finder object: cutoff = 3.5 finder = CutoffNeighborFinder(cutoff, data) # Prefetch the property array containing the particle type information: ptypes = data.particles.particle_types # Loop over all particles: for index in range(data.particles.count): print("Neighbors of particle %i:" % index) # Iterate over the neighbors of the current particle: for neigh in finder.find(index): print(neigh.index, neigh.distance, neigh.delta, neigh.pbc_shift) # The index can be used to access properties of the current neighbor, e.g. type_of_neighbor = ptypes[neigh.index]
Note: In case you rather want to determine the N nearest neighbors of a particle, use the
NearestNeighborFinder
class instead.-
find
(index)¶ Returns an iterator over all neighbors of the given particle.
- Parameters
index (int) – The zero-based index of the central particle whose neighbors should be enumerated.
- Returns
A Python iterator that visits all neighbors of the central particle within the cutoff distance. For each neighbor the iterator returns an object with the following property fields:
index: The zero-based global index of the current neighbor particle.
distance: The distance of the current neighbor from the central particle.
distance_squared: The squared neighbor distance.
delta: The three-dimensional vector connecting the central particle with the current neighbor (taking into account periodicity).
pbc_shift: The periodic shift vector, which specifies how often each periodic boundary of the simulation cell is crossed when going from the central particle to the current neighbor.
The index value returned by the iterator can be used to look up properties of the neighbor particle as demonstrated in the example above.
Note that all periodic images of particles within the cutoff radius are visited. Thus, the same particle index may appear multiple times in the neighbor list of the central particle. In fact, the central particle may be among its own neighbors in a small periodic simulation cell. However, the computed vector (
delta
) and PBC shift (pbc_shift
) will be unique for each visited image of the neighbor particle.
-
find_at
(coords)¶ Returns an iterator over all particles located within the spherical range of the given center position. In contrast to
find()
this method can search for neighbors around arbitrary spatial locations, which don’t have to coincide with any physical particle position.- Parameters
coords – A (x,y,z) coordinate triplet specifying the center location around which to search for particles.
- Returns
A Python iterator enumerating all particles within the cutoff distance. For each neighbor the iterator returns an object with the following properties:
index: The zero-based global index of the current neighbor particle.
distance: The distance of the current particle from the center position.
distance_squared: The squared distance.
delta: The three-dimensional vector from the center to the current neighbor (taking into account periodicity).
pbc_shift: The periodic shift vector, which specifies how often each periodic boundary of the simulation cell is crossed when going from the center point to the current neighbor.
The index value returned by the iterator can be used to look up properties of the neighbor particle as demonstrated in the example above.
Note that all periodic images of particles within the cutoff radius are visited. Thus, the same particle index may appear multiple times in the neighbor list. However, the computed vector (
delta
) and image offset (pbc_shift
) will be unique for each visited image of a neighbor particle.
-
neighbor_distances
(index)¶ Returns the list of distances between some central particle and all its neighbors within the cutoff range.
- Parameters
index (int) – The 0-based index of the central particle whose neighbors should be enumerated.
- Returns
NumPy array containing the radial distances to all neighbor particles within the cutoff range (in arbitrary order).
This method is equivalent to the following code, but performance is typically better:
def neighbor_distances(index): distances = [] for neigh in finder.find(index): distances.append(neigh.distance) return numpy.asarray(distances)
-
neighbor_vectors
(index)¶ Returns the list of vectors from some central particle to all its neighbors within the cutoff range.
- Parameters
index (int) – The 0-based index of the central particle whose neighbors should be enumerated.
- Returns
Two-dimensional NumPy array containing the vectors to all neighbor particles within the cutoff range (in arbitrary order).
The method is equivalent to the following code, but performance is typically better:
def neighbor_vectors(index): vecs = [] for neigh in finder.find(index): vecs.append(neigh.delta) return numpy.asarray(vecs)
-
-
class
ovito.data.
DataCollection
¶ Base:
ovito.data.DataObject
A
DataCollection
is a container class holding together individual data objects, each representing different fragments of a dataset. For example, a dataset loaded from a simulation data file may consist of particles, the simulation cell information and additional auxiliary data such as the current timestep number of the snapshots, etc. All this information is contained in oneDataCollection
, which exposes the individual pieces of information as sub-objects, for example via theDataCollection.particles
,DataCollection.cell
andDataCollection.attributes
fields.Data collections are the elementary entities that get processed within a data
Pipeline
. Each modifier receives a data collection from the preceding modifier, alters it in some way, and passes it on to the next modifier. The output data collection of the last modifier in the pipeline is returned by thePipeline.compute()
method.A data collection essentially consists of a bunch of
DataObjects
, which are all stored in theDataCollection.objects
list. Typically, you don’t access the data objects through this list directly but rather use on of the special accessor fields provided by theDataCollection
class, which give more convenient access to data objects of a particular kind. For example, thesurfaces
dictionary provides key-based access to all theSurfaceMesh
instances currently in the data collection.You can programmatically add or remove data objects from a data collection by manipulating its
objects
list. For instance, to populate a new data collection instance that is initially empty with a newSimulationCell
object:data = DataCollection() cell = SimulationCell() data.objects.append(cell) assert(data.cell is cell)
-
apply
(modifier, frame=None)¶ This method applies a
Modifier
function to the data in this collection in place.It allows modifying a data collection with one of Ovito’s built-in modifiers directly without the need to build up a complete
Pipeline
first. In contrast to a data pipeline, theapply()
method executes the modifier function immediately and changes the input data in place. In other words, the original data in thisDataCollection
will be replaced by the output produced by the invoked modifier function. Note that it is possible to create a copy of the original data collection using theclone()
method if needed. The following code example demonstrates how to use theapply()
method to successively modify the state of a dataset:from ovito.io import import_file from ovito.modifiers import * data = import_file("input/simulation.dump").compute() data.apply(CoordinationAnalysisModifier(cutoff=2.9)) data.apply(ExpressionSelectionModifier(expression="Coordination<9")) data.apply(DeleteSelectedModifier())
Note that it is typically possible to achieve the same result by first populating a
Pipeline
with the modifiers and then calling itscompute()
method at the very end:pipeline = import_file("input/simulation.dump") pipeline.modifiers.append(CoordinationAnalysisModifier(cutoff=2.9)) pipeline.modifiers.append(ExpressionSelectionModifier(expression="Coordination<9")) pipeline.modifiers.append(DeleteSelectedModifier()) data = pipeline.compute()
Also note that
apply()
may be called from a user-defined modifier function (seePythonScriptModifier
) in order to invoke a built-in modifier as a sub-routine:# A user-defined modifier function that calls the built-in ColorCodingModifier # as a sub-routine to assign a color to each atom based on some property # created within the function itself: def modify(frame, data): data.particles_.create_property('idx', data=numpy.arange(data.particles.count)) data.apply(ColorCodingModifier(property='idx'), frame) # Set up a data pipeline that uses the user-defined modifier function: pipeline = import_file("input/simulation.dump") pipeline.modifiers.append(modify) data = pipeline.compute()
- Parameters
modifier – The
Modifier
that will be called by the method to alter the data in place.frame (int) – An optional animation frame number that will be passed to the modifier function, which may use it to implement time-dependent modifications.
-
property
attributes
¶ This field contains a dictionary view with all the global attributes currently associated with this data collection. Global attributes are key-value pairs that represent small tokens of information, typically simple value types such as
int
,float
orstr
. Every attribute has a unique identifier such as'Timestep'
or'ConstructSurfaceMesh.surface_area'
. This identifier serves as look-up key in theattributes
dictionary. Attributes are dynamically generated by modifiers in a data pipeline or come from the data source. For example, if the input simulation file contains timestep information, the timestep number is made available by theFileSource
as the'Timestep'
attribute. It can be retrieved from pipeline’s output data collection:>>> pipeline = import_file('snapshot_140000.dump') >>> pipeline.compute().attributes['Timestep'] 140000
Some modifiers report their calculation results by adding new attributes to the data collection. See each modifier’s reference documentation for the list of attributes it generates. For example, the number of clusters identified by the
ClusterAnalysisModifier
is available in the pipeline output as an attribute namedClusterAnalysis.cluster_count
:pipeline.modifiers.append(ClusterAnalysisModifier(cutoff = 3.1)) data = pipeline.compute() nclusters = data.attributes["ClusterAnalysis.cluster_count"]
The
ovito.io.export_file()
function can be used to output dynamically computed attributes to a text file, possibly as functions of time:export_file(pipeline, "data.txt", "txt/attr", columns = ["Timestep", "ClusterAnalysis.cluster_count"], multiple_frames = True)
If you are writing your own modifier function, you let it add new attributes to a data collection. In the following example, the
CommonNeighborAnalysisModifier
first inserted into the pipeline generates the'CommonNeighborAnalysis.counts.FCC'
attribute to report the number of atoms that have an FCC-like coordination. To compute an atomic fraction from that, we need to divide the count by the total number of atoms in the system. To this end, we append a user-defined modifier function to the pipeline, which computes the fraction and outputs the value as a new attribute named'fcc_fraction'
.pipeline.modifiers.append(CommonNeighborAnalysisModifier()) def compute_fcc_fraction(frame, data): n_fcc = data.attributes['CommonNeighborAnalysis.counts.FCC'] data.attributes['fcc_fraction'] = n_fcc / data.particles.count pipeline.modifiers.append(compute_fcc_fraction) print(pipeline.compute().attributes['fcc_fraction'])
-
property
cell
¶ Returns the
SimulationCell
object, which stores the cell vectors and periodic boundary condition flags, orNone
if there is no cell information associated with this data collection.Note: The
SimulationCell
data object may be read-only if it is currently shared by several data collections. Use thecell_
field instead to request a mutable cell object if you intend to modify it.
-
clone
()¶ Returns a copy of this
DataCollection
containing the same data objects as the original.The method may be used to retain a copy of the original data before modifying a data collection in place, for example using the
apply()
method:original = data.clone() data.apply(ExpressionSelectionModifier(expression="Position.Z < 0")) data.apply(DeleteSelectedModifier()) print("Number of atoms before:", original.particles.count) print("Number of atoms after:", data.particles.count)
Note that the
clone()
method performs an inexpensive, shallow copy, meaning that the newly created collection will still share the data objects with the original collection. Data objects that are shared by two or more data collections are protected against modification by default to avoid unwanted side effects. Thus, in order to subsequently modify the data objects in either the original collection or its copy, you will have to use the underscore notation or theDataObject.make_mutable()
method to explicitly request a deep copy of the particular data object(s) you want to modify. For example:copy = data.clone() # Data objects are shared by original and copy: assert(copy.cell is data.cell) # In order to modify the SimulationCell in the dataset copy, we must request # a mutable version of the SimulationCell using the 'cell_' accessor: copy.cell_.pbc = (False, False, False) # As a result, the cell object in the second data collection has been replaced # with a deep copy and the two data collections no longer share the same # simulation cell object: assert(copy.cell is not data.cell)
-
property
dislocations
¶ Returns the
DislocationNetwork
data object from this data collection orNone
if there isn’t one.
-
property
grids
¶ Returns a dictionary view providing key-based access to all
VoxelGrids
in this data collection. EachVoxelGrid
has a uniqueidentifier
key, which you can use to look it up in this dictionary. To find out which voxel grids exist in the data collection and what their identifier keys are, useprint(data.grids)
Then retrieve the desired
VoxelGrid
from the collection using its identifier key, e.g.charge_density_grid = data.grids['charge-density'] print(charge_density_grid.shape)
-
property
objects
¶ The unordered list of all
DataObjects
stored in this data collection. You can add or remove data objects in this list as needed.Note that typically you don’t have to work with this list directly, because the
DataCollection
class provides several convenience accessor fields for the different flavors of objects in this mixed list. For example, theDataCollection.particles
field returns theParticles
object from this data objects list. Also, dictionary views such asDataCollection.tables
provide key-based access to a particular class of data objects from this list.
-
property
particles
¶ Returns the
Particles
object from this data collection, which stores the particle properties and -as a sub-object- thebonds
between particles.None
is returned if the data collection does not contain any particle data.Note that the
Particles
object may be marked as read-only if it is currently shared by several data collections. If you intend to modify the particles container or its sub-objects in any way, e.g. add, remove or modify particle properties, you must request a mutable version of the particles object using theparticles_
accessor instead.
-
property
surfaces
¶ Returns a dictionary view providing key-based access to all
SurfaceMesh
objects in this data collection. EachSurfaceMesh
has a uniqueidentifier
key, which can be used to look it up in the dictionary. See the documentation of the modifier producing the surface mesh to find out what the right key is, or useprint(data.surfaces)
to see which identifier keys exist. Then retrieve the desired
SurfaceMesh
object from the collection using its identifier key, e.g.surface = data.surfaces['surface'] print(surface.get_vertices())
-
property
tables
¶ A dictionary view of all
DataTable
objects in this data collection. EachDataTable
has a uniqueidentifier
key, which can be used to look it up in this dictionary. You can useprint(data.tables)
to find out which table identifiers exist in the dictionary. Modifiers that generate a data table typically assign a predefined identifier, which can be found in their documentation. Use the key string to retrieve the desired
DataTable
from the dictionary, e.g.rdf = data.tables['coordination-rdf'] print(rdf.xy())
-
property
trajectories
¶ Returns the
TrajectoryLines
object, which holds the continuous particle trajectories traced by theGenerateTrajectoryLinesModifier
.None
is returned if the data collection does not contain aTrajectoryLines
object.
-
-
class
ovito.data.
DataObject
¶ Abstract base class for all data object types in OVITO.
A
DataObject
represents a fragment of data processed in or by a data pipeline. See theovito.data
module for a list of different concrete data object types in OVITO. Data objects are typically contained in aDataCollection
, which represents a whole data set. Furthermore, data objects can be nested into hierarchy. For example, theBonds
data object is part of the parentParticles
data object.Data objects by themselves are non-visual objects. Visualizing the information stored in a data object in images is the responsibility of so-called visual elements. A data object may be associated with a
DataVis
element by assigning it to the data object’svis
field. Each type of visual element exposes a set of parameters that allow you to configure the appearance of the data visualization in rendered images and animations.-
property
identifier
¶ The unique identifier string of the data object. It serves as look-up key in object dictionaries, for example the
DataCollection.tables
collection, or as a target name in various places where a data object needs to be referenced by name, e.g. theTimeAveragingModifier.operate_on
field.Data objects generated by modifiers in a pipeline typically have an automatically assigned identifier, as documented in the description of the respective modifier. When writing your own modifier function, you are responsible for giving new data objects created by your modifier function a meaningful identifier, so that subsequent modifiers in the pipeline can refer to these data objects.
-
make_mutable
(subobj)¶ This helper method requests a deep copy of subobj, which must be a child
DataObject
of this parentDataObject
. A copy will only be made in case the sub-object is currently referenced by at least one more parent object. If, however, the sub-object is exclusive owned by thisDataObject
, no copy is made and the original sub-object is returned as is. The returned object is safe to modify without unexpected side effects, because any shared ownership is converted an exclusive ownership by the method.Please see the section Announcing object modification for a discussion of object ownership and typical use-cases for this method.
- Parameters
subobj (DataObject) – A existing sub-object of this parent data object, for which exclusive ownership is requested.
- Returns
A copy of subobj if its ownership was previously shared with some other parent. Otherwise the original object is returned.
-
property
vis
¶ The
DataVis
element currently associated with this data object, which is responsible for visually rendering the stored data. If set toNone
, the data object remains non-visual and doesn’t appear in rendered images or the viewports. Furthermore, note that the sameDataVis
element may be assigned to multiple data objects in order to synchronize their visual appearance.
-
property
-
class
ovito.data.
DataTable
¶ Base:
ovito.data.PropertyContainer
This data object type represents a series of data points and is used for generating histogram plots and other 2d data graphs. A data table consists of an array of y-values and, optionally, an array of corresponding x-values, one for each data point. Both arrays are standard
Property
objects managed by the data table, which is aPropertyContainer
subclass.If the
x
data array is not present, the x-coordinates of the data points are implicitly determined by the table’sinterval
, which specifies a value range along the x-axis over which the data points are evenly distributed. This is used, for example, for histograms having equisized bins.Data tables generated by modifiers such as
CoordinationAnalysisModifier
andHistogramModifier
are accessible through theDataCollection.tables
dictionary of the data collection. Please see the documentation of a modifier to find out what data table(s) it produces.Examples:
The following code examples demonstrate how to create a new
DataTable
and fill it with data values in Python. You can use these techniques, for instance, to write custom modifier functions that output their results as data plots in OVITO.To create a simple x-y scatter point plot:
# Create a DataTable object, setting the title and plot type: table = DataTable(title='My Scatter Plot', plot_mode=DataTable.PlotMode.Scatter) # Set the x- and y-coordinates of the data points: table.x = table.create_property('X coordinates', data=numpy.linspace(0.0, 10.0, 50)) table.y = table.create_property('Y coordinates', data=numpy.cos(table.x)) # Add the DataTable to the output DataCollection: data.objects.append(table)
Note how the
create_property()
method is used here to create twoProperty
objects holding the x- and y-coordinates of the data points. These property objects are then set asx
andy
arrays of theDataTable
.Generating a graph with several line plots requires creating a vector property object for the
y
array of theDataTable
:table = DataTable(title='Trigonometric functions', plot_mode=DataTable.PlotMode.Line) table.x = table.create_property('Parameter x', data=numpy.linspace(0.0, 14.0, 100)) # Use the x-coords to compute two y-coords per data point: y(x) = (cos(x), sin(x)) y1y2 = numpy.stack((numpy.cos(table.x), numpy.sin(table.x)), axis=1) table.y = table.create_property('f(x)', data=y1y2, components=['cos(x)', 'sin(x)']) data.objects.append(table)
To generate a bar chart, the
x
property must be filled with numeric IDs 0,1,2,3,… denoting the individual bars. Each bar is then assigned a text label by creating aElementType
instance and adding it to theProperty.types
list:table = DataTable(title='My Bar Chart', plot_mode=DataTable.PlotMode.BarChart) table.x = table.create_property('Structure Type', data=[0, 1, 2, 3]) table.x.types.append(ElementType(id=0, name='Other')) table.x.types.append(ElementType(id=1, name='FCC')) table.x.types.append(ElementType(id=2, name='HCP')) table.x.types.append(ElementType(id=3, name='BCC')) table.y = table.create_property('Count', data=[65, 97, 10, 75]) data.objects.append(table)
For histogram plots, one can specify the complete range of values covered by the histogram by setting the table’s
interval
property. The bin counts must be assigned to the table’sy
property. The number of elements in they
property array, together with theinterval
, determine the number of histogram bins and their width:table = DataTable(title='My Histogram', plot_mode=DataTable.PlotMode.Histogram) table.y = table.create_property('Counts', data=[65, 97, 10, 75]) table.interval = (0.0, 2.0) # Four histogram bins of width 0.5 each. table.axis_label_x = 'Values' # Set the x-axis label of the plot. data.objects.append(table)
If you are going to refer to the data table after it has been inserted into the
DataCollection.objects
list, you should give it a uniqueidentifier
at construction time, as shown in the following example:def modify(frame, data): table = DataTable(identifier='trig-func', title='My Plot', plot_mode=DataTable.PlotMode.Line) table.x = table.create_property('X coords', data=numpy.linspace(0.0, 10.0, 50)) table.y = table.create_property('Y coords', data=numpy.cos(frame * table.x)) data.objects.append(table) pipeline.modifiers.append(modify) export_file(pipeline, 'output/data.*.txt', 'txt/table', key='trig-func', multiple_frames=True)
Here, the user-defined modifier function outputs a table to the data pipeline having the
identifier
'trig-func'
. Whenexport_file()
is called later with the export file format'txt/table'
, the table identifier must be specified as extra argumentkey
to select the right data table to export. Furthermore, the table’s identifier may be used as lookup key to retrieve the table from thetables
dictionary view of aDataCollection
.-
property
axis_label_x
¶ The text label of the x-axis. This string is only used for a data plot if the
x
property of the data table isNone
and the x-coordinates of the data points are implicitly defined by the table’sinterval
property. Otherwise thename
of thex
property is used as axis label.
-
property
interval
¶ A pair of float values specifying the x-axis interval covered by the data points in this table. This interval is only used by the table if the data points do not possess explicit x-coordinates (i.e. if the table’s
x
property isNone
). In the absence of explicit x-coordinates, the interval specifies the range of equispaced x-coordinates implicitly generated by the data table.Implicit x-coordinates are typically used in data tables representing histograms, which consist of equally-sized bins covering a certain value range along the x-axis. The bin size is then given by the interval width divided by the number of data points (see
PropertyContainer.count
property). The implicit x-coordinates of data points are placed in the centers of the bins. You can call the table’sxy()
method to let it explicitly calculate the x-coordinates from the value interval for every data point.- Default
(0.0, 0.0)
-
property
plot_mode
¶ The type of graphical plot for rendering the data in this
DataTable
. Must be one of the following predefined constants:DataTable.PlotMode.NoPlot
DataTable.PlotMode.Line
DataTable.PlotMode.Histogram
DataTable.PlotMode.BarChart
DataTable.PlotMode.Scatter
- Default
DataTable.PlotMode.Line
-
property
x
¶ The
Property
containing the x-coordinates of the data points. The data points may not have explicit x-coordinates, so this property may beNone
for a data table. In such a case, the x-coordinates of the data points are implicitly determined by the table’sinterval
.
-
property
-
class
ovito.data.
DislocationNetwork
¶ Base:
ovito.data.DataObject
This data object stores the network of dislocation lines extracted by a
DislocationAnalysisModifier
. You can access it through theDataCollection.dislocations
field.The dislocation network is associated with a
DislocationVis
element controlling the visual appearance of the dislocation lines. It can be accessed through thevis
attribute of theDataObject
base class.Example:
from ovito.io import import_file, export_file from ovito.modifiers import DislocationAnalysisModifier from ovito.data import DislocationNetwork pipeline = import_file("input/simulation.dump") # Extract dislocation lines from a crystal with diamond structure: modifier = DislocationAnalysisModifier() modifier.input_crystal_structure = DislocationAnalysisModifier.Lattice.CubicDiamond pipeline.modifiers.append(modifier) data = pipeline.compute() total_line_length = data.attributes['DislocationAnalysis.total_line_length'] cell_volume = data.attributes['DislocationAnalysis.cell_volume'] print("Dislocation density: %f" % (total_line_length / cell_volume)) # Print list of dislocation lines: print("Found %i dislocation segments" % len(data.dislocations.segments)) for segment in data.dislocations.segments: print("Segment %i: length=%f, Burgers vector=%s" % (segment.id, segment.length, segment.true_burgers_vector)) print(segment.points) # Export dislocation lines to a CA file: export_file(pipeline, "output/dislocations.ca", "ca") # Or export dislocations to a ParaView VTK file: export_file(pipeline, "output/dislocations.vtk", "vtk/disloc")
File export
A dislocation network can be written to a data file in the form of polylines using the
ovito.io.export_file()
function (select thevtk/disloc
output format). During export, a non-periodic version is produced by clipping dislocation lines at the domain boundaries.-
property
segments
¶ The list of dislocation segments in this dislocation network. This list-like object is read-only and contains
DislocationSegment
objects.
-
set_segment
(index, true_burgers_vector=None, cluster_id=None, points=None, custom_color=None)¶ This method allows you to change the data fields of individual dislocation lines. Fields for which no new value is specified will keep their current values.
- Parameters
index – The zero-based index of the dislocation line in the
segments
array to be modified.true_burgers_vector – The lattice-space Burgers vector (
true_burgers_vector
) to be assigned to the dislocation line.cluster_id – The numeric ID of the crystallite cluster the dislocation line is embedded in.
points – A N x 3 NumPy array with the Cartesian coordinates of the dislocation line vertices.
custom_color – RGB color to be used for rendering the line instead of the automatically determined color.
Example of a user-defined modifier function manipulating the dislocation line data:
import numpy as np def modify(frame, data): # Flip Burgers vector and line sense of each dislocation: for index, seg in enumerate(data.dislocations.segments): data.dislocations_.set_segment(index, true_burgers_vector = np.negative(seg.true_burgers_vector), points = np.flipud(seg.points)) # Highlight all 1/6[121] dislocations in a red color: for index, seg in enumerate(data.dislocations.segments): if np.allclose(seg.true_burgers_vector, (1/6, 2/6, 1/6)): data.dislocations_.set_segment(index, custom_color=(1, 0, 0))
-
property
-
class
ovito.data.
DislocationSegment
¶ A single dislocation line from a
DislocationNetwork
.The list of dislocation segments is returned by the
DislocationNetwork.segments
attribute.-
property
cluster_id
¶ The numeric identifier of the crystal cluster of atoms containing this dislocation segment.
The true Burgers vector of the segment is expressed in the local coordinate system of this crystal cluster.
-
property
custom_color
¶ The RGB color value to be used for visualizing this particular dislocation line, overriding the default coloring scheme imposed by the
DislocationVis.indicate_character
mode setting. The custom color is only used if its RGB components are non-negative (i.e. in the range 0-1); otherwise the line will be rendered using the computed color depending on the line’s Burgers vector.- Default
(-1.0, -1.0, -1.0)
-
property
id
¶ The unique identifier of this dislocation segment.
-
property
is_infinite_line
¶ This property indicates whether this segment is an infinite line passing through a periodic simulation box boundary. A segment is considered infinite if it is a closed loop and its start and end points do not coincide.
See also the
is_loop
property.
-
property
is_loop
¶ This property indicates whether this segment forms a closed dislocation loop. Note that an infinite dislocation line passing through a periodic boundary is also considered a loop.
See also the
is_infinite_line
property.
-
property
length
¶ Returns the length of this dislocation segment.
-
property
points
¶ The list of space points that define the shape of this dislocation segment. This is a N x 3 Numpy array, where N is the number of points along the segment. For closed loops, the first and the last point coincide.
-
property
spatial_burgers_vector
¶ The Burgers vector of the segment, expressed in the global coordinate system of the simulation. This vector is calculated by transforming the true Burgers vector from the local lattice coordinate system to the global simulation coordinate system using the average orientation matrix of the crystal cluster the dislocation segment is embedded in.
-
property
true_burgers_vector
¶ The Burgers vector of the segment, expressed in the local coordinate system of the crystal. Also known as the True Burgers vector.
-
property
-
class
ovito.data.
ElementType
¶ Base:
ovito.data.DataObject
This class represents a single type of data elements, for example a particular atom or bond type. It serves as common base class for the
ParticleType
andBondType
classes, which represent these more specific types.Each type has a unique numeric
id
, which is used when looking up types given some numeric value in a typedProperty
array. An example for a typed property is the particle property'Particle Type'
. The property object manages all its associated types in itsProperty.types
list.The
Property.type_by_id()
andProperty.type_by_name()
methods allow looking up a certainElementType
based on its numeric identifier or name string.-
property
color
¶ The color used when rendering elements of this type. This is a RGB tuple with components in the range 0.0 – 1.0.
- Default
(1.0, 1.0, 1.0)
-
property
enabled
¶ Controls whether this type is currently active or inactive. This flag currently has a meaning only in the context of atomic structure identification. Some analysis modifiers manage a list of the structure types they can identify (e.g. FCC, BCC, etc.). The identification of individual structure types can be turned on or off by the user by changing their
enabled
flag. SeeStructureIdentificationModifier.structures
for further information.- Default
True
-
property
id
¶ The unique numeric identifier of the type (typically some positive
int
). The identifier is and must be unique among all element types in thetypes
list of a typedProperty
. Thus, if you create a new element type, make sure you give it a unique id before inserting it into thetypes
list of a typed property.- Default
0
-
property
-
class
ovito.data.
NearestNeighborFinder
(N, data_collection)¶ A utility class that finds the N nearest neighbors of a particle or around a spatial location.
The constructor takes the (maximum) number of requested nearest neighbors, N, and a
DataCollection
containing the input particles and the cell geometry (including periodic boundary flags). N must be a positive integer not greater than 30 (which is the built-in maximum supported by this class).Once the
NearestNeighborFinder
has been constructed, you can call itsfind()
method to iterate over the sorted list of nearest neighbors of a specific particle, for example:from ovito.io import import_file from ovito.data import NearestNeighborFinder # Load input simulation file. pipeline = import_file("input/simulation.dump") data = pipeline.compute() # Initialize neighbor finder object. # Visit the 12 nearest neighbors of each particle. N = 12 finder = NearestNeighborFinder(N, data) # Prefetch the property array containing the particle type information: ptypes = data.particles.particle_types # Loop over all input particles: for index in range(data.particles.count): print("Nearest neighbors of particle %i:" % index) # Iterate over the neighbors of the current particle, starting with the closest: for neigh in finder.find(index): print(neigh.index, neigh.distance, neigh.delta) # The index can be used to access properties of the current neighbor, e.g. type_of_neighbor = ptypes[neigh.index]
Furthermore, the class offers the
find_at()
method, which lets you determine the N nearest particles around an arbitrary spatial location:# Find particles closest to some spatial point (x,y,z): coords = (0, 0, 0) for neigh in finder.find_at(coords): print(neigh.index, neigh.distance, neigh.delta)
Note: In case you rather want to find all neighbor particles within a certain cutoff range of a particle, use the
CutoffNeighborFinder
class instead.-
find
(index)¶ Returns an iterator that visits the N nearest neighbors of the given particle in order of ascending distance.
- Parameters
index (int) – The index of the central particle whose neighbors should be iterated. Particle indices start at 0.
- Returns
A Python iterator that visits the N nearest neighbors of the central particle in order of ascending distance. For each visited neighbor the iterator returns an object with the following attributes:
index: The global index of the current neighbor particle.
distance: The distance of the current neighbor from the central particle.
distance_squared: The squared neighbor distance.
delta: The three-dimensional vector connecting the central particle with the current neighbor (correctly taking into account periodic boundary conditions).
The global index returned by the iterator can be used to look up properties of the neighbor as demonstrated in the first example code above.
Note that several periodic images of the same particle may be visited if the periodic simulation cell is sufficiently small. Then the same particle index may appear more than once in the neighbor list. In fact, the central particle may be among its own neighbors in a sufficiently small periodic simulation cell. However, the computed neighbor vector (
delta
) will be unique for each image of a neighboring particle.The number of neighbors actually visited may be smaller than the requested number, N, if the system contains too few particles and has no periodic boundary conditions.
Note that the
find()
method will not find other particles located exactly at the same spatial position as the central particle for technical reasons. To find such particles too, which are positioned exactly on top of each other, make use of thefind_at()
method instead.
-
find_at
(coords)¶ Returns an iterator that visits the N nearest particles around a spatial point given by coords in order of ascending distance. Unlike the
find()
method, which queries the nearest neighbors of a physical particle, thefind_at()
method allows searching for nearby particles at arbitrary locations in space.- Parameters
coords – A (x,y,z) coordinate triplet specifying the spatial location where the N nearest particles should be queried.
- Returns
A Python iterator that visits the N nearest neighbors in order of ascending distance. For each visited particle the iterator returns an object with the following attributes:
index: The index of the current particle (starting at 0).
distance: The distance of the current neighbor from the query location.
distance_squared: The squared distance to the query location.
delta: The three-dimensional vector from the query point to the current particle (correctly taking into account periodic boundary conditions).
If there exists a particle that is exactly located at the query location given by coords, then it will be returned by this function. This is in contrast to the
find()
function, which does not visit the central particle itself.The number of neighbors actually visited may be smaller than the requested number, N, if the system contains too few particles and has no periodic boundary conditions.
-
-
class
ovito.data.
ParticleType
¶ Base:
ovito.data.ElementType
This data object describes one particle or atom type. In atomistic simulations, each chemical element is typically represented by an instance of the
ParticleType
class. The property fields of the class control how the particles of that type get visualized in terms of e.g. color, particle radius, shape, etc.The
ParticleType
class inherits several general data fields from its base classElementType
, e.g. thecolor
,name
andid
fields. It adds specific fields for particles:radius
andshape
. Furthermore, the class has additional fields controlling the visual appearance of particles with user-defined shapes.The
ParticleType
instances are all stored in theProperty
object with the name'Particle Type'
, which also stores for each particle what its type is. The association between particles and particle types is established via a unique numeric typeid
. The following code shows how to iterate over all particle types in a dataset, which are listed in theProperty.types
field of theparticle_types
property:# Access the property with the name 'Particle Type': prop = data.particles.particle_types # Print list of particle types (their numeric IDs and names) for t in prop.types: print('ID {} -> {}'.format(t.id, t.name)) # Print the numeric type ID of each particle: for tid in prop[...]: print(tid)
The order in which the particle types are stored in the
Property.types
list is arbitrary, and the unique numeric IDs of particle types have no specific meaning in general. A common operation is to find theParticleType
in the list corresponding to a given numeric ID. For this look up operation, theProperty
class provides thetype_by_id()
method:# Look up the particle type with unique ID 2: t = prop.type_by_id(2) print(t.name, t.color, t.radius) # Iterate over all particles and print their type's name: for tid in prop[...]: print(prop.type_by_id(tid).name)
Another common operation is to look up a particle type by name, for example the type representing a certain chemical element. For this kind of look up operation, the
type_by_name()
method may be used, which assumes that each type has a uniquename
(which may not always be true):# Print numeric ID of particle type 'Si': print(prop.type_by_name('Si').id)
-
property
backface_culling
¶ Activates back-face culling for the user-defined particle shape mesh to speed up rendering. If turned on, polygonal sides of the shape mesh facing away from the viewer will not be rendered. You can turn this option off if the particle’s shape is not closed and two-sided rendering is required. This option only has an effect if a user-defined shape has been assigned to the particle type using the
load_shape()
method.- Default
True
-
property
highlight_edges
¶ Activates the highlighting of the polygonal edges of the user-defined particle shape during rendering. This option only has an effect if a user-defined shape has been assigned to the particle type using the
load_shape()
method.- Default
False
-
load_defaults
()¶ Given the type’s chemical
name
, which must have been set before calling this method, initializes the type’scolor
andradius
fields with default values from OVITO’s internal database of chemical elements. This method is useful when creating new atom types while building up a molecule structure.
-
load_shape
(filepath)¶ Assigns a user-defined shape to the particle type. Particles of this type will subsequently be rendered using the polyhedral mesh loaded from the given file. The method will automatically detect the format of the geometry file and supports standard file formats such as OBJ, STL and VTK that contain triangle meshes, see this table.
The shape loaded from the geometry file will be scaled with the
radius
value set for this particle type or the per-particle value stored in theRadius
particle property if present. The shape of each particle will be rendered such that its origin is located at the coordinates of the particle.The following example script demonstrates how to load a user-defined shape for the first particle type (index 0) loaded from a LAMMPS dump file, which can be accessed through the
Property.types
list of theParticle Type
particle property.pipeline = import_file("input/simulation.dump") pipeline.add_to_scene() types = pipeline.source.data.particles_.particle_types_ types.type_by_id_(1).load_shape("input/tetrahedron.vtk") types.type_by_id_(1).highlight_edges = True
-
property
mass
¶ The mass of this particle type.
- Default
0.0
-
property
radius
¶ This property controls the display radius of the particles of this type.
When set to zero, particles of this type will be rendered using the standard size specified by the
ParticlesVis.radius
parameter. Furthermore, precedence is given to any per-particle sizes assigned to theRadius
particle property if that property has been defined.- Default
0.0
The following example script demonstrates how to set the display radii of two particle types loaded from a simulation file, which can be accessed through the
Property.types
list of theParticle Type
particle property.pipeline = import_file("input/simulation.dump") pipeline.add_to_scene() def setup_particle_types(frame, data): types = data.particles_.particle_types_ types.type_by_id_(1).name = "Cu" types.type_by_id_(1).radius = 1.35 types.type_by_id_(2).name = "Zr" types.type_by_id_(2).radius = 1.55 pipeline.modifiers.append(setup_particle_types)
-
property
shape
¶ Selects the geometric shape used when rendering particles of this type. Supported modes are:
ParticlesVis.Shape.Unspecified
(default)ParticlesVis.Shape.Sphere
ParticlesVis.Shape.Box
ParticlesVis.Shape.Circle
ParticlesVis.Shape.Square
ParticlesVis.Shape.Cylinder
ParticlesVis.Shape.Spherocylinder
ParticlesVis.Shape.Mesh
By default, the standard particle shape that is set in the
ParticlesVis
visual element is used to render particles of this type. Parameter values other thanUnspecified
allow you to control the rendering shape on a per-type basis. ModeSphere
includes ellipsoid and superquadric particle shapes, which are enabled by the presence of theAspherical Shape
andSuperquadric Roundness
particle properties.The
load_shape()
method lets you specify a user-defined mesh geometry for this particle type. Calling this method automatically switches the shape parameter to modeMesh
.Setting the shapes of particle types permanently, i.e., for all frames of a loaded simulation trajectory, typically requires a user-defined modifier function. This function is inserted into the
Pipeline
to make the necessary changes to theParticleType
objects associated with theProperty
namedParticle Type
:from ovito.io import import_file from ovito.vis import * # Load a simulation file containing numeric particle types 1, 2, 3, ... pipeline = import_file("input/nylon.data") pipeline.add_to_scene() # Set the default particle shape in the ParticlesVis visual element, # which will be used by all particle types for which we do not specify a different shape below. pipeline.compute().particles.vis.shape = ParticlesVis.Shape.Box pipeline.compute().particles.vis.radius = 1.0 # A user-defined modifier function that configures the shapes of particle types 1 and 2: def setup_particle_types(frame, data): # Write access to property 'Particle Type': types = data.particles_.particle_types_ # Write access to numeric ParticleTypes, which are sub-objects of the Property object: types.type_by_id_(1).radius = 0.5 types.type_by_id_(1).shape = ParticlesVis.Shape.Cylinder types.type_by_id_(2).radius = 1.2 types.type_by_id_(2).shape = ParticlesVis.Shape.Sphere pipeline.modifiers.append(setup_particle_types) # Render a picture of the 3d scene: vp = Viewport(camera_dir = (-2,1,-1)) vp.zoom_all() vp.render_image(filename='output/particles.png', size=(320,240), renderer=TachyonRenderer())
-
property
use_mesh_color
¶ Use the intrinsic mesh color(s) instead of the particle color when rendering particles of this type. This option only has an effect if a user-defined shape has been assigned to the particle type using the
load_shape()
method.- Default
False
-
property
vdw_radius
¶ The van der Waals radius of the particle type. This value is used by the
CreateBondsModifier
to decide which pairs of particles are close enough to be connected by a bond. In contrast to theradius
parameter, the van der Waals radius does not affect the visual appearance of the particles of this type.- Default
0.0
-
property
-
class
ovito.data.
Particles
¶ Base:
ovito.data.PropertyContainer
This container object stores the information associated with a system of particles. It is typically accessed through the
DataCollection.particles
field of a data collection. The current number of particles is given by thecount
attribute that is inherited from thePropertyContainer
base class. The particles may be associated with a set of properties. Each property is represented by aProperty
data object, that is stored in this property container and is basically an array of numeric values of length N, where N is the number of particles in the system. Each property array has a unique name, by which it can be looked up through the dictionary interface of thePropertyContainer
base class. While the user is free to define arbitrary particle properties, OVITO predefines a set of standard properties that each have a fixed data layout, meaning and name. They are listed in the table below.Standard property name
Data type
Component names
Angular Momentum
float
X, Y, Z
Angular Velocity
float
X, Y, Z
Aspherical Shape
float
X, Y, Z
Centrosymmetry
float
Charge
float
Cluster
int64
Color
float
R, G, B
Coordination
int
Deformation Gradient
float
XX, YX, ZX, XY, YY, ZY, XZ, YZ, ZZ
Dipole Magnitude
float
Dipole Orientation
float
X, Y, Z
Displacement Magnitude
float
Displacement
float
X, Y, Z
DNA Strand
int
Elastic Deformation Gradient
float
XX, YX, ZX, XY, YY, ZY, XZ, YZ, ZZ
Elastic Strain
float
XX, YY, ZZ, XY, XZ, YZ
Force
float
X, Y, Z
Kinetic Energy
float
Mass
float
Molecule Identifier
int64
Molecule Type
int
Nucleobase
int
Nucleotide Axis
float
X, Y, Z
Nucleotide Normal
float
X, Y, Z
Orientation
float
X, Y, Z, W
Particle Identifier
int64
Particle Type
int
Periodic Image
int
X, Y, Z
Position
float
X, Y, Z
Potential Energy
float
Radius
float
Rotation
float
X, Y, Z, W
Selection
int
Spin
float
Strain Tensor
float
XX, YY, ZZ, XY, XZ, YZ
Stress Tensor
float
XX, YY, ZZ, XY, XZ, YZ
Stretch Tensor
float
XX, YY, ZZ, XY, XZ, YZ
Structure Type
int
Superquadric Roundness
float
Phi, Theta
Torque
float
X, Y, Z
Total Energy
float
Transparency
float
Vector Color
float
R, G, B
Velocity Magnitude
float
Velocity
float
X, Y, Z
-
property
angles
¶ A
PropertyContainer
storing the list of angles defined for the molecular model (may beNone
).
-
property
colors
¶ The
Property
data array for theColor
standard particle property; orNone
if that property is undefined.
-
create_particle
(position)¶ Adds a new particle to the particle system.
- Parameters
position (tuple) – The xyz coordinates for the new particle.
- Returns
The index of the newly created particle, i.e.
(Particles.count-1)
.
-
delta_vector
(a, b, cell, return_pbcvec=False)¶ Computes the vector connecting two particles a and b in a periodic simulation cell by applying the minimum image convention.
This is a convenience wrapper for the
SimulationCell.delta_vector()
method, which computes the vector between two arbitrary spatial locations \(r_a\) and \(r_b\) taking into account periodic boundary conditions. The version of the method described here takes two particle indices a and b as input, computing the shortest vector \({\Delta} = (r_b - r_a)\) between them using the minimum image convention. Please see theSimulationCell.delta_vector()
method for further information.- Parameters
a – Zero-based index of the first input particle. This may also be an array of particle indices.
b – Zero-based index of the second input particle. This may also be an array of particle indices with the same length as a.
cell – The
SimulationCell
object that defines the periodic domain. Typically,DataCollection.cell
is used here.return_pbcvec (bool) – If True, also returns the vector \(n\), which specifies how often the computed particle-to-particle vector crosses the cell’s face.
- Returns
The delta vector and, optionally, the vector \(n\).
-
property
dihedrals
¶ A
PropertyContainer
storing the list of dihedrals defined for the molecular model (may beNone
).
-
property
forces
¶ The
Property
data array for theForce
standard particle property; orNone
if that property is undefined.
-
property
identifiers
¶ The
Property
data array for theParticle Identifier
standard particle property; orNone
if that property is undefined.
-
property
impropers
¶ A
PropertyContainer
storing the list of impropers defined for the molecular model (may beNone
).
-
property
masses
¶ The
Property
data array for theMass
standard particle property; orNone
if that property is undefined.
-
property
particle_types
¶ The
Property
data array for theParticle Type
standard particle property; orNone
if that property is undefined.
-
property
positions
¶ The
Property
data array for thePosition
standard particle property; orNone
if that property is undefined.
-
property
selection
¶ The
Property
data array for theSelection
standard particle property; orNone
if that property is undefined.
-
property
-
class
ovito.data.
Property
¶ Base:
ovito.data.DataObject
Stores the property values for an array of data elements (e.g. particles, bonds or voxels).
Each property, for example a particle property, is represented by one
Property
object storing the property values for all particles. Thus, aProperty
object is basically an array of values whose length matches the number of data elements.All
Property
objects belonging to the same class of data elements, for example all particle properties, are managed by aPropertyContainer
. In the case of particle properties, the corresponding container class is theParticles
class, which is a specialization of the genericPropertyContainer
base class.Data access
A
Property
object behaves like a Numpy array. For example, you can access the property value for the i-th data element using indexing:positions = data.particles['Position'] print('Position of first particle:', positions[0]) print('Z-coordinate of second particle:', positions[1,2]) for xyz in positions: print(xyz)
Element indices start at zero. Properties can be either vectorial (e.g. velocity vectors are stored as an N x 3 array) or scalar (1-d array of length N). The length of the first array dimension is in both cases equal to the number of data elements (number of particles in the example above). Array elements can either be of data type
float64
,int32
orint64
.If you want to modify the per-element values in a property array, make sure you are working with a modifiable version of the array by employing the underscore notation, e.g.:
modifiable_positions = data.particles_['Position_'] modifiable_positions[0] = (0.1, 2.3, 0.5)
Typed properties
The standard particle property
'Particle Type'
stores the types of particles encoded as integer values, e.g.:>>> data = node.compute() >>> tprop = data.particles['Particle Type'] >>> print(tprop[...]) [2 1 3 ..., 2 1 2]
Here, each number in the property array refers to one of the particle types (e.g. 1=Cu, 2=Ni, 3=Fe, etc.). The defined particle types, each one represented by an instance of the
ParticleType
auxiliary class, are stored in thetypes
array of theProperty
. Each type has a uniqueid
, a human-readablename
and other attributes likecolor
andradius
that control the visual appearance of particles belonging to the type:>>> for type in tprop.types: ... print(type.id, type.name, type.color, type.radius) ... 1 Cu (0.188 0.313 0.972) 0.74 2 Ni (0.564 0.564 0.564) 0.77 3 Fe (1 0.050 0.050) 0.74
IDs of types typically start at 1 and form a consecutive sequence as in the example above. Note, however, that the
types
list may store theParticleType
objects in an arbitrary order. Thus, in general, it is not valid to directly use a type ID as an index into thetypes
array. Instead, thetype_by_id()
method should be used to look up theParticleType
:>>> for i,t in enumerate(tprop): # (loop over the type ID of each particle) ... print('Atom', i, 'is of type', tprop.type_by_id(t).name) ... Atom 0 is of type Ni Atom 1 is of type Cu Atom 2 is of type Fe Atom 3 is of type Cu
Similarly, a
type_by_name()
method exists that looks up aParticleType
by name. For example, to count the number of Fe atoms in a system:>>> Fe_type_id = tprop.type_by_name('Fe').id # Determine ID of the 'Fe' type >>> numpy.count_nonzero(tprop == Fe_type_id) # Count particles having that type ID 957
Note that OVITO supports multiple type classifications. For example, in addition to the
'Particle Type'
standard particle property, which stores the chemical types of atoms (e.g. C, H, Fe, …), the'Structure Type'
property may hold the structural types computed for atoms (e.g. FCC, BCC, …) maintaining its own list of known structure types in thetypes
array.-
property
component_count
¶ The number of vector components if this is a vector property; or 1 if this is a scalar property.
-
property
component_names
¶ The list of names of the vector components if this is a vector property. For example, for the
Position
particle property this field contains['X', 'Y', 'Z']
.
-
property
name
¶ The name of the property.
-
type_by_id
(id)¶ Looks up the
ElementType
with the given numeric ID in thetypes
list. Raises aKeyError
if no type with the numeric ID exists.
-
type_by_name
(name)¶ Looks up the
ElementType
with the given name in the property’stypes
list. If multiple types exists with the same name, the first type is returned. Raises aKeyError
if there is no type with such a name.
-
property
types
¶ The list of
ElementType
instances attached to this property.Note that the element types may be stored in arbitrary order in this list. Thus, it is not valid to use a numeric type ID as an index into this list.
-
property
-
class
ovito.data.
PropertyContainer
¶ Base:
ovito.data.DataObject
A dictionary-like object storing a set of
Property
objects.It implements the
collections.abc.Mapping
interface. That means it can be used like a standard read-only Pythondict
object to access the properties by name, e.g.:data = pipeline.compute() positions = data.particles['Position'] has_selection = 'Selection' in data.particles name_list = data.particles.keys()
New properties can be added with the
create_property()
method as described here.OVITO provides several concrete implementations of the abstract
PropertyContainer
base class:-
property
count
¶ The number of data elements in this container, e.g. the number of particles. This value is always equal to the lengths of the
Property
arrays managed by this container.
-
create_property
(name, dtype=None, components=None, data=None)¶ Adds a new property to the container and optionally initializes it with the per-element data provided by the data parameter. The method returns the new
Property
instance.The method allows to create standard as well as user-defined properties. To create a standard property, one of the standard property names must be provided as name argument:
colors = numpy.random.random_sample(size = (data.particles.count, 3)) data.particles_.create_property('Color', data=colors)
The length of the provided data array must match the number of existing elements in the container, which is given by the
count
attribute. You can alternatively assign the per-element values to the property after its construction:prop = data.particles_.create_property('Color') prop[...] = numpy.random.random_sample(size = prop.shape)
To create a user-defined property, use a non-standard property name:
values = numpy.arange(0, data.particles.count, dtype=int) data.particles_.create_property('myint', data=values)
In this case the data type and the number of vector components of the new property are inferred from the provided data Numpy array. Providing a one-dimensional array creates a scalar property while a two-dimensional array creates a vectorial property. Alternatively, the dtype and components parameters can be specified explicitly if you are going to assign the property values at a later time:
prop = data.particles_.create_property('myvector', dtype=float, components=3) prop[...] = numpy.random.random_sample(size = prop.shape)
If the property to be created already exists in the container, it is replaced with a new one. The existing per-element data from the old property is however retained if data is
None
.Note: If the container contains no properties yet, then the number of elements (e.g. particles or bonds) is still undefined. In this case the
create_property()
method lets you define the number of elements when inserting the very first property by specifying a data array of the desired length. For example, to create a newParticles
container from scratch with 10 particles, a Numpy array of length 10 is used to initialize thePosition
particle property:# An empty Particles container to begin with: particles = Particles() # Create 10 particles with random xyz coordinates: xyz = numpy.random.random_sample(size = (10,3)) particles.create_property('Position', data=xyz)
After the initial
Positions
property has been created, the number of particles in the container is now determined and any subsequently added properties must have the exact same length.- Parameters
name – Either a standard property type constant or a name string.
data – An optional data array with per-element values for initializing the new property. The size of the array must match the element
count
of the container and the shape must be consistent with the number of components of the property to be created.dtype – The element data type when creating a user-defined property. Must be either
int
orfloat
.components (int) – The number of vector components when creating a user-defined property.
- Returns
The newly created
Property
object.
-
delete_elements
(mask)¶ Deletes a subset of the elements from this container. The elements to be deleted must be specified in terms of a 1-dimensional mask array having the same length as the container (see
count
). The method will delete those elements whose corresponding mask value is non-zero, i.e., thei
-th element will be deleted ifmask[i]!=0
.For example, to delete all currently selected particles, i.e., the subset of particles whose
Selection
property is non-zero, one would simply write:data.particles_.delete_elements(data.particles['Selection'])
The effect of this statement is the same as for applying the
DeleteSelectedModifier
to the particles list.
-
delete_indices
(indices)¶ Deletes a subset of the elements from this container. The elements to be deleted must be specified in terms of a sequence of indices, all in the range 0 to
count
-1. The method accepts any type of iterable object, including sequence types and generators.For example, to delete every other particle, one could use Python’s
range()
function to generate all even indices up to the length of the particle container:data.particles_.delete_indices(range(0, data.particles.count, 2))
-
property
title
¶ The title of the data object under which it appears in the user interface of OVITO.
-
property
-
class
ovito.data.
SimulationCell
¶ Base:
ovito.data.DataObject
Stores the geometric shape and the boundary conditions of the simulation cell. A
SimulationCell
data object is typically part of aDataCollection
and can be retrieved through itscell
field:from ovito.io import import_file pipeline = import_file("input/simulation.dump") data = pipeline.compute() # Print cell matrix to the console. [...] is for casting to Numpy array. print(data.cell[...])
The simulation cell geometry is stored as a 3x4 matrix (with column-major ordering). The first three columns of the matrix represent the three cell vectors and the last column is the position of the cell’s origin. For two-dimensional datasets, the
is2D
flag is set. In this case the third cell vector and the z-coordinate of the cell origin are ignored by OVITO in many computations.# Compute simulation box volume by taking the determinant of the # left 3x3 submatrix of the cell matrix: vol = abs(numpy.linalg.det(data.cell[0:3,0:3])) # The SimulationCell.volume property yields the same value. assert numpy.isclose(data.cell.volume, vol)
The
SimulationCell
object behaves like a standard Numpy array of shape (3,4). When modifying the values of the cell matrix, make sure you use the underscore notation to request a modifiable version of theSimulationCell
object:# Make cell twice as large along the Y direction by scaling the second cell vector: data.cell_[:,1] *= 2.0
A
SimulationCell
is always associated with aSimulationCellVis
element, which controls the visual appearance of the simulation box in rendered images. It can be accessed asvis
attribute of the cell, inherited from theDataObject
base class:# Change display color of simulation cell to red: data.cell.vis.rendering_color = (1.0, 0.0, 0.0)
-
delta_vector
(ra, rb, return_pbcvec=False)¶ Computes the correct vector connecting points \(r_a\) and \(r_b\) in a periodic simulation cell by applying the minimum image convention.
The method starts by computing the 3d vector \({\Delta} = r_b - r_a\) for two points \(r_a\) and \(r_b\), which may be located in different images of the periodic simulation cell. The minimum image convention is then applied to obtain the new vector \({\Delta'} = r_b' - r_a\), where the original point \(r_b\) has been replaced by the periodic image \(r_b'\) that is closest to \(r_a\), making the vector \({\Delta'}\) as short as possible (in reduced coordinate space). \(r_b'\) is obtained by translating \(r_b\) an integer number of times along each of the three cell directions: \(r_b' = r_b - H*n\), with \(H\) being the 3x3 cell matrix and \(n\) being a vector of three integers that are chosen by the method such that \(r_b'\) is as close to \(r_a\) as possible.
Note that the periodic image convention is applied only along those cell directions for which periodic boundary conditions are enabled (see
pbc
property). For other directions no shifting is performed, i.e., the corresponding components of \(n = (n_x,n_y,n_z)\) will always be zero.The method is able to compute the results for either an individual pair of input points or for two arrays of input points. In the latter case, i.e. if the input parameters ra and rb are both 2-D arrays of shape Nx3, the method returns a 2-D array containing N output vectors. This allows applying the minimum image convention to a large number of point pairs in one function call.
The optional return_pbcvec flag makes the method return as an additional output the vector \(n\) introduced above. The components of this vector specify the number of times the image point \(r_b'\) needs to be shifted along each of the three cell directions in order to bring it onto the original input point \(r_b\). In other words, it specifies the number of times the computed vector \({\Delta} = r_b - r_a\) crosses a periodic boundary of the cell (either in positive or negative direction). For example, the PBC shift vector \(n = (1,0,-2)\) would indicate that, in order to get from input point \(r_a\) to input point \(r_b\), one has to cross the cell boundaries once in the positive x-direction and twice in the negative z-direction. If return_pbcvec is True, the method returns the tuple (\({\Delta'}\), \(n\)); otherwise it returns just \({\Delta'}\). Note that the vector \(n\) computed by this method can be used, for instance, to correctly initialize the
Bonds.pbc_vectors
property for newly created bonds that cross a periodic cell boundary.- Parameters
ra – The Cartesian xyz coordinates of the first input point(s). Either a 1-D array of length 3 or a 2-D array of shape (N,3).
rb – The Cartesian xyz coordinates of the second input point(s). Must have the same shape as ra.
return_pbcvec (bool) – If True, also returns the vector \(n\), which specifies how often the vector \((r_b'r - r_a)\) crosses the periodic cell boundaries.
- Returns
The vector \({\Delta'}\) and, optionally, the vector \(n\).
Note that there exists also a convenience method
Particles.delta_vector()
, which should be used in cases where \(r_a\) and \(r_b\) are both positions of particles in the simulation cell.
-
property
is2D
¶ Specifies whether the system is two-dimensional (instead of three-dimensional). For two-dimensional systems, the PBC flag in the third direction (Z) and the third cell vector will typically be ignored.
- Default
False
-
property
pbc
¶ A tuple of three Boolean flags specifying whether periodic boundary conditions are enabled along the cell’s directions.
- Default
(False, False, False)
-
property
volume
¶ Computes the volume of the three-dimensional simulation cell. The volume is the absolute value of the determinant of the 3x3 submatrix formed by the three cell vectors.
-
-
class
ovito.data.
SurfaceMesh
¶ Base:
ovito.data.DataObject
A triangle mesh representation of a surface or, more precisely, a two-dimensional manifold that is closed and orientable. Typically, surface meshes are produced within OVITO by modifiers such as the
ConstructSurfaceModifier
,CreateIsosurfaceModifier
orCoordinationPolyhedraModifier
. Please see the user manual page on surface meshes for a more in-depth description of this geometric data structure.Periodic domains
Surface meshes may be embedded in periodic domains, i.e. in simulation cells with periodic boundary conditions applied. That means triangles of a surface mesh can be connected to vertices on opposite sides of a simulation box and wrap around correctly. OVITO takes care of computing the intersections of the periodic surface with the box boundaries and automatically produces a non-periodic representation of the mesh when it comes to visualization, which is cut off at the box boundaries.
The spatial domain the surface mesh is embedded in is specified as a
SimulationCell
object attached to theSurfaceMesh
. It is stored in theSurfaceMesh.domain
field. Note that thisSimulationCell
may, in some situations, not be identical to the global simulationcell
set for theDataCollection
.Spatial regions
Being a closed, orientable manifold a surface mesh divides space into two or more spatial regions. For example, if the surface mesh was constructed by the
ConstructSurfaceModifier
from a set of input particles, then the space enclosed by the surface is the “filled” region and the exterior space is the “empty” region containing no particles.In general, the
SurfaceMesh
class manages a variable list ofregions
, each being identified by a numeric, zero-based index. Thelocate_point()
method allows to determine which spatial region some point in space belongs to.It may be that a surface mesh is degenerate. In such a case there is only one spatial region filling entire space. For example, when there are no input particles, the
ConstructSurfaceModifier
cannot construct a regular surface mesh and the “empty” region fills the entire simulation cell. Conversely, if the periodic simulation cell is completely filled with particles, the “filled” region covers the entire simulation domain and the generated surface mesh will consist of no vertices or faces, i.e., it is degenerate. To discriminate between the two situations, theSurfaceMesh
class has aspace_filling_region
field, which specifies the spatial region that fills entire space in cases where the mesh is degenerate.File export
A surface mesh can be exported to a geometry file in the form of a triangle mesh using OVITO’s
export_file()
function. To this end, a non-periodic version is produced by truncating triangles at the domain boundaries and generating “cap polygons” filling the holes that occur at the intersection of the surface with periodic domain boundaries. The following example code writes a VTK geometry file (vtk/trimesh
export format):from ovito.io import import_file, export_file from ovito.data import SurfaceMesh from ovito.modifiers import ConstructSurfaceModifier # Load a particle set and construct the surface mesh: pipeline = import_file("input/simulation.dump") pipeline.modifiers.append(ConstructSurfaceModifier(radius = 2.8)) mesh = pipeline.compute().surfaces['surface'] # Export the mesh to a VTK file for visualization with ParaView. export_file(mesh, 'output/surface_mesh.vtk', 'vtk/trimesh')
Cutting planes”
A set of cutting planes can be assigned to a
SurfaceMesh
to cut away parts of the mesh for visualization purposes. This may be useful to e.g. cut a hole into a closed surface allowing to look inside the enclosed volume. TheSurfaceMesh
objects manages a list of cutting planes, which are accessible through theget_cutting_planes()
andset_cutting_planes()
methods. Note that the cuts are non-destructive and get performed only on the transient, non-periodic version of the mesh generated during image rendering or when exporting the mesh to a file. The original data structures of the mesh are not affected. TheSliceModifier
, which can act on aSurfaceMesh
, performs the slice by simply adding a new entry to theSurfaceMesh
’s list of cutting planes.Mesh data access
The methods
get_vertices()
,get_faces()
andget_face_adjacency()
methods provide access to the internal data of the surface mesh.-
property
domain
¶ The
SimulationCell
describing the (possibly periodic) domain which this surface mesh is embedded in. Note that this cell generally is independent of and may be different from thecell
found in theDataCollection
.
-
property
faces
¶ The
PropertyContainer
storing the properties of the faces of the mesh.
-
get_cutting_planes
()¶ Returns a N x 4 array containing the definitions of the N cutting planes attached to this
SurfaceMesh
.Each plane is defined by its unit normal vector and a signed displacement magnitude, which determines the plane’s distance from the coordinate origin along the normal, giving four numbers per plane in total. Those parts of the surface mesh which are on the positive side of the plane (in the direction the normal vector) are cut away.
Note that the returned Numpy array is a copy of the internal data stored by the
SurfaceMesh
.
-
get_face_adjacency
()¶ Returns a M x 3 array listing the indices of the three faces that are adjacent to each of the M triangle faces in the mesh. This information can be used to traverse the neighbors of triangle faces. Every triangle face has exactly three neighbors, because surface meshes are closed manifolds.
-
get_faces
()¶ Returns a M x 3 array with the vertex indices of the M triangles in the mesh. Note that the returned Numpy array is a copy of the internal data stored by the
SurfaceMesh
. Also keep in mind that triangle faces can cross the domain boundaries if the periodic boundary conditions are used.
-
get_vertices
()¶ Returns a N x 3 array with the xyz coordinates of the N vertices in the mesh. Note that the returned Numpy array is a copy of the internal data stored by the
SurfaceMesh
.
-
locate_point
(pos, eps=1e-6)¶ Determines the index of the spatial region that contains the given location in 3-D space. Note that region index -1 is typically reserved for the empty/exterior region, which doesn’t contain any atoms or particles. Regions starting at index 0 are used for filled or interior regions.
The parameter eps is used as a precision threshold to detect cases where the query point is positioned exactly on the surface itself, i.e. on the boundary between two spatial regions. Such a condition is indicated by the special return value
None
. You can set eps to 0.0 to disable the point-on-boundary test. Then the method will never yieldNone
as a result.- Parameters
pos – The (x,y,z) coordinates of the query point
eps – Numerical precision threshold for point-on-boundary test
- Returns
The ID of the spatial region containing pos; or
None
if pos is exactly on the dividing surface between two regions
-
property
regions
¶ The
PropertyContainer
storing the properties of the spatial regions of the mesh.
-
set_cutting_planes
(planes)¶ Sets the cutting planes to be applied to this
SurfaceMesh
. The array planes must follow the same format as the one returned byget_cutting_planes()
.
-
property
space_filling_region
¶ Indicates the index of the spatial region that fills the entire domain in case the surface is degenerate, i.e. the mesh has zero faces. The invalid index -1 is typically associated with the empty (exterior) region.
-
property
vertices
¶ The
PropertyContainer
storing the vertex properties of the mesh, including the vertex coordinates.
-
property
-
class
ovito.data.
TrajectoryLines
¶ Base:
ovito.data.PropertyContainer
Data object that stores the trajectory lines of a set of particles, which have been traced by the
GenerateTrajectoryLinesModifier
. It is typically part of a pipeline’s output data collection, from where it can be accessed via theDataCollection.trajectories
field.A
TrajectoryLines
object has an associatedTrajectoryVis
element, which controls the visual appearance of the trajectory lines in rendered images. This visual element is accessible through thevis
attribute of the base class.
-
class
ovito.data.
TriangleMesh
¶ Base:
ovito.data.DataObject
This data object type stores a triangle mesh describing general polyhedral geometry, which may be visualized side by side with the simulation data. Typically, triangle meshes are imported from an external geometry file using the
import_file()
function. See also the corresponding page of the user manual for more information on triangle meshes. This class has no parameters or methods.See also the
SurfaceMesh
class, which is another data object type in OVITO representing surface geometries. In contrast to triangle meshes, however, surface meshes can live inside of periodic simulation domains and are true orientable, closed manifolds. Furthermore, surface meshes can store arbitrary per-vertex and per-face property data – something triangle meshes cannot do.The visual appearance of the triangle mesh in rendered images is controlled through the
TriangleMeshVis
element attached to this data object (see base classDataObject
).
-
class
ovito.data.
VoxelGrid
¶ Base:
ovito.data.PropertyContainer
A two- or three-dimensional structured grid. Each cell (voxel) of the grid is of the same size and shape. The geometry of the entire grid, its
domain
, is defined by an attachedSimulationCell
object, which specific a three-dimensional parallelepiped or a two-dimensional parallelogram. See also the corresponding user manual page for more information on this object type.The
shape
property of the grid specifies the number of voxels along each domain cell vector. The size of an individual voxel is given by domain cell size divided by the number of voxels in each spatial direction.Every voxel of the grid may be associated with one or more field values. The data for these voxel properties is stored in standard
Property
objects, similar to particle or bond properties. Voxel properties can be accessed by name through the dictionary interface that theVoxelGrid
class inherits from itsPropertyContainer
base class.Voxel grids can be loaded from input data files, e.g. a CHGCAR file containing the electron density computed by the VASP code, or they can be dynamically generated within OVITO. The
SpatialBinningModifier
lets you project the information associated with the unstructured particle set to a structured voxel grid.Given a voxel grid, the
CreateIsosurfaceModifier
can then generate aSurfaceMesh
representing an isosurface for a field quantity defined on the voxel grid.Example
The following code example demonstrates how to create a new
VoxelGrid
from scratch and initialize it with data from a Numpy array:# Starting with an empty DataCollection: data = DataCollection() # Create a new SimulationCell object defining the outer spatial dimensions # of the grid and the boundary conditions, and add it to the DataCollection: data.cell = SimulationCell(pbc=(True, True, True), vis=SimulationCellVis(line_width=0.03)) data.cell_[:,:3] = [[10,0,0],[0,10,0],[0,0,10]] # Generate a three-dimensional Numpy array containing the grid cell values. nx = 10; ny = 6; nz = 8 field_data = numpy.random.random((nx, ny, nz)) # Create the VoxelGrid object and give it a unique identifier by which it can be referred to later on. # Link the voxel grid to the SimulationCell object created above, which defines its spatial extensions. # Specify the shape of the grid, i.e. the number of cells in each spatial direction. # Finally, assign a VoxelGridVis visual element to the data object to make the grid visible in the scene. grid = VoxelGrid( identifier = 'field', domain = data.cell, shape = field_data.shape, vis = VoxelGridVis(enabled=True, transparency=0.6)) # Associate a new property with the voxel grid cells and initialize it with the data from the Numpy array. # Note that the data must be provided as linear (1-dim.) array with the following type of memory layout: # The first grid dimension (x) is the fasted changing index while the third grid dimension (z) is the slowest varying index. # In this example, this corresponds to the "Fortran" memory layout of Numpy. grid.create_property('Field Value', data=field_data.flatten(order='F')) # Insert the VoxelGrid object into the DataCollection. data.objects.append(grid) # For demonstration purposes, compute an isosurface on the basis of the VoxelGrid created above. data.apply(CreateIsosurfaceModifier(operate_on='voxels:field', property='Field Value', isolevel=0.7))
-
property
domain
¶ The
SimulationCell
describing the (possibly periodic) domain which this grid is embedded in. Note that this cell generally is independent of and may be different from thecell
found in theDataCollection
.
-
property
shape
¶ A tuple with the numbers of grid cells along each of the three cell vectors of the
domain
.For two-dimensional grids, for which the
is2D
property of thedomain
is set to true, the third entry of theshape
tuple is always equal to 1.Assigning new shape dimensions to the grid automatically resizes the one-dimensional data arrays of the
PropertyContainer
base class if necessary and updates itsPropertyContainer.count
property to be equal to the product of the three shape dimensions, i.e. the total number of voxel cells.
-
property