Applications of the HQS Spin Mapper package

Outline

Applications of the HQS Spin Mapper software package

The HQS Spin Mapper software package contains a selection of modules which provide two major functionalities. The first functionality is the identification of orbital bases containing spin-like orbitals, i.e. orbitals with an occupancy of strictly a single electron. The second one derives effective spin-bath model Hamiltonians for the determined spin-like orbitals. The remaining modules handle input and output tasks.

Accepted input formats (hqs_spin_mapper.transformables)

HQS Spin Mapper currently accepts several different input formats for the Hamiltonian descriptions of a system (i.e. crystals and molecules). Some of the input formats store additional information about the system, which is, for example, needed for the determination of spin-like orbital bases or to reduce the computational of the spin-bath model derivation. There is also the possibility for the user to provide their own input format through the use of the supplied protocols.

Currently supported by default are:

  • LatticeBuilder configuration object instances
  • LatticeBuilder configuration data files as *.yml
  • Struqture FermionHamiltonianSystem object instances
  • Data container (Transformable_QuantumChemistry) including the one- and two-particle integrals specifying the Hamiltonian stored as numpy binary *.npy files
  • Data container (Transformable_Matrices) including the single particle picture matrices of the Hamiltonian as numpy array instances

For each default input format, we provide a specific class derived from the Transformable class. The Transformable object serves as the necessary input argument of the essential methods of the HQS Spin Mapper package and are manipulated by these methods. The Transformable classes satisfy the Supports_SW_Transformation protocol. Any object satisfying the Supports_SW_Transformation protocol can be used with the HQS Spin Mapper functions.

LatticeBuilder object instances and data files (Transformable_LatticeBuilder)

Given an instance of a LatticeBuilder object, it is easy to instantiate a Transformable_LatticeBuilder object from, which can then be used as the input for the subsequent functions. To instantiate the Transformable_LatticeBuilder object, use:

from hqs_spin_mapper.transformables import Transformable_LatticeBuilder

system = Transformable_LatticeBuilder(input_builder)

If the LatticeBuilder description of the system is stored in a valid *.yml file, the Transformable_LatticeBuilder object is instantiated in the same way.

from hqs_spin_mapper.transformables import Transformable_LatticeBuilder

system = Transformable_LatticeBuilder(input_lattice_builder_yaml)

Struqture FermionHamiltonianSystem instances (Transformable_Struqture)

A Transformable_Struqture object can be directly instantiated from a struqture_py.fermions.FermionHamiltonianSystem instance.

from hqs_spin_mapper.transformables import Transformable_Struqture

system = Transformable_Struqture(input_struqture)

Hamiltonians as numpy array instances (Transformable_Matrices)

The description of the Hamiltonian of a given system can be provided as a collection of matrices or tensors respectively, in the single particle picture. The matrices and tensors first need to be stored in a data container class called Hamiltonian_Matrices.

from hqs_spin_mapper.spin_mapper_protocols import Hamiltonian_Matrices

hamiltonian_matrices = Hamiltonian_Matrices(H0_uu=H0_uu,
                                            H0_dd=H0_dd,
                                            H0_ud=H0_ud,
                                            HU=HU,
                                            Jz=None,
                                            Jp=None,
                                            Jc=None)

The spin indices of the quadratic fermionic terms H0_uu, H0_dd, H0_ud indicate the spin component of the electrons. If provided as a matrix, HU denotes density-density interaction terms. If HU is provided as a tensor, more complex four index interaction terms can be specified. The matrix Jz specifies interactions, Jp specifies interactions, and Jc specifies interactions.

The quadratic Hamiltonian matrices H0_uu, H0_dd, H0_ud have to follow the convention The matrices H0_uu, H0_dd, H0_ud have to be provided for three different spin direction combinations , i.e. , and .

If provided as a matrix, a HU has to conform to and if provided as a tensor as

Using the Hamiltonian_Matrices instance, one can instantiate a Transformable_Matrices object that satisfies the Supports_SW_Transformation protocol via

from hqs_spin_mapper.transformables import Transformable_Matrices

system = Transformable_Matrices(hamiltonian_matrices)

Electron integrals from electronic structure codes (Transformable_QuantumChemistry)

If the output of a preceding electronic structure calculations contains the one-electron and two-electron integrals specifying the Hamiltonian stored as numpy binary *.npy files, one can construct a Transformable from the corresponding matrices via instantiation of a Transformable_Matrices object as shown in the previous section.

from hqs_spin_mapper.spin_mapper_protocols import Hamiltonian_Matrices
from hqs_spin_mapper.transformables import Transformable_Matrices

H0_uu = np.load('h0_uu.npy')
H0_dd = np.load('h0_dd.npy')
H0_ud = np.load('h0_ud.npy')
HU = np.load('hu.npy')

hamiltonian_matrices = Hamiltonian_Matrices(H0_uu=H0_uu,
                                            H0_dd=H0_dd,
                                            H0_ud=H0_ud,
                                            HU=HU,
                                            Jz=None,
                                            Jp=None,
                                            Jc=None)

system = Transformable_Matrices(hamiltonian_matrices)

If the results from the preceding electronic structure calculations also include the one-electron density matrix (1-RDM) and the two-electron density matrix (2-RDM), one can also instantiate a Transformable_QuantumChemistry object. The Transformable_QuantumChemistry object also satisfies the Supports_Spin_Optimization protocol. This means that the object can be used as the input argument for the spin basis optimization functionality of Spin Mapper.

from hqs_spin_mapper.transformables import Transformable_QuantumChemistry
from hqs_spin_mapper.spin_mapper_protocols import Hamiltonian_Matrices

hamiltonian_matrices = Hamiltonian_Matrices(H0_uu=H0_uu,
                                            H0_dd=H0_dd,
                                            H0_ud=H0_ud,
                                            HU=HU,
                                            Jz=None,
                                            Jp=None,
                                            Jc=None)

rdm1_uu = np.load('rdm1_uu.npy')
rdm1_dd = np.load('rdm1_dd.npy')
rdm1_ud = np.load('rdm1_ud.npy')
rdm2 = np.load('rdm2.npy')

system = Transformable_QuantumChemistry(hamiltonian_matrices,
                                        rdm1=(rdm1_uu, rdm1_dd, rdm1_ud),
                                        rdm2=rdm2,
                                        tolerance=1e-1)

The input reduced density matrices have to adhere to specific conventions. The spin orbital optimization will otherwise yield non-physical results. We expect the following conventions for the reduced density matrices:

  • single 1-RDM:

  • spin-resolved 1-RDMs:

  • single 2-RDM:

Input exclusively for spin-like orbital optimization (SpinFinder)

If the spin-bath model derivation functionality is not required, the user can instantiate a SpinFinder object, which is the valid minimum input for the functions pertaining to the determination of the spin-like orbital basis.

from hqs_spin_mapper.transformables import SpinFinder

rdm1_uu = np.load('rdm1_uu.npy')
rdm1_dd = np.load('rdm1_dd.npy')
rdm1_ud = np.load('rdm1_ud.npy')
rdm2 = np.load('rdm2.npy')

system = SpinFinder(rdm1=(rdm1_uu, rdm1_dd, rdm1_ud),
                    rdm2=rdm2)

Keyworded constructor arguments for Transformable objects

Fields described by the Supports_SW_Transformation protocol:

  • system_size: Optional[int] = None: Number of orbitals in the system. This is automatically derived from the model description. (available for Transformable_Matrices, Transformable_QuantumChemistry)
  • rotation_matrix: Optional[numpy.ndarray] = None: Basis transformation matrix between different orbital bases. The index denotes the orbital in the old basis and the orbital in the new basis. Default is the identity. (available for Transformable_LatticeBuilder, Transformable_Matrices, Transformable_QuantumChemistry)
  • prefactor_cutoff: float = 1e-6: Terms of the Hamiltonian with coupling constant of absolute value smaller than the chosen value are discarded in the derivation of the spin-bath model (Schrieffer-Wolff transformation). (available for Transformable_LatticeBuilder, Transformable_Struqture, Transformable_Matrices, Transformable_QuantumChemistry)
  • site_type: str = "spinful_fermions": Type of particles occupying the orbitals. (available for Transformable_Matrices, Transformable_QuantumChemistry)
  • spin_indices: Optional[Set[int]] = None: User choice of orbitals to be treated as spin-like. The spin-like orbital optimization will add indices to this set. (available for Transformable_LatticeBuilder, Transformable_Struqture, Transformable_Matrices, Transformable_QuantumChemistry)

Fields described by the Supports_Spin_Optimization protocol:

  • tolerance: float = 1e-1: Chosen discrepancy used in for a basis orbital to be considered spin-like, where denotes the local parity. (available for Transformable_QuantumChemistry, SpinFinder)
  • optimization_steps: int = 4: Number of optimization loops to be performed. Parities are typically converged after 4 loops. (available for Transformable_QuantumChemistry, SpinFinder)

Other mutable fields of Transformable objects

Fields described by the Supports_SW_Transformation protocol:

  • generator: Union[ExpressionSpinful, ExpressionSpinful_complex]: Result for the generator of the Schrieffer-Wolff transformation.
  • transformed_hamiltonian: Union[ExpressionSpinful, ExpressionSpinful_complex]: Result for the transformed fermionic Hamiltonian.

Fields described by the Supports_Spin_Optimization protocol:

  • immutable_indices: Set[int]: Indices of orbitals that are exempt from parity optimization.

Custom Transformable objects (hqs_spin_mapper.spin_mapper_protocols)

The user can create their own custom input format by implementing a class that satisfies the protocols Supports_Spin_Optimization and/or Supports_SW_Transformation.

Functionality: Spin-like orbital basis optimization (hqs_spin_mapper.orbital_optimization)

HQS Spin Mapper provides the functionality to determine the optimal single-particle orbital basis, in which a subset of the basis orbitals exhibits the maximum possible spin-like character. We quantify the spin-like character of basis orbitals on the basis of their local parity . For more details about the use of the local parity as a measure for spin-like character, please consult the Theory section. An orbital in the optimized basis is stored as spin-like if its local parity satisfies , where corresponds to the value stored in the tolerance field of the Transformable.

The spin_orbital_optimization function in the orbital_optimization module determines a basis transformation which would optimize the local parities and stores the corresponding (proposed) basis transformation in the rotation_matrix field of the Transformable. The respective indices of the spin-like orbitals in this optimized basis are stored in the spin_indices field of the Transformable. The spin_orbital_optimization function requires Transformable objects that satisfy the Supports_Spin_Optimization protocol.

from hqs_spin_mapper.orbital_optimization import spin_orbital_optimization

spin_orbital_optimization(system)

Keyworded arguments

  • space_reduction: bool = True: If flag is set, the optimization is restricted to the indices not contained in the immutable_indices field of the Transformable.
  • fast_mode: bool = False: If flag is set, the rotations of the 2-RDM are performed approximately to reduce computation time.
  • target: str = "minimum": Type of extremum to be searched for. Options are: minimum, maximum, extremum.
  • _delta_sufficiently_spin_like: float = 1e-3: Exclude orbitals with initial parities and from the optimization procedure.

The function spin_orbital_optimization rotates the rdm1 and rdm2 to the optimized basis (in order to save memory), but not the matrix descriptions of the Hamiltonians.

The results of the optimized system can subsequently be accessed via:

spin_indices = system.spin_indices
rotation_matrix = system.rotation_matrix

rdm1 = system.rdm1
rdm2 = system.rdm2

Preparation of input data for the Schrieffer-Wolff transformation (hqs_spin_mapper.preconditioning)

For the derivation of the spin-bath model, the spin_indices: Set[int] field of the Transformable instance first requires a finite amount of indices. If the spin orbital optimization feature was used, spin_indices of the Transformable have been set to feature the obtained spin orbitals automatically. Otherwise spin_indices needs to be set manually by the user as in:

spin_like_orbital_indices: Set[int] = {4, 9}

system.spin_indices = spin_like_orbital_indices

An intuitive choice for spin-like orbitals are the orbitals with particularly strong repulsive on-site interactions .

The necessary preconditioning of the input data is done automatically by calling the preconditioning function on the Transformable once.

from hqs_spin_mapper.preconditioning import preconditioning

preconditioning(system)

The following tasks are performed by the preconditioning function:

  • Generalization of the input interaction terms HU, Jz, Jp, Jc to a four index tensor description.

  • Rotation of the matrix and tensor descriptions to the optimized basis by rotation with the rotation_matrix .

  • Mean-field decoupling (if Supports_Spin_Optimization) / Removal of bath exclusive interaction terms (if apply_pre_diagonalization = True).

  • Diagonalization of the quadratic bath Hamiltonian and rotation of the remaining Hamiltonians to the diagonal basis (if apply_pre_diagonalization = True).

  • Placing the spin indices in the middle of the system at (if apply_cross_geometry = True), with defined as

x0 = (system.system_size // 2) - (len(system.spin_indices) // 2)

Keyworded arguments

  • apply_pre_diagonalization: bool = True: Perform simplification of bath exclusive interactions and diagonalization of the resulting quadratic bath Hamiltonian.
  • apply_cross_geometry: bool = False: Move the spin indices to the center of the system (speeds up subsequent DMRG calculations).

Functionality: Effective spin-bath model derivation (hqs_spin_mapper.schrieffer_wolff)

With the input data prepared, the effective spin-bath model derivation becomes a simple call of schrieffer_wolff on the Transformable object.

from hqs_spin_mapper.schrieffer_wolff import schrieffer_wolff

schrieffer_wolff(system)

The results of the derivation, namely the transformed Hamiltonian and the generator of the transformation in are stored in the Transformable object instances as:

transformed_hamiltonian = system.transformed_hamiltonian
generator = system.generator

Keyworded arguments

  • _max_length: int = 5: Size limit for the terms included in the generator.
  • _vector_space_cap: int = 1500000: Size limit of the vector space that the generator is determined from.
  • _max_krylov_space_dimension: int = 200: Size restriction of the Krylov space used in the SVD.
  • rdm2: Optional[numpy.ndarray] = None: 2RDM for use in decoupling (currently not activated).

Extract tensor descriptions for the spin-bath model (hqs_spin_mapper.extract_couplings)

Using the functions contained in the extract_couplings module, one can retrieve the matrix/tensor description of the spin-bath model Hamiltonian.

from hqs_spin_mapper.extract_couplings import extract_quadratic, extract_quartic

(H0_uu, H0_dd, H0_ud) = extract_quadratic(system)
HU = extract_quartic(system)

The matrices and tensor only contain the quadratic/quartic terms of the Hamiltonian. The terms containing more fermionic operators are not extracted!

Serialize the spin-bath model Hamiltonian with struqture (hqs_spin_mapper.struqture_interface)

The transformed Hamiltonian can be exported to an instance of struqture_py.mixed_systems.MixedHamiltonianSystem for use in other HQS software or for storage in a file.

from hqs_spin_mapper.struqture_interface import get_struqture_description

struqture_description = get_struqture_description(system)

Best practices for HQS Spin Mapper

Avoiding local extrema in the orbital optimization

Depending on the initial basis, it is possible that a single call to spin_orbital_optimization fails to determine the true parity optimized basis. The pairwise optimization algorithm can get stuck in local parity extrema of the parameter manifold. This behavior can be discouraged by using the following procedure:

  • Call the spin_orbital_optimization function on the input data with target="minimum".
  • Add the orbital indices that became spin-like to the immutable_indices of the Transformable instance.
  • Call the spin_orbital_optimization function a second time on the updated input data with either target="extremum" or target="maximum".
  • Call the spin_orbital_optimization function on the input data a final time with target="minimum".

Choosing a sensible value for prefactor_cutoff (Important to avoid memory overflow!)

To limit the system memory required for the Schrieffer-Wolff transformation, we recommend adjusting the value of prefactor_cutoff. The terms of the input Hamiltonian in the optimized basis with an absolute coupling constant smaller than prefactor_cutoff will be dropped from the Hamiltonian prior to the Schrieffer-Wolff transformation step. This is done because a large number of small coupling constant terms slows down the SVD step of the transformation significantly, while not significantly changing the result. In the evaluation of the commutators that yield the transformed Hamiltonian, a larger number of terms can cause memory usage to exceed the system memory and resulting in a fatal error. We recommend inspecting the size of the coupling constants of the input Hamiltonian and choosing a value of prefactor_cutoff larger than the coupling constants that can be considered insignificant. When in doubt, choose a larger value for prefactor_cutoff first and then rerun the calculation with a slightly smaller value of prefactor_cutoff to check whether the result is sensitive to the change in prefactor_cutoff. If the result remains sensitive, we recommend gradually decreasing prefactor_cutoff and observing the memory usage.

We recommend to run the scripts or notebooks using the HQS Spin Mapper in an environment with a set ulimit. The ulimit can be set in the terminal before running the python interpreter via:

ulimit -v {maximum_desired_memory_in_kilobytes}

This way the process gets safely terminated without potentially crashing the system.

Estimation of the required memory

We provide a method for a rough estimate of the required memory to perform the Schrieffer-Wolff transformation at the currently set prefactor_cutoff. The function can be called on the Transformable object using:

from hqs_spin_mapper.preconditioning import (
    estimate_required_memory
)

estimate_required_memory(data)

The function estimate_required_memory prints an estimate of the necessary memory in GB. If the estimated value exceeds the realistically available system memory, we recommended to raise the value of the field prefactor_cutoff and to check the estimated memory requirement again.

Post-processing of the transformed Hamiltonian

The result of the model Hamiltonian derivation is stored as an ExpressionSpinful object. These objects offer a set of methods that can be used for further manipulation of the transformed Hamiltonian. Here we list some useful methods of the ExpressionSpinful class.

  • normalize(eps = 0.0): Orders the operators in the terms of the Hamiltonian in a standardized fashion and removes terms with an absolute prefactor smaller than eps.
  • normal_order: Rewrites operators as creators and annihilators and normal orders them (creators to the left of annihilators).
  • sort_descending: Sorts the terms in the expression by the absolute value of the coupling constants in descending order.
  • sort_ascending: Sorts the terms in the expression by the absolute value of the coupling constants in ascending order.
  • find(fermion): Returns the prefactor of the term containing the input operator.
  • project_on_local_hilbert_subspace(x, condition): Projects the Hamiltonian on a subspace of the Hilbert space where the condition is satisfied for site x. condition is the integer corresponding to a bit string, where the 1-bit indicates the x-local state, the 2-bit indicates , the 4-bit indicates , and the 8-bit indicates . The argument condition=6=4+2 thus corresponds to , namely the singly occupied site x.
  • convert_to_spin_model(positions): Replace fermionic operators acting on the positions with spin operators assuming the identity .

A simple script to convert the transformed fermionic Hamiltonian to a true spin-bath Hamiltonian description is given by:

spin_bath_system = system.transformed_hamiltonian

for n in cast(List[int], system.spin_indices):
    spin_bath_system.project_on_local_hilbert_subspace(n, 6)

spin_bath_system = spin_bath_system.convert_to_spin_model(system.spin_indices)