Clebsch-Gordan products

class rascaline.utils.DensityCorrelations(max_angular: int, correlation_order: int, angular_cutoff: int | None = None, selected_keys: Labels | List[Labels] | None = None, skip_redundant: bool | List[bool] | None = False, output_selection: bool | List[bool] | None = None, arrays_backend: str | None = None, cg_backend: str | None = None)

Takes iterative Clebsch-Gordan (CG) tensor products of a density descriptor with itself up to the desired correlation order. Returns TensorMap corresponding to the density correlations output from the specified iteration(s).

The input density descriptor necessarily is body order 2 (i.e. correlation order 1), but can be single- or multi-center. The output is a list of density correlations for each iteration specified in output_selection, up to the target order passed in correlation_order. By default only the last correlation (i.e. the correlation of order correlation_order) is returned.

This function is an iterative special case of the more general correlate_tensors(). As a density is being correlated with itself, some redundant CG tensor products can be skipped with the skip_redundant keyword.

Selections on the angular and parity channels at each iteration can also be controlled with arguments angular_cutoff, angular_selection and parity_selection.

Parameters:
  • max_angular – The maximum angular order for which CG coefficients should be computed and stored. This must be large enough to cover the maximum angular order reached in the CG iterations on a density input to the compute() method.

  • correlation_order – The desired correlation order of the output descriptor. Must be >= 1.

  • angular_cutoff – The maximum angular channel to compute at any given CG iteration, applied globally to all iterations until the target correlation order is reached.

  • selected_keysLabels or list of Labels specifying the angular and/or parity channels to output at each iteration. All Labels objects passed here must only contain key names "o3_lambda" and "o3_sigma". If a single Labels object is given, this is applied to the final iteration only. If a list of Labels is given, each is applied to its corresponding iteration. If None is passed, all angular and parity channels are kept at each iteration, with the global angular_cutoff applied if specified.

  • skip_redundant – Whether to skip redundant CG combinations. Defaults to False, which means all combinations are performed. If a list of bool is passed, this is applied to each iteration. If a single bool is passed, this is applied to all iterations.

  • output_selection – A list of bool specifying whether to output a TensorMap for each iteration. If a single bool is passed as True, outputs from all iterations will be returned. If a list of bool is passed, this controls the output at each corresponding iteration. If None is passed, only the final iteration is output.

  • arrays_backend – Determines the array backend, either "numpy" or "torch".

  • cg_backend

    Determines the backend for the CG combination. It can be "python-sparse", or "python-dense". If the CG combination performs on the sparse coefficients, it means that for each (l1, l2, lambda) block the (m1, m2, mu) coefficients are stored in a sparse format only storing the nonzero coefficients. If this is not given, the most optimal choice is determined given available packages and arrays_backend.

    • "python-dense": Uses the python implementation performing the combinations with the dense CG coefficients.

    • "python-sparse": Uses the python implementation performing the combinations with the sparse CG coefficients.

Returns:

A list of TensorMap corresponding to the density correlations output from the specified iterations. If the output from a single iteration is requested, a TensorMap is returned instead.

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

forward(density: TensorMap) TensorMap | List[TensorMap]

Calls the DensityCorrelations.compute() function.

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

compute(density: TensorMap) TensorMap | List[TensorMap]

Computes the density correlations by taking iterative Clebsch-Gordan (CG) tensor products of the input density descriptor with itself.

Parameters:

density – A density descriptor of body order 2 (correlation order 1), in TensorMap format. This may be, for example, a rascaline SphericalExpansion or LodeSphericalExpansion. Alternatively, this could be multi-center descriptor, such as a pair density.

compute_metadata(density: TensorMap) TensorMap | List[TensorMap]

Returns the metadata-only TensorMap that would be output by the function compute() for the same calculator under the same settings, without performing the actual Clebsch-Gordan tensor products.

Parameters:

density – A density descriptor of body order 2 (correlation order 1), in TensorMap format. This may be, for example, a rascaline SphericalExpansion or LodeSphericalExpansion. Alternatively, this could be multi-center descriptor, such as a pair density.

rascaline.utils.cartesian_to_spherical(tensor: TensorMap, components: List[str], keep_l_in_keys: bool | None = None, remove_blocks_threshold: float | None = 1e-09, cg_backend: str | None = None, cg_coefficients: TensorMap | None = None) TensorMap

Transform a tensor of arbitrary rank from cartesian form to a spherical form.

Starting from a tensor on a basis of product of cartesian coordinates, this function computes the same tensor using a basis of spherical harmonics Y^M_L. For example, a rank 1 tensor with a single “xyz” component would be represented as a single L=1 spherical harmonic; while a rank 5 tensor using a product basis ℝ^3 ℝ^3 ℝ^3 ℝ^3 ℝ^3 would require multiple blocks up to L=5 spherical harmonics.

A single TensorBlock in the input might correspond to multiple TensorBlock in the output. The output keys will contain all the dimensions of the input keys, plus o3_lambda (indicating the spherical harmonics degree) and o3_sigma (indicating that this block is a proper- or improper tensor with +1 and -1 respectively). If keep_l_in_keys is True or if the input tensor is a tensor of rank 3 or more, the keys will also contain multiple l_{i} and k_{i} dimensions, which indicate which angular momenta have been coupled together in which order to get this block.

components specifies which ones of the components of the input TensorMap should be transformed from cartesian to spherical. All these components will be replaced in the output by a single o3_mu component, corresponding to the spherical harmonics M.

By default, symmetric tensors will only contain blocks corresponding to o3_sigma=1. This is achieved by checking the norm of the blocks after the full calculation; and dropping any block with a norm below remove_blocks_epsilon. To keep all blocks regardless of their norm, you can set remove_blocks_epsilon=None.

Parameters:
  • tensor – input tensor, using a cartesian product basis

  • components – components of the input tensor to transform into spherical components

  • keep_l_in_keys

    should the output contains the values of angular momenta that where combined together? This defaults to False for rank 1 and 2 tensors, and True for all other tensors.

    Keys named l_{i} correspond to the input components, with l_1 being the last entry in components and l_N the first one. Keys named k_{i} correspond to intermediary spherical components created during the calculation, i.e. a k_{i} used to be o3_lambda.

  • remove_blocks_threshold – Numerical tolerance to use when determining if a block’s norm is zero or not. Blocks with zero norm will be excluded from the output. Set this to None to keep all blocks in the output.

  • cg_backend – Backend to use for Clebsch-Gordan calculations. This can be "python-dense" or "python-sparse" for dense or sparse operations respectively. If None, this is automatically determined.

  • cg_coefficients – Cache containing Clebsch-Gordan coefficients. This is optional except when using this function from TorchScript. The coefficients should be computed with calculate_cg_coefficients(), using the same cg_backend as this function.

Returns:

TensorMap containing spherical components instead of cartesian components.

rascaline.utils.calculate_cg_coefficients(lambda_max: int, cg_backend: str, use_torch: bool = False) TensorMap

Calculates the Clebsch-Gordan coefficients for all possible combination of angular momenta up to lambda_max.

The structure of the returned TensorMap depends on whether the backend used to perform CG tensor products uses sparse or dense operations.

Dense:
  • samples: _, i.e. a dummy sample.

  • components: [m1, m2, mu] on separate components axes, where m1 and m2 are the m component values for the two arrays being combined and mu is the m component value for the resulting array.

  • properties: cg_coefficient = [0]

Sparse:
  • samples: (m1, m2, mu), where m1 and m2 are the m component values for the two arrays being combined and mu is the m component value for the resulting array.

  • components: [], i.e. no components axis.

  • properties: cg_coefficient = [0]

Parameters:
  • lambda_max – maximum angular momentum value to compute CG coefficients for.

  • cg_backend – TODO

  • use_torch – whether torch tensor or numpy arrays should be used to store the coefficients

Returns:

TensorMap of the Clebsch-Gordan coefficients.