ferrodispcalc package

Submodules

ferrodispcalc.compute

ferrodispcalc.compute.calculate_averaged_structure(traj: list[Atoms], select: list[int] | slice | None = None) Atoms[source]

Compute the time-averaged atomic structure from an MD trajectory.

Atomic coordinates are unwrapped with respect to the first selected frame before averaging to avoid artefacts from periodic boundary crossings.

Parameters:
  • traj (list[Atoms]) – Full MD trajectory as a list of ASE Atoms objects.

  • select (list[int] | slice | None, optional) – Frame selection. None selects the last 50 % of frames. Defaults to None.

Returns:

ASE Atoms object with averaged positions and cell. Element symbols and PBC flags are taken from the first frame of the trajectory.

Return type:

Atoms

ferrodispcalc.compute.calculate_dielectric_constant(polarization: ndarray, volume: float, temperature: float = 300.0, atomic: bool = False) dict[str, float | ndarray][source]

Calculate dielectric tensor components from polarization fluctuations.

The dielectric tensor is computed from the fluctuation formula eps_ij = V / (eps0 * kB * T) * (<Pi Pj> - <Pi><Pj>), where the input polarization is assumed to be in C/m². The input volume must be given in ų and is converted internally to m³. When atomic=False, the polarization is first averaged over all unit cells. When atomic=True, the dielectric tensor is computed independently for each unit cell without spatial averaging, and volume is interpreted as the volume per unit cell.

Parameters:
  • polarization (np.ndarray) – Polarization time series. Shape (n_frames, n_cells, 3) for local polarization from calculate_polarization(), or (n_frames, 3) for an already averaged polarization trajectory.

  • volume (float) – Volume in ų. For atomic=False, this is the total volume of the system corresponding to the polarization trajectory. For atomic=True, this must be the volume of a single unit cell.

  • temperature (float, optional) – Temperature in K. Defaults to 300.0.

  • atomic (bool, optional) – Whether to compute local dielectric tensor components for each unit cell. Defaults to False.

Returns:

Dielectric tensor components eps_xx, eps_yy, eps_zz, eps_xy, eps_xz, and eps_yz. Each value is a scalar when atomic=False and an array of shape (n_cells,) when atomic=True.

Return type:

dict[str, float | np.ndarray]

Notes

This function computes only the ionic fluctuation contribution to the dielectric response and neglects the electronic contribution.

No statistical error analysis is performed internally. Users should apply their own block averaging or related analysis to estimate uncertainties. By default, all provided frames are used in the calculation.

Raises:

ValueError – If the input shape is invalid, if atomic=True but the input is not local polarization, or if volume or temperature is not positive.

ferrodispcalc.compute.calculate_displacement(traj: list[Atoms] | Atoms, nl: ndarray, select: list[int] | slice | None = None) ndarray[source]

Calculate ionic displacements from an MD trajectory.

For each center atom in the neighbor list, the displacement is computed as the mean vector from its neighbors to itself. Periodic boundary conditions are handled via the minimum image convention (MIC).

Parameters:
  • traj (list[Atoms] | Atoms) – Full MD trajectory (list of ASE Atoms) or a single Atoms object.

  • nl (np.ndarray) – Neighbor list array of shape (n_centers, n_neighbors + 1) with 1-based indices as returned by build_neighbor_list().

  • select (list[int] | slice | None, optional) – Frame selection. None selects the last 50 % of frames. Defaults to None.

Returns:

Displacement vectors in Ångström. Shape (n_frames, n_centers, 3) for multiple frames, or (n_centers, 3) for a single frame.

Return type:

np.ndarray

ferrodispcalc.compute.calculate_octahedral_tilt(traj: list[Atoms] | Atoms, nl_bo: ndarray, select: list[int] | slice | None = None) ndarray[source]

Calculate octahedral tilt angles in a perovskite structure.

For each octahedral center, the six neighboring anions are sorted into ±x, ±y, ±z pairs. The tilt angle about each Cartesian axis is the angle between the anion–anion bond vector and the corresponding unit vector.

Parameters:
  • traj (list[Atoms] | Atoms) – Full MD trajectory or a single ASE Atoms object.

  • nl_bo (np.ndarray) – Neighbor list for B–O pairs, shape (n_centers, 7) (center + 6 anion neighbors), 1-based indices.

  • select (list[int] | slice | None, optional) – Frame selection. None selects the last 50 % of frames. Defaults to None.

Returns:

Tilt angles in degrees. Shape (n_frames, n_centers, 3) for multiple frames (columns: x-tilt, y-tilt, z-tilt), or (n_centers, 3) for a single frame.

Return type:

np.ndarray

Raises:

AssertionError – If nl_bo does not have exactly 6 neighbors per center.

ferrodispcalc.compute.calculate_polarization(traj, nl_ba: ndarray, nl_bx: ndarray, born_effective_charge: dict, select: list[int] | slice | None = None) ndarray[source]

Calculate local polarization for each unit cell in a perovskite.

The local polarization at each B-site unit cell is computed from the positions of the B-site cation, the surrounding A-site cations (8 neighbors), and the anions (6 neighbors), weighted by Born effective charges. Results are returned in C/m².

Parameters:
  • traj (list[Atoms] | Atoms) – Full MD trajectory or a single ASE Atoms object.

  • nl_ba (np.ndarray) – Neighbor list for B–A pairs (B-site center, 8 A-site neighbors), shape (n_cells, 9), 1-based indices.

  • nl_bx (np.ndarray) – Neighbor list for B–X pairs (B-site center, 6 anion neighbors), shape (n_cells, 7), 1-based indices.

  • born_effective_charge (dict) – Mapping of element symbol to scalar Born effective charge, e.g. {'Pb': 3.44, 'Ti': 5.18, 'O': -2.87}.

  • select (list[int] | slice | None, optional) – Frame selection. None selects the last 50 % of frames. Defaults to None.

Returns:

Local polarization in C/m². Shape (n_frames, n_cells, 3) for multiple frames, or (n_cells, 3) for a single frame.

Return type:

np.ndarray

Raises:

AssertionError – If the B-site center indices in nl_ba and nl_bx do not match, or if the neighbor counts are not 8 (A-site) and 6 (X-site).

ferrodispcalc.neighborlist

ferrodispcalc.neighborlist.build_neighbor_list(atoms: Atoms, center_elements: list[str], neighbor_elements: list[str], cutoff: float, neighbor_num: int, defect: bool = False) ndarray[source]

Build neighbor lists for selected atom types.

Find the neighbor of selected center atoms within a cutoff radius. And return the 1-based neighbor list array.

Parameters:
  • atoms (ase.Atoms) – ASE Atoms object containing atomic positions, cell and element types.

  • center_elements (list[str]) – Element symbols to treat as centers (e.g. [“Ti”, “Ba”]).

  • neighbor_elements (list[str]) – Element symbols considered as neighbors (e.g. [“O”]).

  • cutoff (float) – Cutoff radius used to search for neighbors (Unit: Angstrom).

  • neighbor_num (int) – Number of neighbors to keep per center. If fewer neighbors are found and defect is False a ValueError is raised; if defect is True missing neighbor slots are filled with the center index.

  • defect (bool, optional) – When True, allow centers with fewer than neighbor_num neighbors and fill missing entries with the center index (default: False).

Returns:

Integer array of shape (n_centers, neighbor_num + 1). Each row contains the center atom index followed by the neighbor atom indices. Indices are 1-based.

Return type:

np.ndarray

Raises:

ValueError – If any center atom has fewer than neighbor_num neighbors and defect is False.

Notes

  • Internally the function uses pymatgen.core.Structure.get_neighbor_list().

  • A small-cell check prints a warning if any cell vector length is less than 4.0 Å.

  • Ensure that the provided atoms object has periodic boundary conditions (PBC) and a valid cell when you expect PBC behavior.

Examples

>>> from ase.io import read
>>> from ferrodispcalc.neighborlist import build_neighbor_list
>>> atoms = read("POSCAR")
>>> nl = build_neighbor_list(
...     atoms,
...     center_elements=['Pb', 'Sr'],
...     neighbor_elements=['O'],
...     cutoff=4.0,
...     neighbor_num=12,
...     defect=False,
... )
ferrodispcalc.neighborlist.save_neighbor_list(nl: ndarray, file_name: str, zero_based: bool = False) None[source]

Save a neighbor list array to a text file.

Parameters:
  • nl (np.ndarray) – Neighbor list array as returned by build_neighbor_list(). The function expects an integer array where each row contains a center index followed by neighbor indices.

  • file_name (str) – Output file path. The neighbor list will be written as a text file with fixed-width integer columns.

  • zero_based (bool, optional) – If True, convert indices to zero-based before saving. By default the neighbor list uses 1-based indices and zero_based is False.

Examples

>>> from ferrodispcalc.neighborlist import save_neighbor_list, build_neighbor_list
>>> nl = build_neighbor_list(...)
>>> save_neighbor_list(nl, 'nl.dat', zero_based=False)

ferrodispcalc.io

ferrodispcalc.io.read_lammps_data(filename: str, type_map: list[str]) Atoms[source]

Read a LAMMPS data file and return the structure as an ASE Atoms object.

Parameters:
  • filename (str) – Path to the LAMMPS data file.

  • type_map (list[str]) – Ordered list of element symbols corresponding to LAMMPS atom types, e.g. ['Pb', 'Ti', 'O'].

Returns:

The structure described in the data file.

Return type:

Atoms

ferrodispcalc.io.read_lammps_dump(filename: str, type_map: List[str], select: int | slice | List[int] | None = None, cache: bool = True) Atoms | List[Atoms] | None[source]

Read a LAMMPS dump trajectory file using the high-performance C++ backend.

The reader memory-maps the file and optionally caches the frame index to a sidecar file (.{filename}.idx) for fast subsequent reads.

Note

Atom coordinates are read as stored in the dump file and are not shifted back to the primary unit cell. A constant offset may therefore exist compared to tools such as dpdata, but this does not affect any relative-displacement calculations.

Parameters:
  • filename (str) – Path to the LAMMPS dump file.

  • type_map (list[str]) – Ordered list of element symbols corresponding to LAMMPS atom types, e.g. ['Pb', 'Ti', 'O'].

  • select (int | slice | list[int] | None, optional) –

    Frames to read.

    • None – read all frames (default).

    • int – read a single frame by index (negative indices supported); returns a single Atoms object instead of a list.

    • slice – read a contiguous or strided range of frames.

    • list[int] – read an arbitrary set of frame indices.

  • cache (bool, optional) – If True (default), write/read a sidecar index file so that re-reading the same dump file is faster.

Returns:

A single Atoms when select is an integer, otherwise a list of Atoms. Returns None (integer select) or [] (other select) if the file is empty.

Return type:

Atoms | list[Atoms] | None

Raises:
  • FileNotFoundError – If filename does not exist.

  • IOError – If the C++ backend fails to open or memory-map the file.

  • IndexError – If an integer select is out of range.

ferrodispcalc.io.read_lammps_log(file_name: str) Dict[source]

Parse thermodynamic data from the last run in a LAMMPS log file.

The function locates the final Per MPI rank memory header and reads all thermo lines up to the matching Loop time of footer.

Parameters:

file_name (str) – Path to the LAMMPS log file.

Returns:

Dictionary mapping each thermo keyword (e.g. 'Step', 'Temp', 'Press') to a NumPy array of values. An extra key 'nframes' holds the number of thermo records read.

Return type:

dict

Raises:

AssertionError – If the number of header columns does not match the data columns.

ferrodispcalc.config

Used for common configuration, including type map and born effective charges.

Author: Denan Li Last modified: 2024-07-15 Email: lidenan@westlake.edu.cn

class ferrodispcalc.config.TypeMap(UniPero: Tuple[str] = ('Ba', 'Pb', 'Ca', 'Sr', 'Bi', 'K', 'Na', 'Hf', 'Ti', 'Zr', 'Nb', 'Mg', 'In', 'Zn', 'O'), PSTO: Tuple[str] = ('Sr', 'Pb', 'Ti', 'O'))[source]

Bases: object

PSTO: Tuple[str] = ('Sr', 'Pb', 'Ti', 'O')
UniPero: Tuple[str] = ('Ba', 'Pb', 'Ca', 'Sr', 'Bi', 'K', 'Na', 'Hf', 'Ti', 'Zr', 'Nb', 'Mg', 'In', 'Zn', 'O')

ferrodispcalc.sed

class ferrodispcalc.sed.DipoleModeResult(mode_movie: ndarray, primitive_mode_movie: ndarray, phases: ndarray, evec_basis: ndarray, amplitude_basis: ndarray, qpoint: ndarray, requested_freq_THz: float, actual_freq_THz: float, primitive_shape: tuple[int, int, int], cell_shape: tuple[int, int, int], freq_method: str, gauge: str, normalize: bool, remove_mean: bool)[source]

Bases: object

Result of a basis-resolved local-mode extraction.

mode_movie

Real-space mode movie with shape (nphase, nx, ny, nz, 3).

Type:

np.ndarray

primitive_mode_movie

Representative primitive-block movie with shape (nphase, px, py, pz, 3).

Type:

np.ndarray

phases

Phase samples in radians with shape (nphase,).

Type:

np.ndarray

evec_basis

Complex basis-resolved eigenvector with shape (px, py, pz, 3).

Type:

np.ndarray

amplitude_basis

Complex, unnormalized Fourier amplitude after gauge fixing, with shape (px, py, pz, 3).

Type:

np.ndarray

qpoint

Reduced-coordinate q-point used for extraction.

Type:

np.ndarray

requested_freq_THz

Target frequency requested by the user.

Type:

float

actual_freq_THz

Frequency actually used. This can differ from the requested value for freq_method="nearest_fft".

Type:

float

primitive_shape

Primitive block size used to group the input grid.

Type:

tuple[int, int, int]

cell_shape

Number of primitive cells along each grid direction.

Type:

tuple[int, int, int]

freq_method

Frequency extraction method, either "nearest_fft" or "direct".

Type:

str

gauge

Gauge-fixing method, either "max_real" or "none".

Type:

str

normalize

Whether evec_basis was normalized by sqrt(sum(abs(A)**2)).

Type:

bool

remove_mean

Whether per-grid-point time means were removed before extraction.

Type:

bool

actual_freq_THz: float
amplitude_basis: ndarray
cell_shape: tuple[int, int, int]
evec_basis: ndarray
freq_method: str
gauge: str
mode_movie: ndarray
normalize: bool
phases: ndarray
primitive_mode_movie: ndarray
primitive_shape: tuple[int, int, int]
qpoint: ndarray
remove_mean: bool
requested_freq_THz: float
class ferrodispcalc.sed.DipoleSedResult(freq_THz: ndarray, qpoints: ndarray, sed: ndarray, primitive_shape: tuple[int, int, int], cell_shape: tuple[int, int, int], basis_summation: str = 'incoherent')[source]

Bases: object

Result of a local-mode SED calculation.

freq_THz

Positive-frequency axis in THz.

Type:

np.ndarray

qpoints

Reduced q-points used in the spatial Fourier transform.

Type:

np.ndarray

sed

SED intensity with shape (nfreq, nq, 4). The last axis stores x, y, z, and total.

Type:

np.ndarray

primitive_shape

Primitive block size used to group the input grid.

Type:

tuple[int, int, int]

cell_shape

Number of primitive cells along each grid direction.

Type:

tuple[int, int, int]

basis_summation

Basis summation rule. "incoherent" sums intensities over basis degrees of freedom; "coherent" sums amplitudes before squaring.

Type:

str

basis_summation: str = 'incoherent'
cell_shape: tuple[int, int, int]
freq_THz: ndarray
primitive_shape: tuple[int, int, int]
qpoints: ndarray
sed: ndarray
ferrodispcalc.sed.calculate_sed(field: ndarray, dt_ps: float, qpoints: ndarray, primitive_shape: tuple[int, int, int] = (1, 1, 1), num_splits: int = 1, remove_mean: bool = False, n_jobs: int = 1, parallel_backend: str = 'threading', *, basis_summation: str = 'incoherent') DipoleSedResult[source]

Calculate a local-mode fluctuation spectrum for a gridded trajectory.

The input field is assumed to be fixed on a regular grid of local modes, such as Ti-centered displacement, local dipole, local polarization, or the corresponding time derivative. Users must provide the reduced q-points explicitly; use generate_commensurate_qpath() when a path should be sampled by q-points commensurate with the simulation cell.

Parameters:
  • field (np.ndarray) – Input trajectory with shape (nframe, nx, ny, nz, 3). The last axis stores the Cartesian x, y, and z components. If dt_ps is in ps and the field is from LAMMPS metal units, displacement-like input is typically in Angstrom and velocity-like input is in Angstrom/ps.

  • dt_ps (float) – Time interval between stored frames in ps. The returned frequency axis is in THz because 1 / ps = 1 THz.

  • qpoints (np.ndarray) – Reduced q-points with shape (nq, 3). The spatial phase is exp(+i * 2*pi * dot(q, cell_index + basis_offset)), where basis_offset is the fractional position inside the primitive block. The opposite sign gives the same intensity after taking |A|^2.

  • primitive_shape (tuple[int, int, int], optional) – Primitive block size in grid units. For example, (1, 1, 5) groups five local modes into each primitive cell along z. Defaults to (1, 1, 1).

  • num_splits (int, optional) – Number of equal time blocks used for block averaging. nframe must be divisible by this value. Defaults to 1.

  • remove_mean (bool, optional) – If True, subtract the per-block time average at every grid point and component before the FFT. This is usually appropriate for a displacement-like field to remove the static/DC component. Velocity-like fields are usually used as-is. Defaults to False.

  • n_jobs (int, optional) – Number of q-points to evaluate in parallel with joblib. 1 means serial execution; negative values follow joblib’s convention. Defaults to 1.

  • parallel_backend (str, optional) – Joblib backend. The default "threading" avoids copying the large block array into worker processes.

  • basis_summation ({"incoherent", "coherent"}, optional) – Basis summation rule. "incoherent" computes sum_b |U_b(q, omega)|^2 and is the default local-fluctuation spectrum. "coherent" computes |sum_b U_b(q, omega)|^2 and preserves basis interference, closer to a structure-factor-like representation.

Returns:

Frequency axis, q-points, SED intensity, and grid metadata.

Return type:

DipoleSedResult

Notes

The normalization is 1 / (Nt * Ncell) before averaging over time blocks, where Nt is the number of frames per block and Ncell is the number of primitive cells. No mass weighting is applied, so this is not a strict mass-weighted phonon SED unless the input field and normalization are defined accordingly. The spatial Fourier transform uses the full local-mode position within the current folded-cell representation and, by default, sums basis intensities incoherently.

ferrodispcalc.sed.extract_eigen_vector(field: ndarray, dt_ps: float, qpoint: ndarray, freq_THz: float, primitive_shape: tuple[int, int, int] = (1, 1, 1), nphase: int = 24, remove_mean: bool = True, freq_method: str = 'nearest_fft', gauge: str = 'max_real', normalize: bool = True) DipoleModeResult[source]

Extract a basis-resolved complex mode and rebuild a real-space movie.

The extraction uses the full local-mode position, cell_index + basis_offset, in the spatial Fourier transform. It keeps the complex amplitude before any basis summation or absolute-square operation, then reconstructs the selected (q, frequency) mode over one artificial vibration period.

Parameters:
  • field (np.ndarray) – Input trajectory with shape (nframe, nx, ny, nz, 3).

  • dt_ps (float) – Time interval between stored frames in ps.

  • qpoint (np.ndarray) – Single reduced-coordinate q-point with shape (3,).

  • freq_THz (float) – Target frequency in THz.

  • primitive_shape (tuple[int, int, int], optional) – Primitive block size in grid units. This should match the SED calculation used to identify the selected peak.

  • nphase (int, optional) – Number of phase frames in the returned one-period movie.

  • remove_mean (bool, optional) – If True, remove the per-grid-point time average before extracting the mode. Defaults to True.

  • freq_method ({"nearest_fft", "direct"}, optional) – "nearest_fft" uses the closest non-negative FFT frequency bin. "direct" projects the trajectory directly at freq_THz.

  • gauge ({"max_real", "none"}, optional) – "max_real" rotates the global complex phase so the largest eigenvector component is positive real. "none" leaves the phase as extracted.

  • normalize (bool, optional) – If True, normalize the eigenvector by sqrt(sum(abs(A)**2)). Defaults to True.

Returns:

Mode movies, basis eigenvector, raw basis amplitude, and metadata.

Return type:

DipoleModeResult

Notes

The returned mode_movie and primitive_mode_movie are normalized local-mode patterns when normalize=True. Their amplitudes are therefore not the physical MD vibration amplitude; use amplitude_basis or the SED intensity to judge mode strength.

ferrodispcalc.sed.generate_commensurate_qpath(q_path: ndarray, cell_shape: tuple[int, int, int]) tuple[ndarray, ndarray][source]

Generate q-points along a path that is commensurate with a cell grid.

For cell shape (Nx, Ny, Nz), allowed reduced q-points satisfy q * (Nx, Ny, Nz) being integer-valued.

Returns:

Reduced q-points and cumulative reduced-coordinate distances.

Return type:

qpoints, q_distances

ferrodispcalc.sed.load_eigen_vector(filename: str | Path) DipoleModeResult[source]

Load a DipoleModeResult saved by save_eigen_vector().

ferrodispcalc.sed.load_sed(filename: str | Path) DipoleSedResult[source]

Load a DipoleSedResult saved by save_sed().

ferrodispcalc.sed.plot_sed(result: DipoleSedResult | ndarray, freq_THz: ndarray | None = None, qpoints: ndarray | None = None, q_distances: ndarray | None = None, component: str = 'total', ax: Axes | None = None, q_labels: tuple[str, str] = ('$\\Gamma$', 'X'), freq_max: float | None = None, freq_interval: float = 5.0, cmap: str = 'RdBu_r', colorbar_min: float | None = None, colorbar_max: float | None = None, use_contourf: bool = False, savepath: str | Path | None = None, show: bool = False) tuple[Figure, Axes][source]

Plot one SED component.

Parameters:
  • result (DipoleSedResult | np.ndarray) – A DipoleSedResult, or a raw SED array with shape (nfreq, nq, 4) or (nfreq, nq). If an array is provided, freq_THz and qpoints must also be provided.

  • freq_THz (np.ndarray | None, optional) – Frequency axis in THz. Ignored when result is a DipoleSedResult.

  • qpoints (np.ndarray | None, optional) – Reduced q-points with shape (nq, 3). Ignored when result is a DipoleSedResult.

  • q_distances (np.ndarray | None, optional) – Cumulative q-path distances. If not provided, they are reconstructed from neighboring q-point differences.

  • component (str, optional) – One of "x", "y", "z", or "total". Defaults to "total".

  • ax (matplotlib.axes.Axes | None, optional) – Existing axes. A new figure and axes are created when None.

  • q_labels (tuple[str, str], optional) – Labels for the first and last q-points. Defaults to Gamma-X labels.

  • freq_max (float | None, optional) – Maximum plotted frequency in THz. Defaults to the maximum available frequency.

  • freq_interval (float, optional) – Frequency tick interval in THz. Defaults to 5.0.

  • cmap (str, optional) – Matplotlib colormap. Defaults to "RdBu_r".

  • colorbar_min (float | None, optional) – Colorbar limits on the log-intensity scale.

  • colorbar_max (float | None, optional) – Colorbar limits on the log-intensity scale.

  • use_contourf (bool, optional) – Use contourf instead of the default imshow rendering.

  • savepath (str | Path | None, optional) – Save the figure when provided.

  • show (bool, optional) – Show the figure interactively before returning.

Returns:

Matplotlib figure and axes.

Return type:

fig, ax

ferrodispcalc.sed.plot_sed_1d(result: DipoleSedResult, qpoint: int | Sequence[float] | ndarray, component: str = 'total', ax: Axes | None = None, style: dict | None = None, file_name: str | Path | None = None, data_only: bool = False) Axes | tuple[ndarray, ndarray][source]

Plot the SED spectrum at one exact q-point.

Parameters:
  • result (DipoleSedResult) – A DipoleSedResult containing the SED intensity, frequency axis, and q-points.

  • qpoint (int | sequence of float | np.ndarray) – Q-point selector. An integer is interpreted as a zero-based q-point index. A length-3 array-like value is interpreted as a reduced q-point coordinate and must match exactly one row in qpoints; no nearest-neighbor or tolerance-based matching is applied.

  • component (str, optional) – Component plotted when data_only is False. One of "total", "x", "y", or "z". Defaults to "total".

  • ax (matplotlib.axes.Axes | None, optional) – Existing axes. A new figure and axes are created when None.

  • style (dict | None, optional) – Matplotlib style arguments passed directly to ax.plot.

  • file_name (str | pathlib.Path | None, optional) – Save the figure to this path when provided. The figure is saved with dpi=650 and bbox_inches="tight".

  • data_only (bool, optional) – If True, return the selected q-point data instead of plotting. The returned intensity array has shape (nfreq, 4) with columns total, x, y, and z. This requires full component SED data with shape (nfreq, nq, 4).

Returns:

  • ax (matplotlib.axes.Axes) – Returned when data_only is False.

  • freq_THz, intensity (tuple[np.ndarray, np.ndarray]) – Returned when data_only is True.

ferrodispcalc.sed.save_eigen_vector(result: DipoleModeResult, filename: str | Path) None[source]

Save an extracted mode result to a compressed npz file.

ferrodispcalc.sed.save_sed(result: DipoleSedResult, filename: str | Path) None[source]

Save SED arrays and grid metadata to a compressed npz file.

ferrodispcalc.vis

ferrodispcalc.vis.grid

ferrodispcalc.vis.grid.grid_data(atoms: Atoms, data: ndarray, element: list[str], target_size: tuple[int, int, int] | None = None, tol: float = 1.0, axis: tuple[tuple, tuple, tuple] = ((1, 0, 0), (0, 1, 0), (0, 0, 1)), return_coord: bool = False) ndarray | tuple[ndarray, ndarray][source]

ferrodispcalc.vis.line_plot

ferrodispcalc.vis.line_plot.line_profile(data: ndarray, ax: Axes = None, along: str = 'x', savepath: str = None, field_prefix: str = None, call_back: callable = None)[source]

ferrodispcalc.vis.plane_plot

ferrodispcalc.vis.plane_plot.plane_profile(data: ndarray, save_dir: str | Path = 'plane_plots', relative: bool = True, select: Dict[str, List[int]] | None = None, hot: bool = True, fig_size: str | Tuple[float, float] = 'auto') None[source]

Plot 2D cross-sections of 3D vector data.

Parameters:

data: np.ndarray

The 4D vector data in (nx, ny, nz, 3) format.

save_dir: str | Path

The directory to save the plots.

relative: bool

Whether to plot the vector in the relative scale (default True).

select: dict, optional

The selected layers to plot. Keys should be ‘x’, ‘y’, ‘z’ (for YZ, XZ, XY planes respectively). e.g., {‘z’: [0, 1], ‘y’: [0]} If None, all layers are plotted.

hot: bool

Whether to plot the angle in the hot colormap (default True).

fig_size: str | tuple

The figure size: ‘auto’, ‘default’, or a (width, height) tuple.

ferrodispcalc.vis.plane_plot.plane_profile_multi(data: ndarray, select: List[str | int | slice | List[int]], save_name: str | Path = 'plane_plots.gif', relative: bool = True, hot: bool = True, fig_size: str | Tuple[float, float] = 'auto', fps: int = 10)[source]

Plot 2D cross-sections of 5D vector data (nframe, nx, ny, nz, 3) and save as a GIF using matplotlib.animation. (使用 matplotlib.animation 绘制 5D 矢量数据的 2D 横截面并另存为 GIF)

Parameters:

data: np.ndarray

The 5D vector data in (nframe, nx, ny, nz, 3) format.

select: list

Defines the plot. Format: [plane_name, layer_index, frame_spec] e.g., [‘xz’, 0, slice(0, None, 2)] -> XZ plane, y_index=0, every 2nd frame.

save_name: str | Path

The output filename (e.g., ‘animation.gif’).

relative: bool

Whether to plot the vector in the relative scale.

hot: bool

Whether to plot the angle in the hot colormap.

fig_size: str | tuple

The figure size: ‘auto’, ‘default’, or a (width, height) tuple.

fps: int

Frames per second for the output GIF.

ferrodispcalc.vis.space_plot

ferrodispcalc.vis.space_plot.space_animation(data: np.ndarray, coord: np.ndarray | None = None, *, color_by: str = 'dz', cmap: str = 'coolwarm', factor: float = 25.0, projection: str = 'ortho', clim: tuple[float, float] | None=None, select: dict | None = None, cell: np.ndarray | None = None, title: str = '3D Vector Field Animation', show_bounding_box: bool = True, show_axes: bool = True, show_slider: bool = True, show_frame_text: bool = True, plotter: pv.Plotter | None = None, stride: int | Sequence[int] = 1, frame_indices: Sequence[int] | slice | None = None, frame_step: int = 1, fps: float = 20.0, autoplay: bool = True, loop: bool = True, dtype: type | str | np.dtype = <class 'numpy.float32'>) pv.Plotter[source]

Play a 3-D vector-field animation in a PyVista window.

Two time-dependent input modes are supported:

  • Grid mode - data has shape (nframe, nx, ny, nz, 3). Coordinates are generated automatically from grid indices and coord is ignored.

  • Point mode - data has shape (nframe, npoint, 3) and coord has shape (npoint, 3).

The scene is created once and updated in-place as frames advance. A slider widget is shown by default for direct frame selection. Keyboard shortcuts are available when more than one frame is selected: space toggles playback, and the left/right arrow keys step backward/forward by one displayed frame.

Parameters:
  • data (np.ndarray) – Vector animation data. Expected shape is (nframe, nx, ny, nz, 3) for grid mode or (nframe, npoint, 3) for point mode.

  • coord (np.ndarray or None) – Cartesian coordinates (npoint, 3). Required for point mode, ignored for grid mode.

  • color_by (str) – Coloring strategy: 'magnitude', 'dx', 'dy', 'dz', 'all' (RGB from components), or 'gradient' (grid mode only).

  • cmap (str) – Matplotlib colormap name (ignored when color_by=’all’).

  • factor (float) – Arrow scale factor passed to pv.PolyData.glyph.

  • projection (str) – 'ortho' for parallel projection, 'persp' for perspective.

  • clim (tuple or None) – Colorbar range (vmin, vmax). None updates the scalar range from the current frame.

  • select (dict or None) – Region filter in fractional coordinates. Keys are 'x', 'y', 'z'; values are [lo, hi] or None (no filter). In grid mode fractions are computed from grid indices; in point mode cell is required.

  • cell (np.ndarray or None) – (3, 3) lattice matrix (rows = lattice vectors). Only needed when select is used in point mode.

  • title (str) – Window title.

  • show_bounding_box (bool) – Draw a wireframe box around the full point cloud.

  • show_axes (bool) – Show orientation axes widget.

  • show_slider (bool) – Show a frame slider at the bottom of the window.

  • show_frame_text (bool) – Show the current frame index in the upper-left corner.

  • plotter (pv.Plotter or None) – Reuse an existing plotter. A new one is created when None.

  • stride (int or sequence of int) – Spatial down-sampling. In grid mode, an int applies to all axes and a 3-item sequence applies to x, y and z. In point mode, only an int stride is accepted.

  • frame_indices (sequence of int, slice, or None) – Frames to expose in the animation. None selects all frames after applying frame_step.

  • frame_step (int) – Temporal down-sampling step used when frame_indices is None; also applied to a slice selection when greater than one.

  • fps (float) – Playback rate in displayed frames per second.

  • autoplay (bool) – Start playback immediately when the window opens.

  • loop (bool) – Restart from the first displayed frame after the last one.

  • dtype (type, str, or np.dtype) – Floating dtype used for the per-frame vector arrays passed to PyVista.

Returns:

The plotter instance. When plotter is None, the window is shown before returning.

Return type:

pv.Plotter

ferrodispcalc.vis.space_plot.space_profile(data: np.ndarray, coord: np.ndarray | None = None, *, color_by: str = 'dz', cmap: str = 'coolwarm', factor: float = 25.0, projection: str = 'ortho', clim: tuple[float, float] | None = None, select: dict | None = None, cell: np.ndarray | None = None, title: str = '3D Vector Field', show_bounding_box: bool = True, show_axes: bool = True, plotter: pv.Plotter | None = None) pv.Plotter[source]

Plot a 3D vector field using PyVista.

Two input modes are supported:

  • Grid modedata has shape (nx, ny, nz, 3). Coordinates are generated automatically from grid indices and coord is ignored.

  • Point modedata has shape (npoint, 3) and coord (also (npoint, 3)) must be provided.

Parameters:
  • data (np.ndarray) – Vector data. (nx, ny, nz, 3) for grid mode or (npoint, 3) for point mode.

  • coord (np.ndarray or None) – Cartesian coordinates (npoint, 3). Required for point mode, ignored for grid mode.

  • color_by (str) – Coloring strategy: 'magnitude', 'dx', 'dy', 'dz', 'all' (RGB from components), or 'gradient' (grid mode only).

  • cmap (str) – Matplotlib colormap name (ignored when color_by=’all’).

  • factor (float) – Arrow scale factor passed to pv.PolyData.glyph.

  • projection (str) – 'ortho' for parallel projection, 'persp' for perspective.

  • clim (tuple or None) – Colorbar range (vmin, vmax). None for automatic.

  • select (dict or None) – Region filter in fractional coordinates. Keys are 'x', 'y', 'z'; values are [lo, hi] or None (no filter). In grid mode fractions are computed from grid indices; in point mode cell is required.

  • cell (np.ndarray or None) – (3, 3) lattice matrix (rows = lattice vectors). Only needed when select is used in point mode.

  • title (str) – Window title.

  • show_bounding_box (bool) – Draw a wireframe box around the full point cloud.

  • show_axes (bool) – Show orientation axes widget.

  • plotter (pv.Plotter or None) – Reuse an existing plotter. A new one is created when None.

Returns:

The plotter instance (already shown).

Return type:

pv.Plotter