Class pypfc

Bases: setup_io
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.
Functions¶
__init__ ¶
__init__(
domain_size: Union[List[float], ndarray],
ndiv: Optional[Union[List[int], ndarray]] = None,
config: Optional[Dict[str, Any]] = None,
**kwargs: Any
) -> None
Initialize PFC simulation with domain parameters and configuration.
Sets up the complete simulation environment including grid discretization, device configuration, and all numerical parameters needed for PFC evolution.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
domain_size
|
array_like of float, shape (3,)
|
Physical size of simulation domain [Lx, Ly, Lz] in lattice parameter units. This defines the spatial extent of the simulation box. |
required |
ndiv
|
array_like of int, shape (3,)
|
Number of grid divisions [nx, ny, nz]. All values must be even numbers
for FFT compatibility. If not provided, automatically calculated as
|
None
|
config
|
dict
|
Configuration parameters as key-value pairs. See the pyPFC overview for a complete list of the configuration parameters. |
None
|
**kwargs
|
dict
|
Individual configuration parameters passed as keyword arguments. These will override any corresponding values in the config dictionary. |
{}
|
Raises:
| Type | Description |
|---|---|
ValueError
|
If any element in |
Examples:
>>> # Custom grid and parameters using config dictionary
>>> config = {'dtime': 5e-5, 'struct': 'BCC', 'verbose': True}
>>> sim = setup_simulation([8.0, 8.0, 8.0], ndiv=[64, 64, 64], config=config)
cleanup ¶
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 ¶
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 |
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 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 the current two-point correlation function.
Returns:
| Type | Description |
|---|---|
Tensor
|
Two-point correlation function in Fourier space. |
get_alat ¶
Get the current lattice parameter.
Returns:
| Type | Description |
|---|---|
float
|
Current lattice parameter value. |
get_alpha ¶
Get the current alpha parameters.
Returns:
| Type | Description |
|---|---|
float
|
Current alpha parameter values. |
get_density ¶
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 the current time step size.
Returns:
| Type | Description |
|---|---|
float
|
Current time step size for integration. |
get_energy ¶
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 ¶
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 ¶
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 the current integration scheme parameters.
Returns:
| Type | Description |
|---|---|
array_like
|
Current parameters for the time integration scheme. |
set_C2_d ¶
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 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 the lattice parameter.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
alat
|
float
|
Lattice parameter. |
required |
set_alpha ¶
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 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 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
|
set_phase_field_smoothing_kernel ¶
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
|
set_update_scheme ¶
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 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 |