tnpy.finite_dmrg.FiniteDMRG
2.1. tnpy.finite_dmrg.FiniteDMRG#
- class tnpy.finite_dmrg.FiniteDMRG(mpo, bond_dim, block_size=1, mps=None, exact_solver_dim=200)[source]#
Bases:
objectThe Density Matrix Renormalization Group (DMRG) algorithm of a finite size system.
- Parameters
mpo (tnpy.operators.MatrixProductOperator) – Matrix product operator.
bond_dim (int) – Bond dimensions.
block_size (int) – Block size of each local update.
mps (tnpy.matrix_product_state.MatrixProductState) – (Optional) Initial guess of
MatrixProductState.exact_solver_dim (int) – A switch point at which the exact eigensolver should be used below certain matrix size.
Examples
To start the algorithm, one can call
run(),>>> fdmrg = FiniteDMRG(mpo, bond_dim=20) >>> fdmrg.run(tol=1e-8)
The optimized
MatrixProductStatecan then be retrieved frommps.- __init__(mpo, bond_dim, block_size=1, mps=None, exact_solver_dim=200)[source]#
The Density Matrix Renormalization Group (DMRG) algorithm of a finite size system.
- Parameters
mpo (tnpy.operators.MatrixProductOperator) – Matrix product operator.
bond_dim (int) – Bond dimensions.
block_size (int) – Block size of each local update.
mps (Optional[tnpy.matrix_product_state.MatrixProductState]) – (Optional) Initial guess of
MatrixProductState.exact_solver_dim (int) – A switch point at which the exact eigensolver should be used below certain matrix size.
Examples
To start the algorithm, one can call
run(),>>> fdmrg = FiniteDMRG(mpo, bond_dim=20) >>> fdmrg.run(tol=1e-8)
The optimized
MatrixProductStatecan then be retrieved frommps.
Methods
__init__(mpo, bond_dim[, block_size, mps, ...])The Density Matrix Renormalization Group (DMRG) algorithm of a finite size system.
one_site_solver(site[, tol])- param site
perturb_wave_function(site[, alpha])Perturb the local tensor of wave function by an amount of
alpha.run([tol, max_sweep, metric])By calling this method,
FiniteDMRGwill start the sweeping procedure until the given tolerance is reached or touching the maximally allowed number of sweeps.sweep([direction, tol])Perform a single sweep on the given
direction.two_site_solver(site[, tol])variance()Attributes
- __init__(mpo, bond_dim, block_size=1, mps=None, exact_solver_dim=200)[source]#
The Density Matrix Renormalization Group (DMRG) algorithm of a finite size system.
- Parameters
mpo (tnpy.operators.MatrixProductOperator) – Matrix product operator.
bond_dim (int) – Bond dimensions.
block_size (int) – Block size of each local update.
mps (Optional[tnpy.matrix_product_state.MatrixProductState]) – (Optional) Initial guess of
MatrixProductState.exact_solver_dim (int) – A switch point at which the exact eigensolver should be used below certain matrix size.
Examples
To start the algorithm, one can call
run(),>>> fdmrg = FiniteDMRG(mpo, bond_dim=20) >>> fdmrg.run(tol=1e-8)
The optimized
MatrixProductStatecan then be retrieved frommps.
- property bond_dim: int#
- property mps: tnpy.matrix_product_state.MatrixProductState#
- property n_sites: int#
- property phys_dim: int#
- one_site_solver(site, tol=1e-08, **kwargs)[source]#
- Parameters
site (int) –
tol (float) – Tolerance to the eigensolver.
**kwargs – Keyword arguments to primme eigensolver.
- Return type
Tuple[float, numpy.ndarray]
Returns:
- perturb_wave_function(site, alpha=1e-05)[source]#
Perturb the local tensor of wave function by an amount of
alpha. For the single-site DMRGone_site_solver(), the modified density matrix can be used against trapping in local minimums. This is an in-place operation.- Parameters
site (int) –
alpha (float) – The small parameter for perturbation.
References
- sweep(direction=Direction.RIGHTWARD, tol=1e-08, **kwargs)[source]#
Perform a single sweep on the given
direction.- Parameters
direction (tnpy.matrix_product_state.Direction) – The left or right direction on which the sweep to perform.
tol (float) – Tolerance to the eigensolver.
**kwargs – Keyword arguments to primme eigensolver.
- Returns
Variationally optimized energy after this sweep.
- Return type
float | None
- run(tol=1e-08, max_sweep=100, metric=Metric.ENERGY, **kwargs)[source]#
By calling this method,
FiniteDMRGwill start the sweeping procedure until the given tolerance is reached or touching the maximally allowed number of sweeps.- Parameters
tol (float) – Required precision to the variationally optimized energy.
max_sweep (int) – Maximum number of sweeps.
metric (tnpy.finite_dmrg.Metric) – By which value to examine the convergence.
**kwargs – Keyword arguments to primme eigensolver.
- Returns
A record of energies computed on each sweep.
- Return type
List[float]
- property measurements: tnpy.matrix_product_state.MatrixProductStateMeasurements#