Пример #1
0
def test_canonical_form(bc, method):
    psi = random_MPS(8, 2, 6, form=None, bc=bc)
    psi2 = psi.copy()
    norm = np.sqrt(psi2.overlap(psi2, ignore_form=True))
    print("norm =", norm)
    psi2.norm /= norm  # normalize psi2
    norm2 = psi.overlap(psi2, ignore_form=True)
    print("norm2 =", norm2)
    assert abs(norm2 - norm) < 1.e-14 * norm
    meth = getattr(psi, method)
    meth(renormalize=False)  # psi.canonical_form_[infinite[2]]()
    psi.test_sanity()
    assert abs(psi.norm - norm) < 1.e-14 * norm
    psi.norm = 1.  # normalized psi
    ov = psi.overlap(psi2, ignore_form=True)
    print("normalized states: overlap <psi_canonical|psi> = 1.-", 1. - ov)
    assert abs(ov - 1.) < 1.e-14
    print("norm_test")
    print(psi.norm_test())
    assert np.max(psi.norm_test()) < 1.e-14  # SB = AS to good precision
    psi3 = psi.copy()
    # call canonical_form again, it shouldn't do anything now
    meth(renormalize=True)
    psi.test_sanity()
    ov = psi.overlap(psi3)
    assert abs(ov - 1.) < 1.e-14
    if method in ['canonical_form_finite', 'canonical_form_infinite2']:
        # check that A = SB S^-1 is orthonormal
        for i in range(psi.L):
            A = psi.get_B(i, 'A')
            c = npc.tensordot(A, A.conj(), axes=[['vL', 'p'], ['vL*', 'p*']])
            A_err = (c - npc.diag(1., c.legs[0])).norm()
            print(A_err)
            assert A_err < 1.e-13
Пример #2
0
def make_environment(psi, H, N=5):

    IdL = H.get_IdL(0)
    assert IdL is not None
    IdR = H.get_IdR(-1)
    assert IdR is not None
    vL = psi.get_B(0, None).get_leg('vL')
    vR = psi.get_B(0, None).get_leg('vR')
    wL = H.get_W(0).get_leg('wL')
    wR = wL.conj()
    dtype = psi.dtype

    LP = npc.diag(1., vR, dtype=dtype, labels=['vR*', 'vR'])
    LP = LP.add_leg(wR, IdL, axis=1, label='wR')
    RP = npc.diag(1., vL, dtype=dtype, labels=['vL*', 'vL'])
    RP = RP.add_leg(wL, IdR, axis=1, label='wL')
    A = psi.get_B(0, form='A')
    B = psi.get_B(0, form='B')
    W = H.get_W(0)

    TA = npc.tensordot(A, W, axes=('p', 'p'))
    TA = npc.tensordot(A.conj(), TA, axes=('p*', 'p*')).itranspose(
        ['vL*', 'vR*', 'wL', 'wR', 'vL', 'vR'])
    TB = npc.tensordot(B, W, axes=('p', 'p'))
    TB = npc.tensordot(B.conj(), TB, axes=('p*', 'p*')).itranspose(
        ['vL*', 'vR*', 'wL', 'wR', 'vL', 'vR'])

    for i in range(N):
        LP = npc.tensordot(LP,
                           TA,
                           axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
        RP = npc.tensordot(RP,
                           TB,
                           axes=(['vL', 'wL', 'vL*'], ['vR', 'wR', 'vR*']))

    return LP, RP
Пример #3
0
def test_lanczos_evolve(n, N_cache, tol=5.e-15):
    # generate Hermitian test array
    leg = gen_random_legcharge(ch, n)
    H = npc.Array.from_func_square(rmat.GUE, leg) - npc.diag(1., leg)
    H_flat = H.to_ndarray()
    H_Op = H  # use `matvec` of the array
    qtotal = leg.to_qflat()[0]
    psi_init = npc.Array.from_func(np.random.random, [leg], qtotal=qtotal)
    psi_init /= npc.norm(psi_init)
    psi_init_flat = psi_init.to_ndarray()
    lanc = lanczos.LanczosEvolution(H_Op, psi_init, {
        'verbose': 1,
        'N_cache': N_cache
    })
    for delta in [-0.1j, 0.1j, 1.j]:  #, 0.1, 1.]:
        psi_final_flat = expm(H_flat * delta).dot(psi_init_flat)
        psi_final, N = lanc.run(delta)
        ov = np.inner(psi_final.to_ndarray().conj(), psi_final_flat)
        ov /= np.linalg.norm(psi_final_flat)
        print("<psi1|psi1_flat>/norm=", ov)
        assert (abs(1. - abs(ov)) < tol)
Пример #4
0
def test_lanczos_evolve(n, N_cache, tol=5.e-14):
    # generate Hermitian test array
    leg = gen_random_legcharge(ch, n)
    H = npc.Array.from_func_square(rmat.GUE, leg) - npc.diag(1., leg)
    H_flat = H.to_ndarray()
    H_Op = H  # use `matvec` of the array
    qtotal = leg.to_qflat()[0]
    psi_init = npc.Array.from_func(np.random.random, [leg], qtotal=qtotal)
    #psi_init /= npc.norm(psi_init) # not necessary
    psi_init_flat = psi_init.to_ndarray()
    lanc = lanczos.LanczosEvolution(H_Op, psi_init, {'N_cache': N_cache})
    for delta in [-0.1j, 0.1j, 1.j, 0.1, 1.]:
        psi_final_flat = expm(H_flat * delta).dot(psi_init_flat)
        norm = np.linalg.norm(psi_final_flat)
        psi_final, N = lanc.run(delta, normalize=False)
        diff = np.linalg.norm(psi_final.to_ndarray() - psi_final_flat)
        print("norm(|psi_final> - |psi_final_flat>)/norm = ",
              diff / norm)  # should be 1.
        assert diff / norm < tol
        psi_final2, N = lanc.run(delta, normalize=True)
        assert npc.norm(psi_final / norm - psi_final2) < tol
Пример #5
0
W.iset_leg_labels(['wL', 'wR', 'p',
                   'p*'])  # wL/wR = virtual left/right of the MPO
Ws = [W] * L

print("3) define 'environments' left and right")

#  .---->- vR     vL ->----.
#  |                       |
#  envL->- wR     wL ->-envR
#  |                       |
#  .---->- vR*    vL*->----.

envL = npc.zeros(
    [W.get_leg('wL').conj(), Bs[0].get_leg('vL').conj(), Bs[0].get_leg('vL')])
envL.iset_leg_labels(['wR', 'vR', 'vR*'])
envL[0, :, :] = npc.diag(1., envL.legs[1])
envR = npc.zeros([
    W.get_leg('wR').conj(), Bs[-1].get_leg('vR').conj(), Bs[-1].get_leg('vR')
])
envR.iset_leg_labels(['wL', 'vL', 'vL*'])
envR[-1, :, :] = npc.diag(1., envR.legs[1])

print("4) contract MPS and MPO to calculate the energy <psi|H|psi>")
contr = envL
for i in range(L):
    # contr labels: wR, vR, vR*
    contr = npc.tensordot(contr, Bs[i], axes=('vR', 'vL'))
    # wR, vR*, vR, p
    contr = npc.tensordot(contr, Ws[i], axes=(['p', 'wR'], ['p*', 'wL']))
    # vR*, vR, wR, p
    contr = npc.tensordot(contr,