Skip to content

Complete API Reference

Listing of all methods across all pyPFC classes.

Class pypfc_grid

copy_from

copy_from(grid: setup_grid) -> None

Copy grid parameters from another grid object.

This method copies all grid configuration parameters from another setup_grid instance, including domain size, grid divisions, and derived parameters.

Parameters:

Name Type Description Default
grid setup_grid

Another setup_grid instance to copy parameters from.

required

get_ddiv

get_ddiv() -> np.ndarray

Get the grid spacing in each direction.

Returns:

Type Description
ndarray

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

get_domain_size

get_domain_size() -> np.ndarray

Get the physical domain size in each direction.

Returns:

Type Description
ndarray

Physical domain size [Lx, Ly, Lz] in lattice parameters.

get_ndiv

get_ndiv() -> np.ndarray

Get the number of grid divisions in each direction.

Returns:

Type Description
ndarray

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

set_ddiv

set_ddiv(ddiv: Union[List[float], ndarray]) -> None

Set the grid spacing in each direction.

Parameters:

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

Grid spacing in each direction [dx, dy, dz].

required

set_ndiv

set_ndiv(ndiv: Union[List[int], ndarray]) -> None

Set the number of grid divisions in each direction.

Updates the grid division parameters and related grid point counts. All divisions must be even numbers for FFT compatibility.

Parameters:

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

Number of grid divisions in each direction [nx, ny, nz]. Must be even numbers.

required

Raises:

Type Description
ValueError

If any value in ndiv is not an even number.

Class pypfc_base

evaluate_C2_d

evaluate_C2_d() -> torch.Tensor

Establish the two-point correlation function for a crystal structure.

Computes the two-point pair correlation function in Fourier space for the specified crystal structure using Gaussian peaks at reciprocal lattice positions.

Returns:

Name Type Description
C2_d (Tensor, shape(nx, ny, nz // 2 + 1))

Two-point pair correlation function on the computational device.

Raises:

Type Description
ValueError

If C20_alpha is negative when C20_amplitude is non-zero.

evaluate_directional_correlation_kernel

evaluate_directional_correlation_kernel(
    H0: ndarray, Rot: ndarray
) -> torch.Tensor

Establish directional correlation kernel for a crystal structure.

Computes directional correlation kernels used in extended PFC models to introduce orientational dependence.

Parameters:

Name Type Description Default
H0 float

Constant modulation of the peak height.

required
Rot ndarray of float, shape (3, 3) or None

Lattice rotation matrix. If None, uses identity matrix.

required

Returns:

Name Type Description
H_d (Tensor, shape(nx, ny, nz // 2 + 1))

Directional correlation kernel on the computational device.

evaluate_k2_d

evaluate_k2_d() -> torch.Tensor

Evaluate the sum of squared wave vectors for FFT operations.

Computes \(k^2 = k_x^2 + k_y^2 + k_z^2\) on the computational device using PyTorch FFT frequency grids. This is fundamental for Fourier-space operations in PFC simulations.

Returns:

Name Type Description
k2_d (Tensor, shape(nx, ny, nz_half))

Sum of squared wave vectors on the device. The z-dimension is reduced due to real FFT symmetry (nz_half = nz//2 + 1).

evaluate_reciprocal_planes

evaluate_reciprocal_planes() -> torch.Tensor

Establish reciprocal vectors/planes for a crystal structure.

Computes reciprocal lattice plane spacing (d-spacing) and wave vectors for crystallographic planes. For cubic systems: d = a / sqrt(h² + k² + l²) where a is the lattice parameter, and reciprocal spacing is k = 2π/d.

Returns:

Name Type Description
k_plane ndarray of float

Reciprocal lattice plane spacings (wave vector magnitudes).

n_plane ndarray of int

Number of symmetrical planes in each family.

den_plane ndarray of float

Atomic density within each plane family.

Raises:

Type Description
ValueError

If the crystal structure is not supported.

ValueError

If there are not enough peaks defined for the requested number of peaks.

Notes

For any family of lattice planes separated by distance d, there are reciprocal lattice points at intervals of 2π/d in reciprocal space.

get_alpha

get_alpha() -> np.ndarray

Get the Gaussian peak widths for the two-point correlation function.

Returns:

Name Type Description
alpha ndarray of float

Gaussian peak widths (α_i) for each peak in the correlation function.

get_csp

get_csp(
    pos: ndarray, normalize_csp: bool = False
) -> np.ndarray

Calculate the centro-symmetry parameter (CSP) for atoms.

Computes CSP values for atoms in a 3D periodic domain to identify crystal defects and disorder. CSP quantifies deviation from centro-symmetric local environments.

Parameters:

Name Type Description Default
pos ndarray of float, shape (natoms, 3)

3D coordinates of atoms.

required
normalize_csp bool

If True, normalizes CSP values to range [0,1].

False

Returns:

Name Type Description
csp ndarray of float, shape (natoms,)

Centro-symmetry parameter for each atom.

References

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

get_device_number

get_device_number() -> int

Get the current GPU device number.

Returns:

Type Description
int

Current GPU device index.

get_device_type

get_device_type() -> str

Get the current computation device type.

Returns:

Type Description
str

Current device type ('CPU' or 'GPU').

get_dtype_cpu

get_dtype_cpu() -> type

Get the current CPU data type.

Returns:

Type Description
dtype

Current NumPy data type used for CPU arrays.

get_dtype_gpu

get_dtype_gpu() -> torch.dtype

Get the current GPU data type.

Returns:

Type Description
dtype

Current PyTorch data type used for GPU tensors.

get_field_average_along_axis

get_field_average_along_axis(
    field: ndarray, axis: int
) -> np.ndarray

Evaluate the mean value of a field variable along a specified axis.

Computes the spatial average of a 3D field along one axis, reducing the dimensionality by averaging over the other two axes.

Parameters:

Name Type Description Default
field ndarray of float, shape (nx, ny, nz)

3D field variable to be averaged.

required
axis str

Axis to average along: x, y or z (case insensitive).

required

Returns:

Name Type Description
result ndarray of float

1D array containing mean values along the specified axis. Shape depends on the axis: (nx,), (ny,), or (nz,).

Raises:

Type Description
ValueError

If axis is not x, y or z.

get_integrated_field_along_axis

get_integrated_field_along_axis(
    field: ndarray, axis: int
) -> np.ndarray

Integrate a field variable along a specified axis.

Performs numerical integration of a 3D field variable along one axis, integrating over the two orthogonal directions.

Parameters:

Name Type Description Default
field ndarray of float, shape (nx, ny, nz)

3D field variable to be integrated.

required
axis str

Axis to integrate along: x, y or z (case insensitive).

required

Returns:

Name Type Description
result ndarray of float

1D array containing integrated values along the specified axis. Shape depends on the axis: (nx,), (ny,), or (nz,).

Raises:

Type Description
ValueError

If axis is not x, y or z.

get_integrated_field_in_volume

get_integrated_field_in_volume(
    field: ndarray, limits: Union[List[float], ndarray]
) -> float

Integrate a field variable within a defined volume.

Performs numerical integration of a field variable over a specified 3D volume on a Cartesian grid.

Parameters:

Name Type Description Default
field ndarray of float, shape (nx, ny, nz)

Field to be integrated over the specified volume.

required
limits array_like of float, length 6

Spatial integration limits: [xmin, xmax, ymin, ymax, zmin, zmax].

required

Returns:

Name Type Description
result float

Result of the volume integration.

get_k

get_k(npoints: int, dspacing: float) -> np.ndarray

Define a 1D wave vector for Fourier space operations.

Parameters:

Name Type Description Default
npoints int

Number of grid points. Must be even.

required
dspacing float

Grid spacing in real space.

required

Returns:

Name Type Description
k ndarray of float

1D wave vector array with proper frequency ordering for FFTs.

Raises:

Type Description
ValueError

If npoints is not an even number.

get_k2_d

get_k2_d() -> torch.Tensor

Get the wave vector magnitude squared tensor.

Returns:

Type Description
Tensor

Wave vector magnitude squared (k²) tensor in Fourier space.

get_phase_field_contour

get_phase_field_contour(
    pf: Union[ndarray, Tensor],
    pf_zoom: float = 1.0,
    evaluate_volume: bool = True,
) -> Union[Tuple[np.ndarray, float], np.ndarray]

Find the iso-contour surface of a 3D phase field using marching cubes.

Extracts iso-surfaces from 3D phase field data using the marching cubes algorithm, with optional volume calculation for enclosed regions.

Parameters:

Name Type Description Default
pf ndarray of float, shape (nx, ny, nz)

3D phase field data for iso-surface extraction.

required
pf_zoom float

Zoom factor for spatial coarsening/refinement.

1.0
evaluate_volume bool

If True, calculates the volume enclosed by the iso-surface.

True

Returns:

Name Type Description
verts ndarray of float, shape (n_vertices, 3)

Vertices of the iso-surface triangulation.

faces ndarray of int, shape (n_faces, 3)

Surface triangulation topology (vertex indices).

volume (float, optional)

Volume enclosed by the iso-surface (only if evaluate_volume=True).

get_rlv

get_rlv(struct: str, alat: float) -> np.ndarray

Get the reciprocal lattice vectors for a crystal structure.

Computes reciprocal lattice vectors for common crystal structures used in phase field crystal modeling.

Parameters:

Name Type Description Default
struct str

Crystal structure type. Options: 'SC', 'BCC', 'FCC', 'DC'.

required
alat float

Lattice parameter.

required

Returns:

Name Type Description
rlv ndarray of float, shape (nrlv, 3)

Reciprocal lattice vectors for the specified crystal structure.

Raises:

Type Description
ValueError

If the crystal structure is not supported.

get_time_stamp

get_time_stamp() -> str

Get current timestamp string.

Returns:

Name Type Description
timestamp str

Current date and time in format: YYYY-MM-DD HH:MM.

get_torch_threads

get_torch_threads() -> Tuple[int, int]

Get the current PyTorch thread configuration.

Returns:

Type Description
tuple of int

(num_threads, num_interop_threads) for PyTorch operations.

get_verbose

get_verbose() -> bool

Get the current verbose output setting.

Returns:

Type Description
bool

Current verbose mode setting.

get_xtal_nearest_neighbors

get_xtal_nearest_neighbors() -> (
    Tuple[np.ndarray, np.ndarray]
)

Get nearest neighbor information for crystal structures.

Computes nearest neighbor distances and coordination numbers for common crystal structures used in phase field crystal modeling.

Returns:

Name Type Description
nnb ndarray of int

Number of nearest and next-nearest neighbors.

nnb_dist ndarray of float

Distances to the nearest and next-nearest neighbors.

Raises:

Type Description
ValueError

If the crystal structure is not supported.

interpolate_atoms

interpolate_atoms(
    intrp_pos: ndarray,
    pos: ndarray,
    values: ndarray,
    num_nnb: int = 8,
    power: int = 2,
) -> np.ndarray

Interpolate values at given positions using inverse distance weighting.

Performs 3D interpolation in a periodic domain using inverse distance weighting: interpolated_value = Σ(wi × vi) / Σ(wi), where wi = 1 / (di^power).

Parameters:

Name Type Description Default
intrp_pos ndarray of float, shape (n_intrp, 3)

3D coordinates of positions where values should be interpolated.

required
pos ndarray of float, shape (n_particles, 3)

3D coordinates of particles with known values.

required
values ndarray of float, shape (n_particles,)

Values at the particle positions to be interpolated.

required
num_nnb int

Number of nearest neighbors to use for interpolation.

8
power float

Power for inverse distance weighting.

2

Returns:

Name Type Description
interp_val ndarray of float, shape (n_intrp,)

Interpolated values at the specified positions.

interpolate_density_maxima

interpolate_density_maxima(
    den: Union[ndarray, Tensor],
    ene: Optional[Union[ndarray, Tensor]] = None,
    pf: Optional[Union[ndarray, Tensor]] = None,
) -> Tuple[
    np.ndarray, Optional[np.ndarray], Optional[np.ndarray]
]

Find density field maxima and interpolate atomic positions and properties.

Identifies local maxima in the density field as atomic positions and performs high-order interpolation to obtain sub-grid precision coordinates. Also interpolates associated field values (density, energy, phase fields) at the atomic positions.

Parameters:

Name Type Description Default
den ndarray of float, shape (nx,ny,nz)

Density field from PFC simulation.

required
ene ndarray of float, shape (nx,ny,nz)

Energy field for interpolation at atomic positions.

None
pf list of ndarray

List of phase fields for interpolation at atomic positions. Each array should have shape (nx,ny,nz).

None

Returns:

Name Type Description
atom_coord ndarray of float, shape (n_maxima,3)

Interpolated coordinates of density maxima (atomic positions).

atom_data ndarray of float, shape (n_maxima, 2+n_phase_fields)

Interpolated field values at atomic positions. Columns: [density, energy, pf1, pf2, ..., pfN]

Notes

The method uses scipy.ndimage for high-order interpolation and applies density thresholding and merging of nearby maxima to remove spurious peaks. The interpolation order is controlled by the _density_interp_order attribute.

set_alpha

set_alpha(alpha: Union[List[float], ndarray]) -> None

Set the Gaussian peak widths for the two-point correlation function.

Parameters:

Name Type Description Default
alpha array_like of float

Gaussian peak widths (α_i) for each peak in the correlation function.

required

set_device_number

set_device_number(device_number: int) -> None

Set the GPU device number for multi-GPU systems.

Parameters:

Name Type Description Default
device_number int

GPU device index (0, 1, 2, ...) for CUDA computations.

required

set_device_type

set_device_type(device_type: str) -> None

Set the computation device type.

Parameters:

Name Type Description Default
device_type str

Device type for computations. Options: 'CPU', 'GPU'.

required

set_dtype_cpu

set_dtype_cpu(dtype: type) -> None

Set the CPU data type for numpy arrays.

Parameters:

Name Type Description Default
dtype dtype

NumPy data type for CPU computations (e.g., np.float32, np.float64).

required

set_dtype_gpu

set_dtype_gpu(dtype: dtype) -> None

Set the GPU data type for PyTorch tensors.

Parameters:

Name Type Description Default
dtype dtype

PyTorch data type for GPU computations (e.g., torch.float32, torch.float64).

required

set_k2_d

set_k2_d(k2_d: Tensor) -> None

Set the wave vector magnitude squared tensor.

Parameters:

Name Type Description Default
k2_d Tensor

Wave vector magnitude squared (k²) tensor in Fourier space. Used for FFT-based operations and differential operators.

required

set_torch_threads

set_torch_threads(
    nthreads: int, nthreads_interop: int
) -> None

Set PyTorch thread configuration for CPU operations.

Parameters:

Name Type Description Default
nthreads int

Number of threads for intra-op parallelism.

required
nthreads_interop int

Number of threads for inter-op parallelism.

required

set_verbose

set_verbose(verbose: bool) -> None

Set verbose output mode for debugging and monitoring.

Parameters:

Name Type Description Default
verbose bool

If True, enables detailed timing and progress output.

required

Class pypfc_pre

do_bicrystal

do_bicrystal(
    xtal_rot: ndarray,
    params: Optional[List[float]] = None,
    liq_width: float = 0.0,
    model: int = 0,
) -> np.ndarray

Define a bicrystal with two different crystal orientations.

Bicrystal Example

Parameters:

Name Type Description Default
xtal_rot ndarray of float, shape (3,3,2)

Crystal orientations (rotation matrices) for the two grains.

required
params list

List containing parameters for the bicrystal model.

  • model=0: [r] - cylindrical crystal radius
  • model=1: [r] - spherical crystal radius
  • model=2: [gb_x1, gb_x2] - grain boundary positions along x
None
liq_width float

Width of the liquid band along the grain boundary.

0.0
model int

Density field layout.

  • 0: Cylindrical crystal, extending through z
  • 1: Spherical crystal
  • 2: Bicrystal with two planar grain boundaries, normal to x
0

Returns:

Name Type Description
density ndarray of float, shape (nx,ny,nz)

Density field.

Raises:

Type Description
ValueError

If the value of model is not supported (should be 0, 1 or 2).

do_polycrystal

do_polycrystal(
    xtal_rot: ndarray,
    params: Optional[List[float]] = None,
    liq_width: float = 0.0,
    model: int = 0,
) -> np.ndarray

Define a polycrystal in a periodic 3D domain.

Polycrystal Example

Parameters:

Name Type Description Default
xtal_rot ndarray of float, shape (3,3,n_xtal)

Crystal orientations (rotation matrices) for n_xtal crystals.

required
params list

List containing parameters for the polycrystal model.

  • model=0: No parameter needed. The number of crystal seeds is determined from the number of provided orientations.
None
liq_width float

Width of the liquid band along the grain boundaries.

0.0
model int

Density field layout.

  • 0: A row of cylindrical seeds along y, with cylinders extending through z
0

Returns:

Name Type Description
density ndarray of float, shape (nx,ny,nz)

Polycrystal density field.

Raises:

Type Description
ValueError

If the value of model is not supported (should be 0).

do_single_crystal

do_single_crystal(
    xtal_rot: Optional[ndarray] = None,
    params: Optional[List[float]] = None,
    model: int = 0,
) -> np.ndarray

Define a single crystal in a periodic 3D domain.

Single Crystal Example

Parameters:

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

Crystal orientation (rotation matrix). Default is an identity matrix.

None
params list

List containing parameters for the single crystal model:

  • model=0: [r] - spherical crystal radius
  • model=1: [x1, x2] - crystal extent in x direction
  • model=2: [r] - cylindrical crystal radius
None
model (int, optional)

Density field layout.

  • 0: Spherical crystal
  • 1: Crystal extending throughout y and z, covering interval in x
  • 2: Cylindrical crystal, extending through z
0

Returns:

Name Type Description
density ndarray of float, shape (nx,ny,nz)

Density field.

Raises:

Type Description
ValueError

If the value of model is not supported (should be 0 or 1).

evaluate_ampl_dens

evaluate_ampl_dens() -> Tuple[np.ndarray, np.ndarray]

Get density field amplitudes and phase densities for PFC simulations.

Amplitude and Density Example

Returns:

Name Type Description
ampl ndarray of float, shape (npeaks,)

Density field amplitudes for the specified crystal structure and number of peaks.

nLnS ndarray of float, shape (2,)

Densities in the liquid (nL) and solid (nS) phases.

Raises:

Type Description
ValueError

If npeaks is not supported for the current crystal structure.

ValueError

If sigma value is not supported for the current configuration.

ValueError

If struct is not 'SC', 'BCC' or 'FCC' or if amplitudes and densities are not available for the specified structure.

Notes

This method provides pre-calculated density field amplitudes and phase densities for different crystal structures (SC, BCC, FCC) and numbers of Fourier peaks in the two-point correlation function. The values depend on the temperature-like parameter sigma.

Fore efficiency, the method uses lookup tables of pre-computed values.

Supported configurations:

struct npeaks sigma range
SC 1 [0.0, 0.25]
2 [0.0, 0.22]
3 [0.0, 0.22]
BCC 1 [0.0, 0.33]
2 [0.0, 0.30]
3 [0.0, 0.25]
FCC 1 [0.0, 0.25]
2 [0.0, 0.25]
3 [0.0, 0.25]

generate_density_field

generate_density_field(
    crd: ndarray, g: ndarray
) -> np.ndarray

Define a 3D density field for (X)PFC modeling.

Parameters:

Name Type Description Default
crd ndarray of float, shape (3,...)

Grid point coordinates [x,y,z].

required
g ndarray of float, shape (3,3)

Rotation matrix for crystal orientation.

required

Returns:

Name Type Description
density ndarray of float

Density field for the specified crystal structure with appropriate Fourier modes and amplitudes.

Raises:

Type Description
ValueError

If struct is not one of the supported crystal structures ('SC', 'BCC', 'FCC', 'DC').

Notes

The density field is generated based on the current crystal structure (struct) and density field amplitudes (ampl) settings.

get_ampl

get_ampl() -> np.ndarray

Get the amplitudes in the density approximation.

Returns:

Name Type Description
ampl ndarray of float, shape (N,)

Amplitudes.

get_density

get_density() -> np.ndarray

Get the density field.

Returns:

Name Type Description
den ndarray of float, shape (nx,ny,nz)

Density field.

get_nlns

get_nlns() -> np.ndarray

Get the liquid and solid phase densities.

Returns:

Name Type Description
nlns ndarray of float, shape (2,)

\([n_{l},n_{s}]\) where \(n_{l}\) is liquid phase density and \(n_{s}\) is solid phase density.

get_npeaks

get_npeaks() -> int

Get the number of peaks in the density field approximation.

Returns:

Name Type Description
npeaks int

Number of peaks in the density field approximation.

get_sigma

get_sigma() -> float

Get the temperature-like parameter sigma.

Returns:

Name Type Description
sigma float

Temperature-like parameter sigma

get_struct

get_struct() -> str

Get the crystal structure.

Returns:

Name Type Description
struct str

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

set_ampl

set_ampl(ampl: Union[List[float], ndarray]) -> None

Set the amplitudes in the density approximation.

Parameters:

Name Type Description Default
ampl array_like of float, shape (N,)

Amplitudes.

required

set_density

set_density(den: ndarray) -> None

Set the density field.

Parameters:

Name Type Description Default
den ndarray of float, shape (nx,ny,nz)

Density field.

required

set_energy

set_energy(ene: ndarray) -> None

Set the PFC energy field.

Parameters:

Name Type Description Default
ene ndarray of float, shape (nx,ny,nz)

PFC energy field.

required

set_nlns

set_nlns(nlns: Union[List[float], ndarray]) -> None

Set the liquid and solid phase densities.

Parameters:

Name Type Description Default
nlns array_like of float, shape (2,)

\([n_{l},n_{s}]\) where \(n_{l}\) is liquid phase density and \(n_{s}\) is solid phase density.

required

set_npeaks

set_npeaks(npeaks: int) -> None

Set the number of peaks in the density field approximation.

Parameters:

Name Type Description Default
npeaks int

Number of peaks in the density field approximation.

required

set_sigma

set_sigma(sigma: float) -> None

Set the temperature-like parameter sigma.

Parameters:

Name Type Description Default
sigma float

Temperature-like parameter sigma.

required

set_struct

set_struct(struct: str) -> None

Set the crystal structure.

Parameters:

Name Type Description Default
struct (FCC, BCC)

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

'FCC'

Class pypfc_io

I/O operations for Phase Field Crystal simulation data.

This class provides comprehensive file I/O functionality for PFC simulations including:

  • Extended XYZ format for atomic positions and properties
  • VTK output for visualization in ParaView/VisIt
  • Binary pickle serialization for Python objects
  • Structured grid data export for continuous fields
  • Point data export for atomic/particle systems
  • Simulation metadata and setup information files

The class supports both text and binary formats, with optional compression for large datasets. All methods are designed to handle typical PFC simulation outputs including atomic positions, density fields, energy fields and phase field data.

Notes

The extended XYZ format follows the convention used in molecular dynamics and materials simulation communities, allowing storage of arbitrary per-atom properties alongside coordinates. VTK output enables direct visualization in scientific visualization software.

append_to_info_file

append_to_info_file(
    info: Union[str, List[str]],
    filename: str = "pypfc_simulation.txt",
    output_path: Optional[str] = None,
) -> None

Append information to a text file.

Adds new content to an existing text file, useful for logging simulation progress or adding additional information to setup files.

Parameters:

Name Type Description Default
info str or list of str

String or list of strings to append to the file.

required
filename str

Name of the output file.

'pypfc_simulation.txt'
output_path str

Path to the output directory. Uses the current working directory as default.

None

load_pickle

load_pickle(filename: str, ndata: int) -> List[Any]

Load data objects from a binary pickle file.

Deserializes Python objects from a pickle file created with save_pickle. Reads a specified number of objects from the binary file.

Parameters:

Name Type Description Default
filename str

Path to input pickle file (without .pickle extension).

required
ndata int

Number of data objects to read from the file.

required

Returns:

Name Type Description
data list

List of Python objects loaded from filename.pickle.

Warning

Only load pickle files from trusted sources, as they can execute arbitrary code during deserialization.

read_extended_xyz

read_extended_xyz(
    filename: str, nfields: int = 0
) -> Tuple[
    np.ndarray, List[float], float, List[np.ndarray]
]

Read PFC data from extended XYZ format file.

Reads atomic positions and associated properties from an extended XYZ file, which may be compressed with gzip. Automatically detects file format and handles both .xyz and .xyz.gz extensions.

Parameters:

Name Type Description Default
filename str

Name of input XYZ file (with or without .xyz/.xyz.gz extension).

required
nfields int

Number of data fields per atom beyond coordinates [x,y,z].

0

Returns:

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

Atomic coordinates.

domain_size ndarray of float, shape (3,)

Domain size [Lx, Ly, Lz] from file header.

time float

Simulation time from file header.

partDen ndarray of float, shape (natoms,)

Particle density values (if available).

partEne ndarray of float, shape (natoms,)

Particle energy values (if available).

partPf ndarray of float, shape (natoms,)

Particle phase field values (if available).

save_pickle

save_pickle(filename: str, data: List[Any]) -> None

Save data objects to a binary pickle file.

Serializes a list of Python objects to a binary pickle file for efficient storage and later retrieval. This is useful for saving simulation state, configuration parameters or processed results.

Parameters:

Name Type Description Default
filename str

Path to the output file (without .pickle extension).

required
data list

List of Python objects to serialize and save.

required
Notes

The output file will be named filename.pickle. Pickle files are Python-specific binary format that preserves object structure and types.

Warning

Only load pickle files from trusted sources, as they can execute arbitrary code during deserialization.

write_extended_xyz

write_extended_xyz(
    filename: str,
    coord: ndarray,
    atom_data: List[ndarray],
    atom_data_labels: List[str],
    simulation_time: float = 0.0,
    gz: bool = True,
) -> None

Save PFC atomic data in extended XYZ format.

Writes atomic positions and associated properties to an extended XYZ file, which is a standard format in molecular dynamics and materials simulation. The format includes atomic coordinates plus arbitrary per-atom properties such as density, energy, phase field values, etc.

Parameters:

Name Type Description Default
filename str

Base name of the output XYZ file (without extension).

required
coord ndarray of float, shape (natoms, 3)

Atomic coordinates [x, y, z] for each atom.

required
atom_data list of ndarray

List of arrays containing per-atom data. Each array must have shape (natoms,) and represent a property for each atom.

required
atom_data_labels list of str

Labels for each data array in atom_data. Must have same length as atom_data list.

required
simulation_time float

Simulation time to include in file header.

0.0
gz bool

If True, compress output using gzip.

True
Notes

The output file will be named filename.xyz or filename.xyz.gz if compression is enabled. The extended XYZ format includes a header with the number of atoms, simulation time, and property labels, followed by atomic coordinates and properties.

Examples:

>>> write_extended_xyz('output', coord, [density, energy], 
...                   ['density', 'energy'], simulation_time=100.0)

write_info_file

write_info_file(
    filename: str = "pypfc_simulation.txt",
    output_path: Optional[str] = None,
) -> None

Write simulation setup information to a file.

Creates a text file containing simulation parameters, grid configuration, and other setup information for documentation and reproducibility.

Parameters:

Name Type Description Default
filename str

Name of the output file.

'pypfc_simulation.txt'
output_path str

Path to the output directory. Uses the current working directory as default.

None

write_vtk_points

write_vtk_points(
    filename: str,
    coord: ndarray,
    scalar_data: List[ndarray],
    scalar_data_name: List[str],
    vector_data: Optional[List[ndarray]] = None,
    vector_data_name: Optional[List[str]] = None,
    tensor_data: Optional[List[ndarray]] = None,
    tensor_data_name: Optional[List[str]] = None,
) -> None

Save 3D point data to VTK file for visualization.

Exports atomic positions and associated scalar, vector and tensor data to a VTK (Visualization Toolkit) file in binary XML format. This format is compatible with ParaView, VisIt and other scientific visualization software.

Parameters:

Name Type Description Default
filename str

Output file name (with or without .vtu extension).

required
coord ndarray of float, shape (natoms, 3)

3D coordinates of points/atoms.

required
scalar_data list of ndarray

List of scalar data arrays. Each array should have shape (natoms,).

required
scalar_data_name list of str

Names/labels for each scalar data array. Must match length of scalarData.

required
vector_data list of ndarray

List of vector data arrays. Each array should have shape (natoms, 3).

None
vector_data_name list of str

Names/labels for each vector data array. Required if vectorData is provided.

None
tensor_data list of ndarray

List of tensor data arrays. Each array should have shape (natoms, 3, 3). Tensors are automatically reshaped to VTK format (natoms, 9).

None
tensor_data_name list of str

Names/labels for each tensor data array. Required if tensorData is provided.

None
Notes

Tensor data is reshaped from (3,3) matrices to 9-component vectors following VTK convention:

\[T = \begin{bmatrix} T_{11} & T_{12} & T_{13} \\ T_{21} & T_{22} & T_{23} \\ T_{31} & T_{32} & T_{33} \end{bmatrix}\]

becomes: \([T_{11}, T_{12}, T_{13}, T_{21}, T_{22}, T_{23}, T_{31}, T_{32}, T_{33}]\)

write_vtk_structured_grid

write_vtk_structured_grid(
    filename: str,
    array_data: List[ndarray],
    array_name: List[str],
) -> None

Save 3D field data to VTK structured grid file.

Exports 3D field data (such as density, energy, or phase fields) to a VTK structured grid file in binary XML format. This format is ideal for visualizing continuous field data in ParaView, VisIt and other scientific visualization software.

Parameters:

Name Type Description Default
filename str

Output file name (with or without .vts extension).

required
array_data list of ndarray

List of 3D numpy arrays containing field data. Each array should have shape (nx, ny, nz) matching the simulation grid.

required
array_name list of str

Names/labels for each data array. Must match length of arrayData.

required
Notes

The output file will be named filename.vts and uses VTK's structured grid format with binary encoding for efficient storage. The grid dimensions and spacing are automatically determined from the inherited grid setup.

Class pypfc

This is the primary class for conducting PFC simulations, providing complete functionality for time evolution, energy evaluation, and phase field analysis. It combines all inherited capabilities from the class hierarchy:

  • Grid setup and domain discretization (pypfc_grid)
  • Mathematical operations and device management (pypfc_base)
  • Density field generation and crystal structures (pypfc_pre)
  • File I/O operations and data export (pypfc_io)

The class uses a configuration-driven approach where all simulation parameters are managed through a DEFAULTS dictionary. User configurations merge with defaults, to provide parameter handling with fallback values.

Notes

All grid divisions (ndiv) must be even numbers for FFT compatibility. The class automatically validates configuration parameters and provides informative error messages for invalid inputs.

cleanup

cleanup() -> None

Clean up allocated tensors and free device memory.

Explicitly deletes PyTorch tensors to free GPU/CPU memory. This is particularly important for GPU simulations to prevent memory leaks and ensure proper resource management.

Notes

This method should be called at the end of simulations, especially when running multiple simulations sequentially or when GPU memory is limited. The method automatically detects which tensors exist based on the update scheme and cleans up accordingly.

do_step_update

do_step_update() -> None

Update the PFC density field using the selected time integration scheme.

Performs one time step of the PFC evolution equation using the time integration method specified by _update_scheme. The method operates in Fourier space for computational efficiency and automatically handles the FFT/iFFT transformations.

Returns:

Name Type Description
f_den_d Tensor

Updated density field in Fourier space. The real-space density field is automatically updated in self._den_d via inverse FFT.

Raises:

Type Description
ValueError

If the specified update_scheme is not supported.

Notes

This method uses precomputed linear terms for efficiency. The appropriate update method is selected based on self._update_scheme and called with the corresponding precomputed Fourier-space operators.

evaluate_energy

evaluate_energy() -> float

Evaluate the PFC energy.

Computes the total free energy of the system using the phase field crystal energy functional.

Returns:

Name Type Description
ene Tensor

Energy field with shape [nx, ny, nz].

eneAv float

Average free energy density.

Notes

Energy is computed in Fourier space for efficiency and transformed back to real space for local energy density visualization.

get_C2_d

get_C2_d() -> torch.Tensor

Get the current two-point correlation function.

Returns:

Type Description
Tensor

Two-point correlation function in Fourier space.

get_alat

get_alat() -> float

Get the current lattice parameter.

Returns:

Type Description
float

Current lattice parameter value.

get_alpha

get_alpha() -> np.ndarray

Get the current alpha parameters.

Returns:

Type Description
float

Current alpha parameter values.

get_density

get_density() -> np.ndarray

Get the PFC density field and its mean value.

Returns the current density field and its spatial average, transferring data from GPU to CPU if necessary.

Returns:

Name Type Description
den ndarray of float, shape (nx,ny,nz)

Density field on CPU.

mean_den float

Spatially averaged density.

get_dtime

get_dtime() -> float

Get the current time step size.

Returns:

Type Description
float

Current time step size for integration.

get_energy

get_energy() -> np.ndarray

Get the PFC energy field and its mean value.

Computes the local energy density field and its spatial average using the current density field configuration.

Returns:

Name Type Description
ene ndarray of float, shape (nx,ny,nz)

Local energy density field on CPU.

mean_ene float

Spatially averaged energy density.

get_phase_field

get_phase_field() -> Union[np.ndarray, List[np.ndarray]]

Evaluate phase field using wavelet filtering.

Computes the phase field by applying wavelet convolution followed by Gaussian smoothing. The phase field identifies crystalline regions and their orientations within the PFC density field.

Returns:

Name Type Description
pf ndarray of float, shape (nx,ny,nz) or list of such arrays

Phase field(s) on CPU. Returns a single array for isotropic kernels or a list of arrays for directional analysis with multiple kernels.

Notes

The method automatically handles both single and multiple wavelet kernels for comprehensive grain orientation analysis.

get_update_scheme

get_update_scheme() -> str

Establish the PFC time integration scheme and return method handle.

Configures the selected time integration scheme and returns a function handle to the appropriate update method. This method sets up scheme-specific parameters and precomputed terms for efficient time stepping.

Returns:

Name Type Description
update_density callable

Function handle to the selected time integration method.

Raises:

Type Description
ValueError

If alpha, beta, gamma parameters are not provided for '2nd_order' scheme.

ValueError

If f_denOld_d is not provided for '2nd_order' scheme.

ValueError

If the specified update_scheme is not supported.

Notes

The method automatically precomputes linear terms in Fourier space to optimize performance during time stepping. Different schemes require different numbers of precomputed terms and have varying stability properties and computational costs.

get_update_scheme_params

get_update_scheme_params() -> np.ndarray

Get the current integration scheme parameters.

Returns:

Type Description
array_like

Current parameters for the time integration scheme.

set_C2_d

set_C2_d(C2_d: Tensor) -> None

Set the two-point correlation function in Fourier space.

Parameters:

Name Type Description Default
C2_d Tensor

Two-point correlation function tensor in Fourier space. Must match the grid dimensions and be on the correct device.

required

set_H2

set_H2(H0: float, Rot: ndarray) -> None

Set the directional correlation kernel for extended PFC models.

Configures the H2 kernel used in extended phase field crystal models to introduce directional correlations.

Parameters:

Name Type Description Default
H0 float

Amplitude of the directional correlation kernel.

required
Rot array_like or None

Rotation matrix for orienting the kernel. If None, uses identity orientation.

required

set_alat

set_alat(alat: float) -> None

Set the lattice parameter.

Parameters:

Name Type Description Default
alat float

Lattice parameter.

required

set_alpha

set_alpha(alpha: Union[List[float], ndarray]) -> None

Set the alpha parameters, controlling the widths of the Gaussian peaks of the pair correlation function.

Parameters:

Name Type Description Default
alpha float

Pair correlation peak widths.

required

set_density

set_density(density: ndarray) -> None

Set the PFC density field.

Updates the density field and automatically computes its Fourier transform for subsequent calculations.

Parameters:

Name Type Description Default
density ndarray of float, shape (nx,ny,nz)

New density field configuration.

required

set_dtime

set_dtime(dtime: float) -> None

Set the time step.

Parameters:

Name Type Description Default
dtime float

Time step size. Must be positive and small enough for stability.

required

set_phase_field_kernel

set_phase_field_kernel(
    H0: float = 1.0,
    Rot: Optional[Union[ndarray, List[ndarray]]] = None,
) -> None

Set phase field kernel for directional analysis.

Configures the correlation kernel used for phase field evaluation, allowing for directional filtering or isotropic analysis.

Parameters:

Name Type Description Default
H0 float

Kernel amplitude.

1.0
Rot ndarray of float, shape (3,3)

Rotation matrix for directional kernel. If None, uses isotropic two-point correlation function.

None

set_phase_field_smoothing_kernel

set_phase_field_smoothing_kernel(
    pf_gauss_var: Optional[float] = None,
) -> None

Set phase field smoothing kernel.

Configures Gaussian smoothing parameters for phase field calculations.

Parameters:

Name Type Description Default
pf_gauss_var float

Gaussian variance for smoothing kernel. If None, uses current value.

None

set_update_scheme

set_update_scheme(update_scheme: str) -> None

Set the time integration scheme for density evolution.

Parameters:

Name Type Description Default
update_scheme str

Time integration method. Options: '1st_order', '2nd_order', 'exponential'.

required

set_update_scheme_params

set_update_scheme_params(
    params: Union[List[float], ndarray],
) -> None

Set parameters for the time integration scheme.

Parameters:

Name Type Description Default
params array_like

Parameters specific to the chosen integration scheme. Format depends on the selected scheme type.

required

Class pypfc_ovito

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