FlatLinearOperator

  • full name: tenpy.linalg.sparse.FlatLinearOperator

  • parent module: tenpy.linalg.sparse

  • type: class

Inheritance Diagram

Inheritance diagram of tenpy.linalg.sparse.FlatLinearOperator

Methods

FlatLinearOperator.__init__(npc_matvec, leg, …)

Initialize this LinearOperator.

FlatLinearOperator.adjoint()

Hermitian adjoint.

FlatLinearOperator.dot(x)

Matrix-matrix or matrix-vector multiplication.

FlatLinearOperator.eigenvectors([num_ev, …])

Find (dominant) eigenvector(s) of self using scipy.sparse.linalg.eigs().

FlatLinearOperator.flat_to_npc(vec)

Convert flat numpy vector of selected charge sector into npc Array.

FlatLinearOperator.flat_to_npc_None_sector(vec)

Convert flat vector of undetermined charge sectors into npc Array.

FlatLinearOperator.flat_to_npc_all_sectors(vec)

Convert flat vector of all charge sectors into npc Array with extra “charge” leg.

FlatLinearOperator.from_NpcArray(mat[, …])

Create a FlatLinearOperator from a square Array.

FlatLinearOperator.from_guess_with_pipe(…)

Create a FlatLinearOperator` from a matvec function acting on multiple legs.

FlatLinearOperator.matmat(X)

Matrix-matrix multiplication.

FlatLinearOperator.matvec(x)

Matrix-vector multiplication.

FlatLinearOperator.npc_to_flat(npc_vec)

Convert npc Array into a 1D ndarray, inverse of flat_to_npc().

FlatLinearOperator.npc_to_flat_all_sectors(npc_vec)

Convert npc Array with qtotal = self.charge_sector into ndarray.

FlatLinearOperator.rmatmat(X)

Adjoint matrix-matrix multiplication.

FlatLinearOperator.rmatvec(x)

Adjoint matrix-vector multiplication.

FlatLinearOperator.transpose()

Transpose this linear operator.

Class Attributes and Properties

FlatLinearOperator.H

Hermitian adjoint.

FlatLinearOperator.T

Transpose this linear operator.

FlatLinearOperator.charge_sector

Charge sector of the vector which is acted on.

FlatLinearOperator.ndim

class tenpy.linalg.sparse.FlatLinearOperator(*args, **kwargs)[source]

Bases: scipy.sparse.linalg.interface.LinearOperator

Square Linear operator acting on numpy arrays based on a matvec acting on npc Arrays.

Note that this class represents a square linear operator. In terms of charges, this means it has legs [self.leg.conj(), self.leg] and trivial (zero) qtotal.

Parameters
  • npc_matvec (function) – Function to calculate the action of the linear operator on an npc vector (with the specified leg). Has to return an npc vector with the same leg.

  • leg (LegCharge) – Leg of the vector on which npc_matvec can act on.

  • dtype (np.dtype) – The data type of the arrays.

  • charge_sector (None | charges | 0) – Selects the charge sector of the vector onto which the Linear operator acts. None stands for all sectors, 0 stands for the zero-charge sector. Defaults to 0, i.e., assumes the dominant eigenvector is in charge sector 0.

  • vec_label (None | str) – Label to be set to the npc vector before acting on it with npc_matvec. Ignored if None.

possible_charge_sectors

Each row corresponds to one possible choice for charge_sector.

Type

ndarray[QTYPE, ndim=2]

shape

The dimensions for the selected charge sector.

Type

(int, int)

dtype

The data type of the arrays.

Type

np.dtype

leg

Leg of the vector on which npc_matvec can act on.

Type

LegCharge

vec_label

Label to be set to the npc vector before acting on it with npc_matvec. Ignored if None.

Type

None | str

npc_matvec

Function to calculate the action of the linear operator on an npc vector (with one leg).

Type

function

matvec_count

The number of times npc_matvec was called.

Type

int

_mask

The indices of leg corresponding to the charge_sector to be diagonalized.

Type

ndarray[ndim=1, bool]

_npc_matvec_multileg

Only set if initalized with from_guess_with_pipe(). The npc_matvec function to be wrapped around. Takes the npc Array in multidimensional form and returns it that way.

Type

function | None

_labels_split

Only set if initalized with from_guess_with_pipe(). Labels of the guess before combining them into a pipe (stored as leg).

Type

list of str

classmethod from_NpcArray(mat, charge_sector=0)[source]

Create a FlatLinearOperator from a square Array.

Parameters
  • mat (Array) – A square matrix, with contractable legs.

  • charge_sector (None | charges | 0) – Selects the charge sector of the vector onto which the Linear operator acts. None stands for all sectors, 0 stands for the zero-charge sector. Defaults to 0, i.e., assumes the dominant eigenvector is in charge sector 0.

classmethod from_guess_with_pipe(npc_matvec, v0_guess, labels_split=None, dtype=None)[source]

Create a FlatLinearOperator` from a matvec function acting on multiple legs.

This function creates a wrapper matvec function to allow acting on a “vector” with multiple legs. The wrapper combines the legs into a LegPipe before calling the actual matvec function, and splits them again in the end.

Parameters
  • npc_matvec (function) – Function to calculate the action of the linear operator on an npc vector with the given split labels labels_split. Has to return an npc vector with the same legs.

  • v0_guess (Array) – Initial guess/starting vector which can be applied to npc_matvec.

  • labels_split (None | list of str) – Labels of v0_guess in the order in which they are to be combined into a LegPipe. None defaults to v0_guess.get_leg_labels().

  • dtype (np.dtype | None) – The data type of the arrays. None defaults to dtype of v0_guess (!).

Returns

  • lin_op (cls) – Instance of the class to be used as linear operator

  • guess_flat (np.ndarray) – Numpy vector representing the guess v0_guess.

property charge_sector

Charge sector of the vector which is acted on.

flat_to_npc(vec)[source]

Convert flat numpy vector of selected charge sector into npc Array.

If charge_sector is not None, convert to a 1D npc vector with leg self.leg. Otherwise convert vec, which can be non-zero in all charge sectors, to a npc matrix with an additional 'charge' leg to allow representing the full vector at once.

Parameters

vec (1D ndarray) – Numpy vector to be converted. Should have the entries according to self.charge_sector.

Returns

npc_vec – Same as vec, but converted into a npc array.

Return type

Array

npc_to_flat(npc_vec)[source]

Convert npc Array into a 1D ndarray, inverse of flat_to_npc().

Parameters

npc_vec (Array) – Npc Array to be converted. If self.charge_sector is not None, this should be a 1D array with that qtotal. If self.charge_sector is not None, it should have an additional "charge" leg, (as returned by flat_to_npc()).

Returns

vec – Same entries as npc_vec, but converted into a flat Numpy array.

Return type

1D ndarray

flat_to_npc_all_sectors(vec)[source]

Convert flat vector of all charge sectors into npc Array with extra “charge” leg.

Deprecated since version 0.7.3: This is merged into flat_to_npc() with self.charge_sector = None.

Parameters

vec (1D ndarray) – Numpy vector to be converted.

Returns

npc_vec – Same as vec, but converted into a npc array.

Return type

Array

flat_to_npc_None_sector(vec, cutoff=1e-10)[source]

Convert flat vector of undetermined charge sectors into npc Array.

The charge sector to be used is chosen as the block with the maximal norm, not by self.charge_sector (which might be None).

Parameters

vec (1D ndarray) – Numpy vector to be converted.

Returns

npc_vec – Same as vec, but converted into a npc array.

Return type

Array

npc_to_flat_all_sectors(npc_vec)[source]

Convert npc Array with qtotal = self.charge_sector into ndarray.

Deprecated since version 0.7.3: This is merged into npc_to_flat() with self.charge_sector = None.

Parameters

npc_vec (Array) – Npc Array to be converted. Should only have entries in self.charge_sector.

Returns

vec – Same as npc_vec, but converted into a flat Numpy array.

Return type

1D ndarray

eigenvectors(num_ev=1, max_num_ev=None, max_tol=1e-12, which='LM', v0=None, v0_npc=None, cutoff=1e-12, hermitian=False, **kwargs)[source]

Find (dominant) eigenvector(s) of self using scipy.sparse.linalg.eigs().

If no charge_sector was selected, we look in all charge sectors.

Parameters
  • num_ev (int) – Number of eigenvalues/vectors to look for.

  • max_num_ev (int) – scipy.sparse.linalg.speigs() somtimes raises a NoConvergenceError for small num_ev, which might be avoided by increasing num_ev. As a work-around, we try it again in the case of an error, just with larger num_ev up to max_num_ev. None defaults to num_ev + 2.

  • max_tol (float) – After the first NoConvergenceError we increase the tol argument to that value.

  • which (str) – Which eigenvalues to look for, see scipy.sparse.linalg.eigs(). More details also in argsort().

  • v0 (Array) – Initial guess as a “flat” numpy array.

  • v0_npc (Array) – Initial guess, to be converted by npc_to_flat().

  • cutoff (float) – Only used if self.charge_sector is None; in that case it determines when entries in a given charge-block are considered nonzero, and what counts as degenerate.

  • hermitian (bool) – If False (default), use scipy.sparse.linalg.eigs() If True, assume that self is hermitian and use scipy.sparse.linalg.eigsh().

  • **kwargs – Further keyword arguments given to scipy.sparse.linalg.eigsh() or scipy.sparse.linalg.eigs(), respectively.

Returns

  • eta (1D ndarray) – The eigenvalues, sorted according to which.

  • w (list of Array) – The eigenvectors corresponding to eta, as npc.Array with LegPipe.

property H

Hermitian adjoint.

Returns the Hermitian adjoint of self, aka the Hermitian conjugate or Hermitian transpose. For a complex matrix, the Hermitian adjoint is equal to the conjugate transpose.

Can be abbreviated self.H instead of self.adjoint().

Returns

A_H – Hermitian adjoint of self.

Return type

LinearOperator

property T

Transpose this linear operator.

Returns a LinearOperator that represents the transpose of this one. Can be abbreviated self.T instead of self.transpose().

adjoint()

Hermitian adjoint.

Returns the Hermitian adjoint of self, aka the Hermitian conjugate or Hermitian transpose. For a complex matrix, the Hermitian adjoint is equal to the conjugate transpose.

Can be abbreviated self.H instead of self.adjoint().

Returns

A_H – Hermitian adjoint of self.

Return type

LinearOperator

dot(x)

Matrix-matrix or matrix-vector multiplication.

Parameters

x (array_like) – 1-d or 2-d array, representing a vector or matrix.

Returns

Ax – 1-d or 2-d array (depending on the shape of x) that represents the result of applying this linear operator on x.

Return type

array

matmat(X)

Matrix-matrix multiplication.

Performs the operation y=A*X where A is an MxN linear operator and X dense N*K matrix or ndarray.

Parameters

X ({matrix, ndarray}) – An array with shape (N,K).

Returns

Y – A matrix or ndarray with shape (M,K) depending on the type of the X argument.

Return type

{matrix, ndarray}

Notes

This matmat wraps any user-specified matmat routine or overridden _matmat method to ensure that y has the correct type.

matvec(x)

Matrix-vector multiplication.

Performs the operation y=A*x where A is an MxN linear operator and x is a column vector or 1-d array.

Parameters

x ({matrix, ndarray}) – An array with shape (N,) or (N,1).

Returns

y – A matrix or ndarray with shape (M,) or (M,1) depending on the type and shape of the x argument.

Return type

{matrix, ndarray}

Notes

This matvec wraps the user-specified matvec routine or overridden _matvec method to ensure that y has the correct shape and type.

rmatmat(X)

Adjoint matrix-matrix multiplication.

Performs the operation y = A^H * x where A is an MxN linear operator and x is a column vector or 1-d array, or 2-d array. The default implementation defers to the adjoint.

Parameters

X ({matrix, ndarray}) – A matrix or 2D array.

Returns

Y – A matrix or 2D array depending on the type of the input.

Return type

{matrix, ndarray}

Notes

This rmatmat wraps the user-specified rmatmat routine.

rmatvec(x)

Adjoint matrix-vector multiplication.

Performs the operation y = A^H * x where A is an MxN linear operator and x is a column vector or 1-d array.

Parameters

x ({matrix, ndarray}) – An array with shape (M,) or (M,1).

Returns

y – A matrix or ndarray with shape (N,) or (N,1) depending on the type and shape of the x argument.

Return type

{matrix, ndarray}

Notes

This rmatvec wraps the user-specified rmatvec routine or overridden _rmatvec method to ensure that y has the correct shape and type.

transpose()

Transpose this linear operator.

Returns a LinearOperator that represents the transpose of this one. Can be abbreviated self.T instead of self.transpose().