PowerSpectrum

class rascaline.torch.utils.PowerSpectrum(calculator_1: CalculatorModule, calculator_2: CalculatorModule | None = None, types: List[int] | None = None)

Bases: Module

General power spectrum of one or of two calculators.

If calculator_2 is provided, the invariants \(p_{nl}\) are generated by taking quadratic combinations of calculator_1’s spherical expansion \(\rho_{nlm}\) and calculator_2’s spherical expansion \(\nu_{nlm}\) according to Bartók et. al.

\[p_{nl} = \rho_{nlm}^\dagger \cdot \nu_{nlm}\]

where we use the Einstein summation convention. If gradients are present the invariants of those are constructed as

\[\nabla p_{nl} = \nabla \rho_{nlm}^\dagger \cdot \nu_{nlm} + \rho_{nlm}^\dagger \cdot \nabla \nu_{nlm}\]

Note

Currently only supports gradients with respect to positions.

If calculator_2=None invariants are generated by combining the coefficients of the spherical expansion of calculator_1. The spherical expansions given as input can only be rascaline.SphericalExpansion or rascaline.LodeSphericalExpansion.

Parameters:
  • calculator_1 – first calculator

  • calculator_1 – second calculator

  • types – List of "neighbor_type" to use in the properties of the output. This option might be useful when running the calculation on subset of a whole dataset and trying to join along the sample dimension after the calculation. If None, blocks are filled with "neighbor_type" found in the systems.

Raises:

Example

As an example we calculate the power spectrum for a short range (sr) spherical expansion and a long-range (lr) LODE spherical expansion for a NaCl crystal.

>>> import rascaline
>>> import ase

Construct the NaCl crystal

>>> atoms = ase.Atoms(
...     symbols="NaCl",
...     positions=[[0, 0, 0], [0.5, 0.5, 0.5]],
...     pbc=True,
...     cell=[1, 1, 1],
... )

Define the hyper parameters for the short-range spherical expansion

>>> sr_hypers = {
...     "cutoff": 1.0,
...     "max_radial": 6,
...     "max_angular": 2,
...     "atomic_gaussian_width": 0.3,
...     "center_atom_weight": 1.0,
...     "radial_basis": {
...         "Gto": {},
...     },
...     "cutoff_function": {
...         "ShiftedCosine": {"width": 0.5},
...     },
... }

Define the hyper parameters for the long-range LODE spherical expansion from the hyper parameters of the short-range spherical expansion

>>> lr_hypers = sr_hypers.copy()
>>> lr_hypers.pop("cutoff_function")
{'ShiftedCosine': {'width': 0.5}}
>>> lr_hypers["potential_exponent"] = 1

Construct the calculators

>>> sr_calculator = rascaline.SphericalExpansion(**sr_hypers)
>>> lr_calculator = rascaline.LodeSphericalExpansion(**lr_hypers)

Construct the power spectrum calculators and compute the spherical expansion

>>> calculator = rascaline.utils.PowerSpectrum(sr_calculator, lr_calculator)
>>> power_spectrum = calculator.compute(atoms)

The resulting invariants are stored as metatensor.TensorMap as for any other calculator

>>> power_spectrum.keys
Labels(
    center_type
        11
        17
)
>>> power_spectrum[0]
TensorBlock
    samples (1): ['system', 'atom']
    components (): []
    properties (432): ['l', 'neighbor_1_type', 'n_1', 'neighbor_2_type', 'n_2']
    gradients: None

See also

If you are interested in the SOAP power spectrum you can the use the faster rascaline.SoapPowerSpectrum.

Initialize internal Module state, shared by both nn.Module and ScriptModule.

property name

Name of this calculator.

compute(systems: System | List[System], gradients: List[str] | None = None, use_native_system: bool = True) TensorMap

Runs a calculation with this calculator on the given systems.

See rascaline.calculators.CalculatorBase.compute() for details on the parameters.

Raises:

NotImplementedError – If a spherical expansions contains a gradient with respect to an unknwon parameter.

forward(systems: System | List[System], gradients: List[str] | None = None, use_native_system: bool = True) TensorMap

Calls the PowerSpectrum.compute() function.

This is intended for torch.nn.Module compatibility, and should be ignored in pure Python mode.