pfaffian

Tools for converting Pfaffian wave functions into matrix product states (MPS).

Representing Nambu correlation matrices

The input of the algorithm is the Nambu correlation matrix that contains \(\langle cc \rangle\) and \(\langle c^\dagger c \rangle\) correlators in addition to the standard \(\langle c^\dagger c \rangle\) ones. Such a correlation matrix may be specified in two different bases:

  • In the standard complex-fermion basis, the Nambu correlation matrix is listed as 2×2 blocks C[2*i : 2*i+2, 2*j : 2*j+2] of

    \[\begin{split}\begin{pmatrix} \langle c_j^\dagger c_i \rangle & \langle c_j c_i \rangle \\ \langle c_j^\dagger c_i^\dagger \rangle & \langle c_j c_i^\dagger \rangle \end{pmatrix},\end{split}\]

    where \(c_i\) and \(c_i^\dagger\) are the standard fermion annihilation and creation operators on site \(0\le i < L.\)

  • We can also use the Majorana basis, defined by

    \[\begin{split}\gamma_{2n} &= (c^\dagger_n + c_n) / \sqrt{2},\\ \gamma_{2n+1} &= i(c^\dagger_n - c_n) / \sqrt{2}\end{split}\]

    for all \(0\le n < L\). In this basis, the correlation matrix is specified as the matrix \(\langle \gamma_j\gamma_i \rangle.\)

Classes

SchmidtModes

Boguliubov excitations that generate the Schmidt vectors of a Nambu mean-field state.

SchmidtVectors

Schmidt vectors of a Nambu mean-field state.

MPSTensorData

Data for computing one MPS tensor of a Pfaffian state.

High-level functions

correlation_matrix

Ground-state correlation matrix of a mean-field Nambu Hamiltonian.

parity

Fermion parity of a Boguliubov vacuum.

C_to_MPS

MPS representation of a Nambu mean-field ground state from its correlation matrix.

C_to_iMPS

iMPS representation of a Nambu mean-field state from correlation matrices.

Helper functions

vector_C2M

Converts mode vectors from complex-fermion to Majorana basis.

vector_M2C

Converts mode vectors from Majorana to complex-fermion basis.

matrix_C2M

Converts Hamiltonian or correlation matrices from complex-fermion to Majorana basis.

matrix_M2C

Converts Hamiltonian or correlation matrices from Majorana to complex-fermion basis.

assert_nambu

Indicates if a matrix is not Nambu symmetric.

assert_nambu_hamiltonian

Indicates if a Hamiltonian is not Nambu symmetric.

assert_nambu_correlation

Indicates if a correlation matrix is not Nambu symmetric.

Static variables

fermion_site

Lattice site prototype for the parity-conserving fermion MPS.

fermion_leg

LegCharge for the single-site Hilbert space of the parity-conserving fermion MPS.

chinfo

ChargeInfo for fermion parity conservation.

Classes

class temfpy.pfaffian.SchmidtModes(nL, nR, e, vL, vR, pL, pR)[source]

Boguliubov excitations that generate the Schmidt vectors of a Nambu mean-field state.

nL: int

Size of the left half of the system.

nR: int

Size of the right half of the system.

e: ndarray

array (n_entangled,) – Entangled eigenvalues \(0 < \lambda\le 1/2\) of the diagonal blocks of the correlation matrix.

vL: ndarray | None

array (2nL, 2nL) – Eigenvectors of the left-left block of the correlation matrix in the complex-fermion basis, if computed.

The eigenvectors are the columns of the matrix in the order of eigenvalues

  • \(0\le \lambda\le 1/2\) in increasing order

  • \(1\ge \lambda\ge 1/2\) in decreasing order,

i.e. the entangled modes are vL[:, nL-n_entangled : nL] (eigenvalues: e) and vL[:, -n_entangled:] (eigenvalues: 1-e).

The eigenvectors obey the Nambu symmetry:

vL[::2, nL:] == vL[1::2, :nL].conj() vL[1::2, nL:] == vL[::2, :nL].conj()

vR: ndarray | None

array (2nR, 2nR) – Eigenvectors of the right-right block of the correlation matrix in the complex-fermion basis, if computed.

The eigenvectors are the columns of the matrix in the order of eigenvalues

  • \(1/2\ge \lambda\ge 0\) in decreasing order

  • \(1/2\le \lambda\le 1\) in increasing order,

i.e. the entangled modes are vR[:, :n_entangled] (eigenvalues: e[::-1]) and vR[:, nR : nR+n_entangled] (eigenvalues: 1-e[::-1]).

They obey the Nambu symmetry:

vR[::2, nR:] == vR[1::2, :nR].conj()
vR[1::2, nR:] == vR[::2, :nR].conj()

If both vL and vR are known, their entangled components the left-right block of the Hamiltonian:

S = sqrt(e * (1-e))
S = concatenate((S, -S))
vL_entangled.T.conj() @ C[:2*nL, 2*nL:] @ vR_entangled[:, ::-1] == diag(S)

However, to better handle fermion anticommutation, the sign of vR is reversed if the left Boguliubov vacuum is parity odd.

pL: int | None

0 or 1 – Parity of the Boguliubov vacuum that annihilates the operators defined by eigenvectors vL[:, nL:], if computed.

pR: int | None

0 or 1 – Parity of the Boguliubov vacuum that annihilates the operators defined by eigenvectors vR[:, nR:], if computed.

parity(which='T')[source]

Parity of the Boguliubov vacuum on the specified half of the system.

Parameters:

which (str) – Either “L” for left or “R” for right side or “T” (default) for the overall parity of the state.

Return type:

int | None

Returns:

pL, pR, or the overall parity, depending on which or None if the requested parity is unknown.

property n_entangled: int

Number of entangled orbitals.

size(which='T')[source]

Size of the specified half or the whole of the system.

Parameters:

which (str) – Either “L” for left or “R” for right side or “T” (default) for total size.

Return type:

int

Returns:

nL, nR, or their sum, depending on which.

property vL_entangled: ndarray | None

Entangled left Schmidt mode orbitals, if computed.

property vR_entangled: ndarray | None

Entangled right Schmidt mode orbitals, if computed.

mode_vectors(which, entangled=False)[source]

Returns the Schmidt mode orbitals on the specified side.

Parameters:
  • which (str) – Either “L” for left or “R” for right side.

  • entangled (bool) – Whether to return the entangled (True) or all (False, default) eigenvectors.

Return type:

ndarray

Returns:

Either vL or vR, depending on which, truncated to entangled modes if entangled.

eigenvalues(which, entangled=False)[source]

Returns the Schmidt mode eigenvalues on the specified side.

Parameters:
  • which (str) – Either “L” for left or “R” for right side.

  • entangled (bool) – Whether to return the entangled (True) or all (False, default) eigenvalues.

Return type:

ndarray | None

Returns:

The eigenvalues corresponding to mode_vectors() with the same parameters.

property singular_values: ndarray | None

Singular values of the offdiagonal correlation matrix blocks.

If vL_entangled and vR_entangled are both known, satisfies

C_LR == vL_entangled @ diag(S) @ vR_entangled[:, ::-1].T.conj()

Otherwise, None.

classmethod from_correlation_matrix(C, x, trunc_par, *, basis, which='LR', diag_tol=1e-08, total_parity=None)[source]

Computes the SchmidtModes of a Nambu mean-field state with correlation matrix C for an entanglement cut between sites x-1 and x (zero-indexed).

We do this by diagonalising the left-left and right-right blocks of C. The eigenvectors give the Schmidt modes, the eigenvalues their relative weight in the entangled mode. We only treat modes with eigenvalues away from 0 or 1 by at least cutoff as entangled.

We enforce Nambu symmetry by redefining the eigenvectors for eigenvalues above 1/2 in terms of the other half.

If both left and right modes are computed, the entangled ones are paired up by ensuring that the eigenvalues sum to 1, and that they SVD-diagonalise the top right block C_LR and their signs are fixed so as to handle anticommutation signs optimally.

Parameters:
  • C (ndarray) – The correlation matrix.

  • x (int) – Position of the entanglement cut.

  • trunc_par (dict | StoppingCondition) –

    Which Schmidt modes should be kept as entangled.

    Must be either a StoppingCondition object or a dictionary with matching keys.

  • basis (str) – whether the correlation matrix is supplied in Majorana ("M", default) or complex-fermion ("C") basis.

  • which (str) – Whether to return Left and/or Right Schmidt modes. Must be a combination of "L", "R".

  • diag_tol (float) – If which=="LR", largest allowed offdiagonal matrix element in diagonalised / SVD correlation submatrices before an error is raised.

  • total_parity (int | None) – Used to infer the parity of the leading Schmidt state on both sides if Schmidt modes are only computed on one side.

Return type:

SchmidtModes

Note

  • If trunc_par.svd_min is not provided, a default of 1e-6 (i.e., a truncation threshold of 1e-12) is used.

  • If trunc_par.degeneracy_tol is not provided, the degeneracy tolerance defaults to 1e-12.

property e_ratio: ndarray

\(\log((1-\lambda)/\lambda\) for all eigenvalues in e.

embed_subsets(sets)[source]

Given an array of \(\gamma^\dagger\) excitations listed in order of the left singular vectors, generates the set of excitations on both sides.

Parameters:

sets (bool ndarray (n, n_entangled)) –

Array of occupation numbers.

Each row specifies one Schmidt state by listing whether each \(\gamma^\dagger\) excitation is “occupied”.

The order of the excitations is consistent with vL.

Return type:

tuple[ndarray | None, ndarray | None]

Returns:

schmidt_values(sets)[source]

Schmidt values of the Schmidt vectors with given excitation numbers.

Parameters:

sets (bool ndarray (n, n_entangled)) –

Array of occupation numbers.

Each row specifies one Schmidt state by listing whether each \(\gamma^\dagger\) excitation is “occupied”.

The order of the excitations is consistent with vL.

Returns:

λ – Schmidt values corresponding to the input Schmidt vectors.

Return type:

ndarray (n,)

class temfpy.pfaffian.SchmidtVectors(modes, left_sets, right_sets, schmidt_values, idx_n, idx_parity)[source]

Schmidt vectors of a Nambu mean-field state.

modes: SchmidtModes

The Boguliubov excitations underlying the Schmidt vectors.

left_sets: ndarray | None

bool (n_schmidt, n_entangled) – Left Schmidt vectors.

Each row contains the occupation of all entangled \(\gamma_L^\dagger\) excitations, described by vL[:, :n_entangled], in one left Schmidt vector, if modes.vL is not None.

right_sets: ndarray | None

bool (n_schmidt, n_entangled) – Right Schmidt vectors.

Each row contains the occupation of all entangled \(\gamma_R^\dagger\) excitations, described by vR[:, nR-n_entangled:nR], in one right Schmidt vector, if vR is not None.

If both left_sets and right_sets are computed, they are related by left_sets = right_sets[:, ::-1].

schmidt_values: ndarray

array (n_schmidt,): Schmidt values \(\lambda_\alpha\) corresponding to each Schmidt vector.

Collated by number of Boguliubov excitations (cf. idx_n); sorted in decreasing order within each excitation-number sector.

idx_n: dict[int, slice]

Dictionary mapping the number of \(\gamma^\dagger\) excitations to slice of sets/singular values corresponding to that excitation number.

That is, each row in left_sets[idx_n[n]] contains n True entries.

Slices follow each other in the order 0, 2,…, 1, 3,…

idx_parity: dict[int, slice]

Dictionary mapping the parity (0 or 1) of \(\gamma^\dagger\) excitations to slice of sets/singular values corresponding to that excitation number.

That is, each row in left_sets[idx_parity[n]] contains an even (odd) number of True entries if n=0 (1).

property n_schmidt: int

Number of Schmidt vectors.

property n_entangled: int

Number of entangled orbitals.

property nL: int

Size of the left half of the system.

property nR: int

Size of the right half of the system.

size(which='T')[source]

Size of the specified half or the whole of the system.

Parameters:

which (str) – Either “L” for left or “R” for right side or “T” (default) for total size.

Return type:

int

Returns:

the appropriate system size.

property vL: ndarray | None

Left Schmidt mode orbitals vL.

property vR: ndarray | None

Right Schmidt mode orbitals vR.

mode_vectors(which, entangled=False)[source]

Returns the Schmidt mode orbitals on the specified side.

Parameters:
  • which (str) – Either “L” for left or “R” for right side.

  • entangled (bool) – Whether to return the entangled (True) or all (False, default) eigenvectors.

Return type:

ndarray

Returns:

Either vL or vR, depending on which, truncated to entangled modes if entangled.

property pL: int | None

Parity of left Boguliubov vacuum pL.

property pR: int | None

Parity of right Boguliubov vacuum pR.

parity(which='T')[source]

Parity of the Boguliubov vacuum on the specified half of the system.

Parameters:

which (str) – Either “L” for left or “R” for right side or “T” (default) for the overall parity of the state.

Return type:

int | None

Returns:

pL, pR, or the overall parity, depending on which.

sets(which)[source]

Returns the sets of Boguliubov excitations on the specified side.

Parameters:

which (str) – Either “L” for left or “R” for right side.

Returns:

Either left_sets or right_sets, depending on which.

classmethod from_schmidt_modes(modes, trunc_par)[source]

Generates the most significant SchmidtVectors from an instance of SchmidtModes.

Parameters:
  • modes (SchmidtModes) – The Schmidt modes.

  • trunc_par (dict | StoppingCondition) –

    Specifies which Schmidt states should be kept.

    Must be either a StoppingCondition object or a dictionary with matching keys.

    The filtering function is_sector is applied to the number of \(\gamma^\dagger\) excitations.

Return type:

SchmidtVectors

classmethod from_correlation_matrix(C, x, trunc_par, *, basis, which='LR', diag_tol=1e-08, total_parity=None)[source]

Computes the most significant SchmidtVectors of a Nambu mean-field state with correlation matrix C for an entanglement cut between sites x-1 and x (zero-indexed).

Calls SchmidtModes.from_correlation_matrix() (see there for details of the parameters) and from_schmidt_modes().

Return type:

SchmidtVectors

class temfpy.pfaffian.MPSTensorData(mode, norm, pfaffian_matrix, labels, qtotal, leg_bra, new_sets_bra, idx_n_bra, leg_idx_bra, leg_ket, new_sets_ket, idx_n_ket)[source]

Data for computing one MPS tensor of a Pfaffian state.

  • If mode is "left", contains an implicit description of the left canonical tensor

    \[A^{n_i}_{\alpha\beta} = (\langle n_i | \otimes_g \langle L^{(i-1)}_{\alpha}|) |L^{(i)}_{\beta} \rangle.\]
  • If mode is "right", contains an implicit description of the right canonical tensor

    \[B^{n_i}_{\beta\alpha} = (\langle R^{(i)}_\alpha| \otimes_g \langle n_i|) |R^{(i-1)}_\beta \rangle.\]
mode: str

Whether the overlap is between "left" or "right" Schmidt vectors.

norm: float

Overlap of the normalised Boguliubov vacua for the two entanglement cuts.

pfaffian_matrix: ndarray

Matrix of Pfaffian entries for all bogulon excitations that appear in at least one Schmidt state.

It is a block matrix, with the first rows/columns corresponding to excitations on the ket side, followed by excitations on the bra side.

labels: list[str, str]

Labels of the bra and ket leg(pipe)s.

qtotal: int

Total fermion parity of the tensor to ensure matching parities in Schmidt vectors.

  • 0 if mode == "left".

  • 0 if mode == "right" and virtual legs are labelled with fermion parity to the right.

  • Otherwise, difference of total fermion parity between the ket and bra Schmidt vectors.

leg_bra: LegCharge

Tensor leg corresponding to the bra Schmidt vectors.

If the tensor contains a physical leg, it is united to the bra leg as an unsorted tenpy.linalg.charges.LegPipe.

new_sets_bra: ndarray

Bra Schmidt vectors in terms of the excitations in pfaffian_matrix.

Given the block structure of pfaffian_matrix, each row starts with a number of False entries corresponding to the ket excitations.

If the tensor contains a physical leg, it is double the length of sets and contains all Schmidt vectors with the on-site degree of freedom once empty, one filled.

If needed for matching the bra and ket parities, the Boguliubov orbital with λ closest to 1/2 is particle-hole flipped.

After all this, the sets are resorted by parity and number of Boguliubov excitations.

idx_n_bra: ndarray

Analogue of idx_n for new_sets_bra.

leg_idx_bra: ndarray

Mapping between leg_bra and new_sets_bra.

Namely, leg_idx_bra[i] gives the index in leg_bra corresponding to new_sets_bra[i].

leg_ket: LegCharge

Tensor leg corresponding to the ket Schmidt vectors.

new_sets_ket: ndarray

Ket Schmidt vectors in terms of the excitations in pfaffian_matrix.

Given the block structure of pfaffian_matrix, each row ends with a number of False entries corresponding to the bra excitations.

idx_n_ket: ndarray

idx_n of the ket Schmidt vector.

classmethod from_schmidt_vectors(Schmidt_bra, Schmidt_ket, mode, *, nambu_tolerance=1e-08, min_SV=1e-06)[source]

Builds an MPSTensorData object from the Schmidt vectors of the two surrounding entanglement cuts.

Parameters:
  • Schmidt_bra (SchmidtVectors) –

    Schmidt vectors on the “bra” side.

    That is, for the entanglement cut to the left if mode=="left" or to the right if mode="right".

    Must contain a description of left (right) Schmidt vectors if mode=="left" ("right").

  • Schmidt_ket (SchmidtVectors) –

    Schmidt vectors on the “ket” side.

    That is, for the entanglement cut to the right if mode=="left" or to the left if mode="right".

    Must contain a description of left (right) Schmidt vectors if mode=="left" ("right").

  • mode (str) – Specifies which Schmidt vectors are used to build the MPS tensor and the canonical form the tensor is computed in.

  • nambu_tolerance (float) – numerical tolerance for checking Nambu symmetry

  • min_SV (float) – smallest acceptable singular value of U before it is considered singular

Return type:

MPSTensorData

to_npc_array()[source]

Computes the MPS tensor as a TeNPy Array object.

If mode is “left”, returns a left canonical tensor; if “right”, a right canonical one.

Return type:

Array

High-level functions

temfpy.pfaffian.correlation_matrix(H, basis=None, *, rtol=0, atol=1e-10)[source]

Ground-state correlation matrix of a mean-field Nambu Hamiltonian.

Parameters:
  • H (ndarray) –

    Real-space Nambu Hamiltonian.

    If basis is given, it must be in the basis indicated by its first character.

  • basis (str | None) –

    Basis used by the input and the output.

    If specified, it is a string of the form "X->Y", where the characters X and Y must be M or C:

    • X controls whether the Hamiltonian is given in the Majorana or the complex-fermion basis.

    • Y controls whether the correlation matrix is returned in the Majorana or the complex-fermion basis.

    If not specified, the output is returned in the same basis as the input, without testing for Nambu symmetry.

  • rtol (float) – Relative tolerance for eigenvalue and matrix-element checks.

  • atol (float) – Absolute tolerance for eigenvalue and matrix-element checks.

Return type:

ndarray

Returns:

The correlation matrix C in the format indicated by the last character of mode.

Notes

In Majorana mode "M->*", the Hamiltonian is expected as the array of coefficients of \(\gamma_i \gamma_j\).

In complex-fermion mode "C->*", the Hamiltonian is expected as 2×2 blocks of coefficients of

\[\begin{split}\begin{pmatrix} c_i^\dagger c_j & c_i^\dagger c_j^\dagger \\ c_i c_j & c_i c_j^\dagger \end{pmatrix}.\end{split}\]

For the format of the correlation matrices, see Representing Nambu correlation matrices.

temfpy.pfaffian.parity(V, *, tol=1e-12)[source]

Fermion parity of a Boguliubov vacuum.

The calculation is based on the Bloch-Messiah decomposition. If the Boguliubov operators \(\gamma\) that define the vacuum are related to the standard fermion operators \(c\) by

\[\begin{split}(\gamma^\dagger, \gamma) = (c^\dagger, c) \begin{pmatrix} U & V^* \\ V & U^* \end{pmatrix},\end{split}\]

the singular values of the matrix \(V\) (and also \(V^*\)) are (in decreasing order) \(1, \dots, 1, \sigma_1, \sigma_1, \dots, \sigma_n, \sigma_n, 0, \dots, 0\). The parity of the vacuum is that of the number of completely filled modes (i.e., singular values 1).

Parameters:
  • V (ndarray) – submatrix of Nambu mode unitary that maps \(c\) to \(\gamma^\dagger\) (or \(c^\dagger\) to \(\gamma\)).

  • tol (float) – tolerance for considering two numbers equal (only needed for edge cases of 1 and 2 modes)

Return type:

int

Returns:

parity as 0 (even) or 1 (odd)

temfpy.pfaffian.C_to_MPS(C, trunc_par, *, basis, diag_tol=1e-08, ortho_center=None)[source]

MPS representation of a Nambu mean-field ground state from its correlation matrix.

Parameters:
  • C (ndarray) – Nambu correlation matrix in the basis indicated by basis.

  • trunc_par (dict | StoppingCondition) –

    Specifies which Schmidt states should be kept.

    Must be either a StoppingCondition object or a dictionary with matching keys.

    Only specify the field sectors if you know what you are doing!!

  • basis (str) – “M” or “C”, indicates whether the correlation matrix is given in the Majorana or the complex-fermion basis.

  • ortho_center (int) – Orthogonality centre of the mixed canonical MPS. Midpoint of the chain by default.

  • diag_tol (float) – Largest allowed offdiagonal matrix element in diagonalised / SVD correlation submatrices before an error is raised.

Return type:

MPS

Returns:

The wave function as a TeNPy MPS object.

Note

  • If trunc_par.svd_min is not provided, the truncation threshold defaults to 1e-6.

  • If trunc_par.degeneracy_tol is not provided, the degeneracy tolerance defaults to 1e-12.

temfpy.pfaffian.C_to_iMPS(C_short, C_long, trunc_par, sites_per_cell, cut, *, basis, diag_tol=1e-08, unitary_tol=1e-06, schmidt_tol=1e-06)[source]

iMPS representation of a Nambu mean-field state from correlation matrices.

The two correlation matrices are expected to represent the ground states of a gapped, translation invariant Hamiltonian on two system sizes that differ by one repeating unit cell.

The method is analogous to iMPS.MPS_to_iMPS(), with two differences:

  • No explicit MPS tensors are computed for the environment of the iMPS unit cell. Instead, the Schmidt vector overlaps needed for gauge fixing are computed using the Pfaffian state overlap formulas implemented in MPSTensorData.

  • The rightmost tensor is computed directly using the right Schmidt vectors of the shorter chain. This means that no separate iMPSErrors are returned for the right side.

Parameters:
  • C_short (ndarray) – Nambu correlation matrix in the basis indicated by basis for the shorter chain.

  • C_long (ndarray) – Nambu correlation matrix in the basis indicated by basis for the longer chain.

  • trunc_par (dict | StoppingCondition) –

    Specifies which Schmidt states should be kept.

    Must be either a StoppingCondition object or a dictionary with matching keys.

    Only specify the field sectors if you know what you are doing!!

  • sites_per_cell (int) – Size of the iMPS unit cell.

  • cut (int) – First site of the repeating unit cell in C_long.

  • basis (str) – “M” or “C”, indicates whether the correlation matrix is given in the Majorana or the complex-fermion basis.

  • diag_tol (float) – Largest allowed offdiagonal matrix element in diagonalised / SVD correlation submatrices before an error is raised.

  • unitary_tol (float) – Maximum deviation of the gauge rotation matrices from unitarity before a warning is raised.

  • schmidt_tol (float) – Maximum mixing of unequal Schmidt values by the gauge rotation matrices before a warning is raised.

Return type:

tuple[MPS, iMPSError]

Returns:

  • iMPS (MPS) – iMPS with unit cell size sites_per_cell, constructed from the additional unit cell of mps_long.

  • validation_metric (iMPSError) – Errors introduced during the conversion.

Note

  • If trunc_par.svd_min is not provided, the truncation threshold defaults to 1e-6.

  • If trunc_par.degeneracy_tol is not provided, the degeneracy tolerance defaults to 1e-12.

Helper functions

temfpy.pfaffian.vector_C2M(v)[source]

Converts mode vectors from complex-fermion to Majorana basis.

Parameters:

v (ndarray) –

A vector or a multidimensional array of vectors in the complex-fermion basis.

In the latter case, the “site” dimension is expected as first (i.e. several vectors should be supplied as columns of a matrix).

The complex-fermion basis should be listed in the order \(c^\dagger_1, c_1, \dots, c^\dagger_L, c_L\).

Return type:

ndarray

Returns:

The given vectors in the Majorana basis.

temfpy.pfaffian.vector_M2C(v)[source]

Converts mode vectors from Majorana to complex-fermion basis.

Parameters:

v (ndarray) –

A vector or a multidimensional array of vectors in the Majorana basis.

In the latter case, the “site” dimension is expected as first (i.e. several vectors should be supplied as columns of a matrix).

Return type:

ndarray

Returns:

The given vectors in the complex-fermion basis.

The complex-fermion basis is listed in the order \(c^\dagger_1, c_1, \dots, c^\dagger_L, c_L\).

temfpy.pfaffian.matrix_C2M(H)[source]

Converts Hamiltonian or correlation matrices from complex-fermion to Majorana basis.

Parameters:

H (ndarray) –

A matrix in the complex-fermion basis.

See the module documentation for the expected layout.

Return type:

ndarray

Returns:

The given matrix in the Majorana basis.

temfpy.pfaffian.matrix_M2C(H)[source]

Converts Hamiltonian or correlation matrices from Majorana to complex-fermion basis.

Parameters:

H (ndarray) – A matrix in the Majorana basis.

Return type:

ndarray

Returns:

The given matrix in the complex-fermion basis.

See the module documentation for the detailed layout.

temfpy.pfaffian.assert_nambu(C, basis=None, offset=None, name='', rtol=0, atol=1e-10)[source]

Indicates if a matrix is not Nambu symmetric.

In Majorana basis, Nambu symmetry requires that the matrix be imaginary and antisymmetric, except for a constant offset / 2 on the diagonal.

In complex-fermion basis, the matrix can be split into blocks by creation and annihilation operators, e.g.,

\[\begin{split}C = \begin{pmatrix} \langle c_j^\dagger c_i \rangle & \langle c_j c_i \rangle \\ \langle c_j^\dagger c_i^\dagger \rangle & \langle c_j c_i^\dagger \rangle \end{pmatrix} =: \begin{pmatrix} C_{11} & C_{12} \\ C_{21} & C_{22} \end{pmatrix},\end{split}\]

which must satisfy

\[\begin{split}C_{11} + C_{22}^* &= \texttt{offset}\times\mathbb{I} \\ C_{12} &= -C_{21}^*.\end{split}\]
Parameters:
  • C (ndarray) – Nambu matrix to be checked.

  • basis (str) –

    Whether the matrix is supplied in Majorana ("M") or complex-fermion ("C") basis.

    If unspecified, only checks if the matrix is Hermitian.

  • offset (float) –

    Constant diagonal offset.

    Should be 0 for Hamiltonians and 1 for correlation matrices.

  • name (str) – Type of matrix (e.g. "Hamiltonian")

  • rtol (float) – Relative tolerance for matrix-element checks.

  • atol (float) – Absolute tolerance for matrix-element checks.

Return type:

ndarray

Returns:

Regularised Nambu matrix.

  • Skew-Hermitian component is removed.

  • In the Majorana basis, small real-part deviations are deleted.

  • In the complex-fermion basis, small imaginary parts are pruned if the matrix is almost real.

Raises:

AssertionError | testing.ComparisonWarning – If C is not Hermitian or Nambu symmetric up to the given tolerance. Whether an exception or a warning is raised is determined by testing.TEST_ACTION.

temfpy.pfaffian.assert_nambu_hamiltonian(C: ndarray, basis: str = None, *, offset: float = 0, name: str = 'Hamiltonian', rtol: float = 0, atol: float = 1e-10) ndarray

Indicates if a Hamiltonian is not Nambu symmetric.

See assert_nambu() for details. (Arguments offset and name are fixed.)

temfpy.pfaffian.assert_nambu_correlation(C: ndarray, basis: str = None, *, offset: float = 1, name: str = 'correlation matrix', rtol: float = 0, atol: float = 1e-10) ndarray

Indicates if a correlation matrix is not Nambu symmetric.

See assert_nambu() for details. (Arguments offset and name are fixed.)