Skip to content

Class pypfc_base

pyPFC logo

Bases: setup_grid

Functions

__init__

__init__(
    domain_size: ndarray,
    ndiv: ndarray,
    config: Dict[str, Any],
) -> None

Initialize the base PFC setup with domain parameters and device configuration.

Parameters:

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

Physical size of the simulation domain [Lx, Ly, Lz] in lattice parameter units.

required
ndiv ndarray of int, shape (3,)

Number of grid divisions [nx, ny, nz]. Must be even numbers for FFT compatibility.

required
config dict

Configuration parameters as key-value pairs. See the pyPFC overview for a complete list of the configuration parameters.

required

Raises:

Type Description
ValueError

If dtype_gpu is not torch.float32 or torch.float64.

ValueError

If GPU is requested but no GPU is available.

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