Skip to main content

Overview

The Elasticity task computes the elastic tensor of a crystal by:
  1. Performing a full cell relaxation with the OPT task.
  2. Generating a set of deformed structures via pymatgen.analysis.elasticity.DeformedStructureSet using the specified normal and shear strains.
  3. Evaluating the stress tensor of each deformed structure with the supplied calculator.
  4. Fitting the elastic constants (C_ij) by linear regression of stress vs. strain.
The task is registered in Prefect as Elasticity with a TASK_SOURCE + INPUTS cache policy.

Function signature

from mlip_arena.tasks.elasticity import run as ELASTICITY

result = ELASTICITY(
    atoms,
    calculator,
    optimizer="BFGSLineSearch",
    optimizer_kwargs=None,
    filter="FrechetCell",
    filter_kwargs=None,
    criterion=None,
    normal_strains=np.linspace(-0.01, 0.01, 4),
    shear_strains=np.linspace(-0.06, 0.06, 4),
    persist_opt=True,
    cache_opt=False,
)

Parameters

body.atoms
ase.Atoms
required
Input crystal structure. A copy is made internally.
body.calculator
ase.calculators.calculator.BaseCalculator
required
ASE-compatible calculator for stress and force evaluations.
body.optimizer
Optimizer | str
default:"BFGSLineSearch"
Optimizer used for the initial relaxation. See OPT for accepted string values.
body.optimizer_kwargs
dict | None
default:"None"
Extra keyword arguments forwarded to the optimizer constructor.
body.filter
Filter | str | None
default:"FrechetCell"
Filter applied during the initial relaxation to allow cell optimization. See OPT for accepted string values.
body.filter_kwargs
dict | None
default:"None"
Extra keyword arguments forwarded to the filter constructor.
body.criterion
dict | None
default:"None"
Convergence criteria passed to optimizer.run(). See OPT for details.
body.normal_strains
list[number] | np.ndarray | None
default:"np.linspace(-0.01, 0.01, 4)"
Normal strain values applied to generate deformed structures. Defaults to four values uniformly spaced from −0.01 to +0.01 (i.e., ±1%).
body.shear_strains
list[number] | np.ndarray | None
default:"np.linspace(-0.06, 0.06, 4)"
Shear strain values applied to generate deformed structures. Defaults to four values uniformly spaced from −0.06 to +0.06 (i.e., ±6%).
body.persist_opt
boolean
default:"true"
When True, the result of the initial OPT sub-task is persisted (saved to Prefect result storage). Useful for inspecting the relaxed structure independently of the elasticity result.
body.cache_opt
boolean
default:"false"
When True, the initial OPT sub-task uses caching so that identical relaxations are not recomputed across calls. When False, the cache is refreshed every time.

Return value

{
    "elastic_tensor": ElasticTensor,
    "residuals_sum":  float,
}
elastic_tensor
pymatgen.analysis.elasticity.ElasticTensor
The fitted elastic tensor (Voigt notation, 6×6) with near-zero components zeroed out at a tolerance of 1e-7. The ElasticTensor object provides properties including:
et = result["elastic_tensor"]

# Bulk and shear moduli (GPa)
print(et.k_voigt, et.g_voigt)  # Voigt average
print(et.k_reuss, et.g_reuss)  # Reuss average
print(et.k_vrh,   et.g_vrh)    # VRH average

# Young's modulus and Poisson ratio
print(et.y_mod)       # eV/ų — convert: × 160.218 to GPa
print(et.universal_anisotropy)
residuals_sum
number
Sum of least-squares fit residuals across all C_ij linear regressions. A lower value indicates a better-quality fit and a more linearly elastic material response in the applied strain range.
If the initial relaxation fails, the function returns a Prefect State object instead of a dict. Check isinstance(result, dict) before accessing keys.

Example

import numpy as np
from ase.build import bulk
from mlip_arena.models import MLIPEnum
from mlip_arena.tasks.elasticity import run as ELASTICITY
from mlip_arena.tasks.utils import get_calculator

# BCC iron
atoms = bulk("Fe", "bcc", a=2.87)
calculator = get_calculator(MLIPEnum["MACE-MP(M)"])

result = ELASTICITY(
    atoms=atoms,
    calculator=calculator,
    optimizer="BFGSLineSearch",
    filter="FrechetCell",
    criterion={"fmax": 0.001},
    normal_strains=np.linspace(-0.01, 0.01, 4),
    shear_strains=np.linspace(-0.06, 0.06, 4),
    persist_opt=True,
    cache_opt=False,
)

et = result["elastic_tensor"]
print(f"K_VRH  = {et.k_vrh:.1f} GPa")
print(f"G_VRH  = {et.g_vrh:.1f} GPa")
print(f"Residuals sum = {result['residuals_sum']:.4f}")

# Full Voigt matrix
print(et.voigt)