symmray.sparse.sparse_fermionic_array

Fermionic symmetric arrays with block sparse backend.

Classes

FermionicArray

A fermionic block symmetry array.

Z2FermionicArray

A fermionic block array with Z2 symmetry.

U1FermionicArray

A fermionic block array with U1 symmetry.

Z2Z2FermionicArray

A fermionic block array with Z2 x Z2 symmetry.

U1U1FermionicArray

A fermionic block array with U1 x U1 symmetry.

Functions

argsort(seq)

Module Contents

symmray.sparse.sparse_fermionic_array.argsort(seq)[source]
class symmray.sparse.sparse_fermionic_array.FermionicArray(indices, charge=None, blocks=(), phases=(), label=None, dummy_modes=None, symmetry=None)[source]

Bases: symmray.fermionic_common.FermionicCommon, symmray.sparse.sparse_array_common.SparseArrayCommon, symmray.sparse.sparse_data_common.BlockCommon, symmray.array_common.ArrayCommon, symmray.common.SymmrayCommon

A fermionic block symmetry array.

Parameters:
  • indices (tuple of Index) – The indices of the array.

  • charge (hashable, optionals) – The total charge of the array, if not given it will be inferred from either the first sector or set to the identity charge, if no sectors are given.

  • blocks (dict, optional) – The blocks of the array, by default empty.

  • phases (dict, optional) – The lazy phases of each block, by default empty.

  • label (hashable, optional) – An optional label for the array, potentially needed for ordering dummy odd fermionic modes.

  • symmetry (str or Symmetry, optional) – The symmetry of the array, if not using a specific symmetry class.

__slots__ = ('_blocks', '_charge', '_dummy_modes', '_indices', '_label', '_phases', '_symmetry')
fermionic = True
static_symmetry = None
_phases
_dummy_modes = ()
to_pytree()[source]

Convert this sparse fermionic array to a pytree purely of non-symmray containers and objects.

classmethod from_pytree(data)[source]

Create a sparse fermionic array from a pytree purely of non-symmray containers and objects.

property parity

The parity of the total charge.

property phases

The lazy phases of each sector. Trivial phases are not necessarily stored.

copy()[source]

Create a copy of this fermionic array.

Returns:

The copied array.

Return type:

FermionicArray

new_with(indices, charge, blocks, label=None)[source]

Create a new block sparse fermionic array of the same class as this one. Unlike copy, this does not copy over any existing data and drops by default label, phases, and dummy_modes.

copy_with(indices=None, blocks=None, charge=None, phases=None)[source]

Create a copy of this fermionic array with some attributes replaced.

Parameters:
  • indices (tuple of Index, optional) – The new indices, if None, the original indices are used.

  • blocks (dict, optional) – The new blocks, if None, the original blocks are used.

  • charge (int, optional) – The new total charge, if None, the original charge is used.

  • phases (dict, optional) – The new phases, if None, the original phases are used.

modify(indices=None, blocks=None, charge=None, phases=None, dummy_modes=None)[source]

Modify this fermionic array inplace. This is for internal use, and does not perform any checks on the updated attributes.

Parameters:
  • indices (tuple of Index, optional) – The new indices, if None, the original indices are used.

  • blocks (dict, optional) – The new blocks, if None, the original blocks are used.

  • charge (int, optional) – The new total charge, if None, the original charge is used.

  • phases (dict, optional) – The new phases, if None, the original phases are used.

  • dummy_modes (object or FermionicOperator, optional) – The new dummy_modes, if None, the original dummy_modes are used.

randomize_phases(seed=None, inplace=False) FermionicArray[source]

Randomize the phases of each sector to either +1 or -1. This is useful for testing.

Parameters:
  • seed (int or numpy.random.Generator, optional) – The random seed or generator, by default None.

  • inplace (bool, optional) – Whether to perform the operation in place.

Returns:

The phase randomized array.

Return type:

FermionicArray

phase_sync(inplace=False) FermionicArray[source]

Multiply all lazy phases into the block arrays.

Parameters:

inplace (bool, optional) – Whether to perform the operation in place.

Returns:

The resolved array, which now has no lazy phases.

Return type:

FermionicArray

phase_flip(*axs, inplace=False) FermionicArray[source]

Flip the phase of all sectors with odd parity at the given axis.

Parameters:
  • ax (int) – The axis along which to flip the phase.

  • inplace (bool, optional) – Whether to perform the operation in place.

Returns:

The phase-flipped array.

Return type:

FermionicArray

phase_transpose(axes=None, inplace=False)[source]

Phase this fermionic array as if it were transposed virtually, i.e. the actual arrays are not transposed. Useful when one wants the actual data layout to differ from the required fermionic mode layout.

Parameters:
  • axes (tuple of int, optional) – The permutation of axes, by default None.

  • inplace (bool, optional) – Whether to perform the operation in place.

Return type:

FermionicArray

phase_sector(sector, inplace=False)[source]

Flip the phase of a specific sector.

Parameters:
  • sector (tuple of hashable) – The sector to flip the phase of.

  • inplace (bool, optional) – Whether to perform the operation in place.

Return type:

FermionicArray

phase_global(inplace=False)[source]

Flip the global phase of the array.

Parameters:

inplace (bool, optional) – Whether to perform the operation in place.

Return type:

FermionicArray

_map_blocks(fn_block=None, fn_sector=None, fn_filter=None)[source]

Map the blocks and their keys (sectors) of the array inplace.

_resolve_dummy_modes_conj(phase_permutation=True)[source]

Assuming we have effectively taken the conjugate of a fermionic array with dummy modes, get their new order and compute any phase changes coming from moving back to the beginning of the index order.

_resolve_dummy_modes_combine(left, right)[source]

Calculate the new combined dummy odd modes and any associated global phases combing from contracting two fermionic arrays a and b. This modifies this array in place.

_resolve_dummy_modes_squeeze(axes_squeeze)[source]

Assuming we are about to squeeze away axes_squeeze, compute the phases associated with moving them to the beginning of the array, and then turn them into dummy modes.

transpose(axes=None, phase=True, inplace=False)[source]

Transpose the fermionic array, by default accounting for the phases accumulated from swapping odd charges.

Parameters:
  • axes (tuple of int, optional) – The permutation of axes, by default None.

  • phase (bool, optional) – Whether to flip the phase of sectors whose odd charges undergo a odd permutation. By default True.

  • inplace (bool, optional) – Whether to perform the operation in place.

Returns:

The transposed array.

Return type:

FermionicArray

conj(phase_permutation=True, phase_dual=False, inplace=False)[source]

Conjugate this fermionic array. By default this include phases from both the virtual flipping of all axes, but not the conjugation of dual indices, such that:

(
    tensordot_fermionic(x.conj(), x, ndim) ==
    tensordot_fermionic(x, x.conj(), ndim)
)

If all indices have matching dualness (i.e. all bra or all ket), or you set phase_dual=True then the above contractions will also be equal to x.norm() ** 2.

Parameters:
  • phase_permutation (bool, optional) – Whether to flip the phase of sectors whose odd charges undergo a odd permutation due to virtually flipping the order of axes, by default True.

  • phase_dual (bool, optional) – Whether to flip the phase of dual indices, by default False. If a FermionicArray has a mix of dual and non-dual indices, and you are explicitly forming the norm, you may want to set this to True. But if it is part of a large tensor network you only need to flip the phase of true ‘outer’ dual indices.

Return type:

FermionicArray

dagger(phase_dual=False, inplace=False)[source]

Fermionic adjoint, implements .H attribute.

Parameters:
  • phase_dual (bool, optional) – Whether to flip the phase of dual indices, by default False. If a FermionicArray has a mix of dual and non-dual indices, and you are explicitly forming the norm, you may want to set this to True. But if it is part of a large tensor network you only need to flip the phase of ‘outer’ dual indices.

  • inplace (bool, optional) – Whether to perform the operation in place.

Returns:

The conjugate transposed array.

Return type:

FermionicArray

eigh(phase_eigenvalues=True, **kwargs) tuple[symmray.sparse.sparse_vector.BlockVector, FermionicArray][source]

Eigenvalue decomposition of this assumed Hermitian block sparse fermionic array.

Parameters:

phase_eigenvalues (bool, optional) – If True, any local phase will be absorbed into the eigenvalues, such that U @ diag(w) @ U.H == a always holds. By default True. If False then one of U or U.H should be phased individually to account for local phasesin the above expression.

Returns:

  • w (BlockVector) – The eigenvalues as a vector. If phase_eigenvalues is True, then eigenvalues corresponding to odd parity sectors will be negated depending of the inner index dualness.

  • U (FermionicArray) – The array of eigenvectors.

class symmray.sparse.sparse_fermionic_array.Z2FermionicArray(indices, charge=None, blocks=(), phases=(), label=None, dummy_modes=None, symmetry=None)[source]

Bases: FermionicArray

A fermionic block array with Z2 symmetry.

static_symmetry
to_pyblock3(flat=False)[source]
class symmray.sparse.sparse_fermionic_array.U1FermionicArray(indices, charge=None, blocks=(), phases=(), label=None, dummy_modes=None, symmetry=None)[source]

Bases: FermionicArray

A fermionic block array with U1 symmetry.

static_symmetry
to_pyblock3(flat=False)[source]
class symmray.sparse.sparse_fermionic_array.Z2Z2FermionicArray(indices, charge=None, blocks=(), phases=(), label=None, dummy_modes=None, symmetry=None)[source]

Bases: FermionicArray

A fermionic block array with Z2 x Z2 symmetry.

static_symmetry
class symmray.sparse.sparse_fermionic_array.U1U1FermionicArray(indices, charge=None, blocks=(), phases=(), label=None, dummy_modes=None, symmetry=None)[source]

Bases: FermionicArray

A fermionic block array with U1 x U1 symmetry.

static_symmetry