10. Calculation Configurations#

class renormalizer.utils.configs.BondDimDistri(value)[source]#

Bases: Enum

Bond dimension distribution

center_gauss = 'center gaussian'#

Guassian distribution which peaks at the center.

uniform = 'uniform'#

uniform distribution.

class renormalizer.utils.configs.CompressConfig(criteria: CompressCriteria = CompressCriteria.threshold, threshold: float = 0.001, bonddim_distri: BondDimDistri = BondDimDistri.uniform, max_bonddim: int = 32, vmethod: str = '2site', vprocedure=None, vrtol=1e-05, vguess_m=(5, 5), dump_matrix_size=inf, dump_matrix_dir='./', ofs: Optional[OFS] = None, ofs_swap_jw: bool = False)[source]#

Bases: object

MPS/MPO Compress Configuration.

Parameters:
  • criteria (CompressCriteria, optional) – The criteria for compression. Default is CompressCriteria.threshold.

  • threshold (float, optional) – The threshold to keep states if criteria is set to CompressCriteria.threshold or CompressCriteria.both. Default is \(10^{-3}\).

  • bonddim_distri (BondDimDistri, optional) – Bond dimension distribution if criteria is set to CompressCriteria.fixed or CompressCriteria.both. Default is BondDimDistri.uniform.

  • max_bonddim (int, optional) – Maximum bond dimension M under various bond dimension distributions. Default is 32.

  • vmethod (str, optional) –

    The algorithm used in variational compression. The default is 2site.

    • 1site: optimize one site at a time during sweep.

    • 2site: optimize two site at a time during sweep.

    Note

    It is recommended to use the 2site algorithm in variational optimization, though 1site algorithm is much cheaper. The reason is that 1site algorithm is easy to stuck into local minima while 2site algorithm is more robust. For complicated Hamiltonian, the problem is severer. The renormalized basis selection rule used here to partly solve the local minimal problem is according to Chan. J. Chem. Phys. 2004, 120, 3172. For efficiency, you could use 2site algorithm for several sweeps and then switch to 1site algorithm to converge the final mps.

  • vprocedure (list, optional) –

    The procedure to do variational compression. In the sublist, the first element is a CompressConfig, which controls the compression in each sweep. For backward compatibility, it can also be an int, which indicates CompressCriteria.fixed with the int as the max_bonddim. The Default is

    if vmethod == "1site":
        vprocedure = [[max_bonddim,1.0],[max_bonddim,0.7],
                      [max_bonddim,0.5],[max_bonddim,0.3],
                      [max_bonddim,0.1]] + [[max_bonddim,0],]*10
    else:
        vprocedure = [[max_bonddim,0.5],[max_bonddim,0.3],
                      [max_bonddim,0.1]] + [[max_bonddim,0],]*10
    

  • vrtol (float, optional) – The threshold of convergence in variational compression. Default is \(10^{-5}\).

  • vguess_m (tuple(int), optional) – The bond dimension of compressed_mpo and compressed_mps to construct the initial guess compressed_mpo @ compressed_mps in MatrixProduct.variational_compress. The default is (5,5).

  • dump_matrix_size (int or float, optional) –

    The total bytes threshold for dumping matrix into disks. Every local site matrix in the MPS/MPO with a size larger than the specified size will be moved from memory to the disk at dump_matrix_dir. The matrix will be moved back to memory upon access. The dumped matrices on the disks will be removed once the MPS/MPO is garbage-collected.

    Note

    If dump_matrix_size is set and the program exits abnormally, it is possible that the dumped matrices are not deleted automatically. A common case when this can happen is that the program is killed by a signal. In such cases it is recommended to check and clean dump_matrix_dir after the exit of the program since the dumped matrices may take a lot of disk space.

  • dump_matrix_dir (str, optional) – The directory to dump matrix when matrix is larger than dump_matrix_size.

  • ofs (OFS, optional) – Whether optimize the DOF ordering by OFS. The default value is None which means does not perform OFS.

  • ofs_swap_jw (bool, optional) – Whether swap the Jordan-Wigner transformation ordering when OFS is enabled. Set to True for and only for ab initio Hamiltonian constructed by the experimental renormalizer.model.h_qc.qc_model. Default is False.

See also

CompressCriteria

Compression criteria

renormalizer.mps.mp.MatrixProduct.variational_compress

Variational compression

renormalizer.mps.Mpo.contract

compress mpo @ mps/mpdm/mpo

property bonddim_should_set#
compute_m_trunc(sigma: ndarray, idx: int, left: bool) int[source]#
copy() CompressConfig[source]#
relax()[source]#
set_bonddim(length, max_value=None)[source]#
property threshold#
update(other: CompressConfig)[source]#
class renormalizer.utils.configs.CompressCriteria(value)[source]#

Bases: Enum

Criteria for compression.

both = 'both'#

Compress combining threshold and fixed. The bond dimension for the two criteria are both calculated and the smaller one is used.

fixed = 'fixed'#

Compress depending on pre-set fixed bond dimension.

threshold = 'threshold'#

Compress depending on pre-set threshold. States with singular value smaller than the threshold are discarded.

class renormalizer.utils.configs.EvolveConfig(method: EvolveMethod = EvolveMethod.prop_and_compress, adaptive=False, guess_dt=0.1, adaptive_rtol=0.0005, taylor_order: Optional[int] = None, rk_solver='C_RK4', reg_epsilon=1e-10, ivp_rtol=1e-05, ivp_atol=1e-08, ivp_solver='krylov', force_ovlp=True)[source]#

Bases: object

check_valid_dt(evolve_dt: complex)[source]#
copy()[source]#
property is_tdvp#
class renormalizer.utils.configs.EvolveMethod(value)[source]#

Bases: Enum

Time evolution methods.

prop_and_compress = 'P&C'#
prop_and_compress_tdrk = 'P&C TD RK'#
prop_and_compress_tdrk4 = 'P&C TD RK4'#
tdvp_mu_cmf = 'TDVP Matrix Unfolding Constant Mean Field'#
tdvp_mu_vmf = 'TDVP Matrix Unfolding Variable Mean Field'#
tdvp_ps = 'TDVP PS one-site'#
tdvp_ps2 = 'TDVP PS two-site'#
tdvp_vmf = 'TDVP Variable Mean Field'#
class renormalizer.utils.configs.OFS(value)[source]#

Bases: Enum

On the fly swapping criteria

ofs_d = 'OFS-D'#
ofs_debug = 'OFS-Debug'#
ofs_ds = 'OFS-D/S'#
ofs_s = 'OFS-S'#
class renormalizer.utils.configs.OptimizeConfig(procedure=None)[source]#

Bases: object

DMRG ground state algorithm optimization configuration

The procedure to do ground state optimization. In the sublist, the first element is a CompressConfig, which controls the compression in each sweep. For backward compatibility, it can also be an int, which indicates CompressCriteria.fixed with the int as the max_bonddim. The second element is the percent to choose the renormalied basis from each symmetry block to avoid trapping into local minimum.

copy()[source]#
renormalizer.utils.configs.parse_memory_limit(x) float[source]#