Skip to content

Class pypfc_ovito

This class provides custom interfaces to selected functionalities in OVITO, useful for post-processing of pyPFC simulation output and crystal structure analysis.

The ovito Python package must be installed.

Bases: setup_grid

Functions

__init__

__init__(
    ndiv: Union[List[int], ndarray],
    ddiv: Union[List[float], ndarray],
    dtype_cpu: dtype = np.double,
    struct: str = "FCC",
    pos: Optional[ndarray] = None,
    verbose: bool = False,
) -> None

Initialize the pypfc_ovito setup with domain size and grid divisions.

Parameters:

Name Type Description Default
ndiv ndarray of int, shape (3,)

Number of grid divisions along each coordinate axis [nx, ny, nz].

required
ddiv ndarray of float, shape (3,)

Grid spacing along each coordinate axis [dx, dy, dz].

required
dtype_cpu dtype

CPU data type for numerical arrays.

double
struct str

Crystal structure type.

  • 'FCC': Face-centered cubic
  • 'BCC': Body-centered cubic
'FCC'
pos ndarray of float, shape (N,3)

Initial atomic positions. If None, positions will be generated based on the crystal structure.

None
verbose bool

Enable verbose output for debugging and monitoring.

False

do_ovito_csp

do_ovito_csp(num_neighbors: int = 12) -> np.ndarray

Evaluate the Centro-Symmetry Parameter (CSP) for each atom using OVITO.

The CSP value is a useful measure for identifying defects in crystal structures. Lower values indicate more perfect crystal structure.

Parameters:

Name Type Description Default
num_neighbors int

Number of nearest neighbors to consider for CSP calculation.

  • FCC structure: 12
  • BCC structure: 8
12

Returns:

Name Type Description
csp_values (ndarray, shape(N))

Array of CSP values for each atom.

References

A. Stukowski (2010), Visualization and analysis of atomistic simulation data with OVITO - the Open Visualization Tool, Modelling Simul. Mater. Sci. Eng. 18, 015012. https://doi.org/10.1088/0965-0393/18/1/015012

C.L. Kelchner, S.J. Plimpton and J.C. Hamilton, Dislocation nucleation and defect structure during surface indentation, Phys. Rev. B, 58(17):11085-11088, 1998. https://doi.org/10.1103/PhysRevB.58.11085

do_ovito_dxa

do_ovito_dxa(
    rep: Union[List[int], ndarray] = [1, 1, 1],
    tol: float = 1e-08,
) -> Tuple[
    np.ndarray,
    np.ndarray,
    np.ndarray,
    np.ndarray,
    np.ndarray,
    List[np.ndarray],
]

Perform dislocation analysis using OVITO's Dislocation eXtraction Algorithm (DXA).

Parameters:

Name Type Description Default
rep array_like of int, shape (3,)

Number of simulation box replications to consider in the periodic directions: [n_rep_x, n_rep_y, n_rep_z]. Default is no replication.

[1, 1, 1]
tol float

Tolerance parameter for dislocation type identification.

1e-08

Returns:

Name Type Description
disl_type_ids ndarray of int, shape (N,)

Array of dislocation type IDs:

  • 0: \(\langle 100\rangle\) Composite
  • 1: \(\frac{1}{2}\langle 110\rangle\) Perfect edge
  • 2: \(\frac{1}{6}\langle 112\rangle\) Shockley screw
  • 3: \(\frac{1}{6}\langle 110\rangle\) Stair-rod
  • 4: \(\frac{1}{3}\langle 100\rangle\) Hirth
  • 5: \(\frac{1}{3}\langle 111\rangle\) Frank
  • 6: other
disl_coord ndarray of float, shape (N,3)

Coordinates of the first point of each dislocation line.

disl_line_len ndarray of float, shape (N,)

Array of dislocation line lengths.

disl_line_dir ndarray of float, shape (N,3)

Array of dislocation line direction unit vectors.

disl_burg_vec ndarray of float, shape (N,3)

Array of Burgers vectors.

disl_segm_pts ndarray of int, shape (N,)

List of dislocation segment points.

References

A. Stukowski, V.V. Bulatov and A. Arsenlis (2012), Automated identification and indexing of dislocations in crystal interfaces, Modelling Simul. Mater. Sci. Eng. 20, 085007. https://doi.org/10.1088/0965-0393/20/8/085007

do_ovito_ptm

do_ovito_ptm(
    ref_rot: Optional[ndarray] = None,
    output_euler_ang: bool = False,
    output_strain: bool = False,
) -> Union[
    Tuple[np.ndarray, np.ndarray],
    Tuple[np.ndarray, np.ndarray, np.ndarray],
]

Evaluate crystal structure, crystallographic orientation and elastic strain tensor using OVITO's Polyhedral Template Matching (PTM) algorithm.

Parameters:

Name Type Description Default
ref_rot array_like of float, shape (3,3)

Reference rotation matrix.

None
output_euler_ang bool

Format of output orientation representation:

  • True: Euler angles (ZXZ convention).
  • False: Quaternions. This is the default.
False
output_strain bool

Whether to evaluate the local elastic Green-Lagrange strain tensors:

  • True: Compute strain tensors for each atom.
  • False: Skip strain calculation.
False

Returns:

Name Type Description
structure_id ndarray of int, shape (N,)

Integer array of structure types per atom:

  • 0: Other/Unknown
  • 1: FCC (face-centered cubic)
  • 2: HCP (hexagonal close-packed)
  • 3: BCC (body-centered cubic)
  • 4: ICO (icosahedral coordination)
  • 5: SC (simple cubic)
  • 6: Cubic diamond
  • 7: Hexagonal diamond
  • 8: Graphene
rot ndarray of float, shape (N,3) or (N,4)

Local crystal orientation expressed as Euler angles or quaternions, saved per atom as \([\varphi_{1}, \Psi, \varphi_{2}]\) or (if output_euler_ang=True).

strain ndarray of float, shape (N,6)

Local Green-Lagrange strain tensor (if output_strain=True), saved per atom as \([\epsilon_{xx}, \epsilon_{yy}, \epsilon_{zz}, \epsilon_{xy}, \epsilon_{xz}, \epsilon_{yz}]\).

References

P.M. Larsen et al. (2016), Robust Structural Identification via Polyhedral Template Matching, Modelling Simul. Mater. Sci. Eng. 24, 055007. https://doi.org/10.1088/0965-0393/24/5/055007

A. Stukowski (2010), Visualization and analysis of atomistic simulation data with OVITO - the Open Visualization Tool, Modelling Simul. Mater. Sci. Eng. 18, 015012. https://doi.org/10.1088/0965-0393/18/1/015012

get_coord

get_coord() -> np.ndarray

Get atom coordinates.

Returns:

Name Type Description
coord ndarray of float, shape (N,3)

Atom (x, y, z) coordinates.

get_struct

get_struct() -> str

Get crystal structure type.

Returns:

Name Type Description
struct str

Crystal structure type: 'FCC', 'BCC'.

set_coord

set_coord(coord: ndarray) -> None

Set atom coordinates.

Parameters:

Name Type Description Default
coord ndarray of float, shape (N,3)

Atom (x, y, z) coordinates.

required

set_struct

set_struct(struct: str) -> None

Set crystal structure type.

Parameters:

Name Type Description Default
struct (FCC, BCC)

Crystal structure type: 'FCC', 'BCC'.

'FCC'

set_verbose

set_verbose(verbose: bool) -> None

Set verbose output mode.

Parameters:

Name Type Description Default
verbose bool

Provide verbose output, or not.

required