def diagNQt(dmrg, fbmps, fname, debug=False): t0 = time.time() nsite = dmrg.nsite prefixL = fname + 'L_site_' prefixR = fname + 'R_site_' # L->R sweeps mpo_dmrg_init.genBmat(dmrg, fname, -1) nii = numpy.zeros(nsite) for isite in range(0, nsite): fL = h5py.File(prefixL + str(isite - 1), "r") fR = h5py.File(prefixR + str(isite + 1), "r") bsite = mpo_dmrg_io.loadSite(fbmps, isite, dmrg.ifQt) tmpl = fL['mat'].value tmpr = fR['mat'].value npmat = mpo_dmrg_opers.genNpMat() # # i---*---j # / | \ # * L *m * R # \ | / # a---*---b # tmp = numpy.einsum('mn,lnr->lmr', npmat, bsite) tmp = numpy.einsum('Ll,lmr->Lmr', tmpl, tmp) tmp = numpy.einsum('lmr,rR->lmR', tmp, tmpr) nii[isite] = numpy.tensordot(tmp, bsite, axes=([0, 1, 2], [0, 1, 2])) fL.close() fR.close() # Final print ' sum of nii =', numpy.sum(nii) print ' nii =', nii t1 = time.time() print ' time for diagNQt = %.2f s' % (t1 - t0), ' rank =', dmrg.comm.rank return nii
def fmpsQtReverse(fmps0, fmps1, status, isym=2): print('\n[qtensor_api.fmpsQtReverse] status=', status) nsite = fmps0['nsite'].value fmps1['nsite'] = nsite # Check symmetry lenSym = len(fmps0['qnum0'].value[0]) if lenSym != isym: print(' error: isym is not consistent!') print(' qnum0:', fmps0['qnum0'].value[0]) print(' lenSym/isym =', lenSym, isym) exit() # Save qnums for isite in range(nsite + 1): ql = fmps0['qnum' + str(isite)].value fmps1['qnum' + str(isite)] = ql # Cast into numpy tensor for isite in range(nsite): qtsite = mpo_dmrg_io.loadSite(fmps0, isite, True) npsite = qtsite.toDenseTensor(qtsite.idlst) print(' isite =', isite, ' shape=', npsite.shape) fmps1['site' + str(isite)] = npsite return 0
def genCAopsNQt(norder, dmrg, fbmps, fkmps, fname, status, debug=False): if dmrg.comm.rank == 0: print '[mpo_dmrg_block.genCAopsNQt] status=', status print ' fname = ', fname t0 = time.time() bnsite = fbmps['nsite'].value knsite = fkmps['nsite'].value assert bnsite == knsite nsite = bnsite sbas = 2 * nsite prefix = fname + '_site_' # # On-site operators # sgnn = numpy.array([1., -1., -1., 1.]) cre = [0] * 2 cre[0] = mpo_dmrg_opers.genElemSpatialMat(0, 0, 1) cre[1] = mpo_dmrg_opers.genElemSpatialMat(1, 0, 1) ann = [0] * 2 ann[0] = mpo_dmrg_opers.genElemSpatialMat(0, 0, 0) ann[1] = mpo_dmrg_opers.genElemSpatialMat(1, 0, 0) cc2 = cre[0].dot(cre[1]) aa2 = ann[0].dot(ann[1]) ca2 = [cre[i].dot(ann[j]) for i in range(2) for j in range(2)] ac2 = [ann[i].dot(cre[j]) for i in range(2) for j in range(2)] # # L->R sweeps # if status == 'L': mpo_dmrg_init.genBmat(dmrg, fname, -1) for isite in range(0, nsite): if debug: print ' isite=', isite, ' of nsite=', nsite ti = time.time() f0 = h5py.File(prefix + str(isite - 1), "r") f1name = prefix + str(isite) f1 = h5py.File(f1name, "w") bsite = mpo_dmrg_io.loadSite(fbmps, isite, dmrg.ifQt) ksite = mpo_dmrg_io.loadSite(fkmps, isite, dmrg.ifQt) # 0. Qnums porbs = 2 * isite nqum = fkmps['qnum' + str(isite)].value[:, 0] sgnl = numpy.power(-1.0, nqum) ksite2 = numpy.einsum('l,lnr->lnr', sgnl, ksite) # 1. Similar to mpo_dmrg_init.genSopsNQt: <li[A]|lj[B]> oplst = ['mat'] renorm_l1(bsite, ksite, f0, f1, oplst) # 2. Part-1: <l'n'|op|ln> = <l'|op|l>*<n'|n> if p < 2k oplst = ['op_C_'+str(p) for p in range(porbs)]\ + ['op_A_'+str(p) for p in range(porbs)] renorm_l1(bsite, ksite, f0, f1, oplst) # Part-2: <l'n'|op|ln> = <l'|l>(-1)^l*<n'|op|n> if p = 2k,2k+1 # <L'|op|L> = A[l'n',L']*<l'|l><n'|op|n>*A[ln,L] tmp0 = renorm_l2_absorb_lop(bsite, f0['mat']) for ispin in range(2): tmp = renorm_l2_transform_nop(tmp0, cre[ispin], ksite2) f1['op_C_' + str(2 * isite + ispin)] = tmp tmp = renorm_l2_transform_nop(tmp0, ann[ispin], ksite2) f1['op_A_' + str(2 * isite + ispin)] = tmp # 3. Part-1: <l'n'|aa(p<q)|ln> = <l'|apq|l><n'|n> if p,q < 2k oplst = ['op_CC_'+str(p)+'_'+str(q) for p in range(porbs) for q in range(p+1,porbs)]\ + ['op_AA_'+str(p)+'_'+str(q) for p in range(porbs) for q in range(p+1,porbs)] renorm_l1(bsite, ksite, f0, f1, oplst) # Part-2: <l'|ap|l>*<n'|aq|n> if p < 2k and q in 2k,2k+1 for p in range(porbs): # CC tmp0 = renorm_l2_absorb_lop(bsite, f0['op_C_' + str(p)]) for ispin in range(2): tmp = renorm_l2_transform_nop(tmp0, cre[ispin], ksite2) f1['op_CC_' + str(p) + '_' + str(2 * isite + ispin)] = tmp # AA tmp0 = renorm_l2_absorb_lop(bsite, f0['op_A_' + str(p)]) for ispin in range(2): tmp = renorm_l2_transform_nop(tmp0, ann[ispin], ksite2) f1['op_AA_' + str(p) + '_' + str(2 * isite + ispin)] = tmp # Part-3: <l'|l>*<n'|apq|n> if p,q in 2k,2k+1 # <L'|opq|L> = A[l'n',L']*<l'|l><n'|op|n>*A[ln,L] tmp0 = renorm_l2_absorb_lop(bsite, f0['mat']) tmp = renorm_l2_transform_nop(tmp0, cc2, ksite) f1['op_CC_' + str(2 * isite + 0) + '_' + str(2 * isite + 1)] = tmp tmp = renorm_l2_transform_nop(tmp0, aa2, ksite) f1['op_AA_' + str(2 * isite + 0) + '_' + str(2 * isite + 1)] = tmp # 4. Part-1: <l'n'|cp*aq|ln> = <l'|cp*aq|l><n'|n> if p,q < 2k oplst = ['op_CA_'+str(p)+'_'+str(q) for p in range(porbs) for q in range(porbs)]\ + ['op_AC_'+str(p)+'_'+str(q) for p in range(porbs) for q in range(porbs)] renorm_l1(bsite, ksite, f0, f1, oplst) # Part-2: cp[l]*aq[n] or cq[n]*ap[l]=-ap[l]*cq[n] # ap[l]*cq[n] or aq[n]*cp[l]=-cp[l]*aq[n] # for p < 2*k and q in 2k,2k+1 for p in range(porbs): # CA tmp0 = renorm_l2_absorb_lop(bsite, f0['op_C_' + str(p)]) for ispin in range(2): tmp = renorm_l2_transform_nop(tmp0, ann[ispin], ksite2) f1['op_CA_' + str(p) + '_' + str(2 * isite + ispin)] = tmp f1['op_AC_' + str(2 * isite + ispin) + '_' + str(p)] = -tmp # AC tmp0 = renorm_l2_absorb_lop(bsite, f0['op_A_' + str(p)]) for ispin in range(2): tmp = renorm_l2_transform_nop(tmp0, cre[ispin], ksite2) f1['op_AC_' + str(p) + '_' + str(2 * isite + ispin)] = tmp f1['op_CA_' + str(2 * isite + ispin) + '_' + str(p)] = -tmp # Part-3: <l'|l>*<n'|cp*aq|n> if p,q in 2k,2k+1 tmp0 = renorm_l2_absorb_lop(bsite, f0['mat']) for ispin in range(2): for jspin in range(2): tmp = renorm_l2_transform_nop(tmp0, ca2[ispin * 2 + jspin], ksite) f1['op_CA_' + str(2 * isite + ispin) + '_' + str(2 * isite + jspin)] = tmp tmp = renorm_l2_transform_nop(tmp0, ac2[ispin * 2 + jspin], ksite) f1['op_AC_' + str(2 * isite + ispin) + '_' + str(2 * isite + jspin)] = tmp # final isite f0.close() f1.close() tf = time.time() if dmrg.comm.rank == 0: print ' isite =', os.path.split( f1name)[-1], ' t = %.2f s' % (tf - ti) # # L<-R sweeps: For properties, left canonical form is assumed even for R sweeps! # elif status == 'R': mpo_dmrg_init.genBmat(dmrg, fname, nsite) for isite in range(nsite - 1, -1, -1): if debug: print ' isite=', isite, ' of nsite=', nsite ti = time.time() f0 = h5py.File(prefix + str(isite + 1), "r") f1name = prefix + str(isite) f1 = h5py.File(f1name, "w") bsite = mpo_dmrg_io.loadSite(fbmps, isite, dmrg.ifQt) ksite = mpo_dmrg_io.loadSite(fkmps, isite, dmrg.ifQt) # 0. Qnums porbs = 2 * (isite + 1) ksite2 = numpy.einsum('n,lnr->lnr', sgnn, ksite) # 1. <ri[A]|rj[B]> oplst = ['mat'] renorm_r1(bsite, ksite, f0, f1, oplst) # 2. Part-1: <n'r'|op|nr> = <n'|n>(-1)^n*<r'|op|r> if p > 2k # A[l',n',r']<n'|n>(-1)^n*<r'|op|r>*A[l,n,r] oplst = ['op_C_'+str(p) for p in range(porbs,sbas)]\ + ['op_A_'+str(p) for p in range(porbs,sbas)] renorm_r1(bsite, ksite2, f0, f1, oplst) # Part-2: <n'r'|op|nr> = <n'|op|n>*<r'|r> if p = 2k,2k+1 # <l'|op|l> = A[l',n'r']*<n'|op|n>*<r'|r>*A[l,nr] tmp0 = renorm_r2_absorb_rop(ksite, f0['mat']) for ispin in range(2): tmp = renorm_r2_transform_nop(tmp0, cre[ispin], bsite) f1['op_C_' + str(2 * isite + ispin)] = tmp tmp = renorm_r2_transform_nop(tmp0, ann[ispin], bsite) f1['op_A_' + str(2 * isite + ispin)] = tmp # 3. Part-1: <n'r'|aa(p<q)|n'r'> = <n'|n><r'|apq|r> if p,q > 2k # <l'|aa(p<q)|l> = A[l',n',r']<n'|n><r'|apq|r>A[l,n,r] oplst = ['op_CC_'+str(p)+'_'+str(q) for p in range(porbs,sbas) for q in range(p+1,porbs)]\ + ['op_AA_'+str(p)+'_'+str(q) for p in range(porbs,sbas) for q in range(p+1,porbs)] renorm_r1(bsite, ksite, f0, f1, oplst) # Part-2: <n'r'|aa(p<q)|n'r'> = <n'|ap|n>*S*<r'|aq|r> if p in 2k,2k+1 and q > 2k # <l'|aa(p<q)|l> = A[l',n',r']*<n'|ap|n>S*<r'|aq|r>*A[l,n,r] for p in range(porbs, sbas): # CC tmp0 = renorm_r2_absorb_rop(ksite2, f0['op_C_' + str(p)]) for ispin in range(2): tmp = renorm_r2_transform_nop(tmp0, cre[ispin], bsite) f1['op_CC_' + str(2 * isite + ispin) + '_' + str(p)] = tmp # AA tmp0 = renorm_r2_absorb_rop(ksite2, f0['op_A_' + str(p)]) for ispin in range(2): tmp = renorm_r2_transform_nop(tmp0, ann[ispin], bsite) f1['op_AA_' + str(2 * isite + ispin) + '_' + str(p)] = tmp # Part-3: <n'r'|aa(p<q)|n'r'> = <n'|apq|n>*<r'|r> if p,q in 2k,2k+1 # <l'|aa(p<q)|l> = A[l',n',r]*<n'|apq|n>*<r'|r>*A[l,n,r] tmp0 = renorm_r2_absorb_rop(ksite, f0['mat']) tmp = renorm_r2_transform_nop(tmp0, cc2, bsite) f1['op_CC_' + str(2 * isite + 0) + '_' + str(2 * isite + 1)] = tmp tmp = renorm_r2_transform_nop(tmp0, aa2, bsite) f1['op_AA_' + str(2 * isite + 0) + '_' + str(2 * isite + 1)] = tmp # 4. Part-1: <n'r'|cp*aq|nr> = <n'|n><r'|cp*aq|r> if p,q > 2k # <l'|cp*aq|l> = A[l',nr']*<r'|cp*aq|r>*A[l,nr] oplst = ['op_CA_'+str(p)+'_'+str(q) for p in range(porbs,sbas) for q in range(porbs,sbas)]\ + ['op_AC_'+str(p)+'_'+str(q) for p in range(porbs,sbas) for q in range(porbs,sbas)] renorm_r1(bsite, ksite, f0, f1, oplst) # Part-2: <n'r'|cp*aq|nr> = <n'|cp|n>S*<r'|aq|r> for p in range(porbs, sbas): # CA tmp0 = renorm_r2_absorb_rop(ksite2, f0['op_A_' + str(p)]) for ispin in range(2): tmp = renorm_r2_transform_nop(tmp0, cre[ispin], bsite) f1['op_CA_' + str(2 * isite + ispin) + '_' + str(p)] = tmp f1['op_AC_' + str(p) + '_' + str(2 * isite + ispin)] = -tmp # AA tmp0 = renorm_r2_absorb_rop(ksite2, f0['op_C_' + str(p)]) for ispin in range(2): tmp = renorm_r2_transform_nop(tmp0, ann[ispin], bsite) f1['op_AC_' + str(2 * isite + ispin) + '_' + str(p)] = tmp f1['op_CA_' + str(p) + '_' + str(2 * isite + ispin)] = -tmp # Part-3: <n'r'|cp*aq|n'r'> = <n'|cp*aq|n>*<r'|r> if p,q in 2k,2k+1 tmp0 = renorm_r2_absorb_rop(ksite, f0['mat']) for ispin in range(2): for jspin in range(2): tmp = renorm_r2_transform_nop(tmp0, ca2[ispin * 2 + jspin], bsite) f1['op_CA_' + str(2 * isite + ispin) + '_' + str(2 * isite + jspin)] = tmp tmp = renorm_r2_transform_nop(tmp0, ac2[ispin * 2 + jspin], bsite) f1['op_AC_' + str(2 * isite + ispin) + '_' + str(2 * isite + jspin)] = tmp # final isite f0.close() f1.close() tf = time.time() if dmrg.comm.rank == 0: print ' isite =', os.path.split( f1name)[-1], ' t = %.2f s' % (tf - ti) # CHECK f = h5py.File(f1name, 'r') rdm1 = numpy.zeros((sbas, sbas)) hdm1 = numpy.zeros((sbas, sbas)) for i in range(sbas): for j in range(sbas): rdm1[i, j] = f['op_CA_' + str(i) + '_' + str(j)].value hdm1[j, i] = f['op_AC_' + str(i) + '_' + str(j)].value sab = f['mat'].value[0, 0] #print rdm1+hdm1 print ' ovlap=', sab print ' P+H-I=', numpy.linalg.norm(rdm1 + hdm1 - numpy.identity(sbas) * sab) print ' skewP=', numpy.linalg.norm(rdm1 - rdm1.T) print ' skewH=', numpy.linalg.norm(hdm1 - hdm1.T) print ' trace=', numpy.trace(rdm1) f.close() t1 = time.time() dmrg.comm.Barrier() print ' time for genHops = %.2f s' % (t1 - t0), ' rank =', dmrg.comm.rank return rdm1
def evalPropsNQt(dmrg,fbmps,fkmps,fop,status,debug=False): t0 = time.time() # sites bnsite = fbmps['nsite'].value knsite = fkmps['nsite'].value assert bnsite == knsite nsite = bnsite nop = fop['nop'].value if debug: print ' nop = ',nop print ' nsite = ',nsite # Create directory - maybe we could use dmrg.path+'/tmpdirProps' in future path = dmrg.path+'/tmpdirProps' sysutil_io.createDIR(path,debug) fname = path+'/trop' prefix = fname+'_site_' # L->R sweeps if status == 'L': mpo_dmrg_init.genBopsNQt(fname,nop,-1) for isite in range(0,nsite): if debug: print ' isite=',isite,' of nsite=',nsite ti = time.time() f0 = h5py.File(prefix+str(isite-1),"r") f1name = prefix+str(isite) f1 = h5py.File(f1name,"w") bsite = mpo_dmrg_io.loadSite(fbmps,isite,False) ksite = mpo_dmrg_io.loadSite(fkmps,isite,False) if isite == nsite-1: exphop = numpy.zeros(nop,dtype=dmrg_dtype) for iop in range(nop): if debug: print ' iop=',iop,' of nop=',nop cop = fop['site'+str(isite)+'/op'+str(iop)].value tmp = f0['opers'+str(iop)].value #--- kernel --- tmp = numpy.tensordot(bsite.conj(),tmp,axes=([0],[1])) tmp = numpy.tensordot(cop,tmp,axes=([0,2],[2,0])) tmp = numpy.tensordot(tmp,ksite,axes=([1,3],[1,0])) #--- kernel --- f1['opers'+str(iop)] = tmp # Get the expectation value <psi|Hx|psi> at the boundary if isite == nsite-1: exphop[iop] = tmp[0,0,0] f0.close() f1.close() # final isite tf = time.time() if debug: print ' isite =',os.path.split(f1name)[-1],\ ' nop =',nop,' t = %.2f s'%(tf-ti) elif status == 'R': mpo_dmrg_init.genBopsNQt(fname,nop,nsite) for isite in range(nsite-1,-1,-1): if debug: print ' isite=',isite,' of nsite=',nsite ti = time.time() f0 = h5py.File(prefix+str(isite+1),"r") f1name = prefix+str(isite) f1 = h5py.File(f1name,"w") bsite = mpo_dmrg_io.loadSite(fbmps,isite,False) ksite = mpo_dmrg_io.loadSite(fkmps,isite,False) if isite == 0: exphop = numpy.zeros(nop,dtype=dmrg_dtype) for iop in range(nop): if debug: print ' iop=',iop,' of nop=',nop cop = fop['site'+str(isite)+'/op'+str(iop)].value tmp = f0['opers'+str(iop)].value #--- kernel --- tmp = numpy.tensordot(bsite.conj(),tmp,axes=([2],[1])) tmp = numpy.tensordot(cop,tmp,axes=([1,2],[2,1])) tmp = numpy.tensordot(tmp,ksite,axes=([1,3],[1,2])) #--- kernel --- f1['opers'+str(iop)] = tmp # Get the expectation value <psi|Hx|psi> at the boundary if isite == 0: exphop[iop] = tmp[0,0,0] f0.close() f1.close() # final isite tf = time.time() if debug: print ' isite =',os.path.split(f1name)[-1],\ ' nop =',nop,' t = %.2f s'%(tf-ti) # final sysutil_io.deleteDIR(path,1,debug) t1=time.time() if debug: print ' time for evalProps = %.2f s'%(t1-t0) return exphop
def evalPropsQt(dmrg,fbmps,fkmps,fop,status,debug=False): t0 = time.time() # sites bnsite = fbmps['nsite'].value knsite = fkmps['nsite'].value assert bnsite == knsite nsite = bnsite nop = fop['nop'].value if debug: print ' nop = ',nop print ' nsite = ',nsite # Create directory path = dmrg.path+'/tmpdirProps' sysutil_io.createDIR(path,debug) fname = path+'/trop' prefix = fname+'_site_' # L->R sweeps if status == 'L': mpo_dmrg_initQt.genBopsQt(dmrg,fname,nop,-1,ifslc=False) lop = qtensor.qtensor() cop = qtensor.qtensor() for isite in range(0,nsite): if debug: print ' isite=',isite,' of nsite=',nsite ti = time.time() f0 = h5py.File(prefix+str(isite-1),"r") f1name = prefix+str(isite) f1 = h5py.File(f1name,"w") bsite = mpo_dmrg_io.loadSite(fbmps,isite,True) ksite = mpo_dmrg_io.loadSite(fkmps,isite,True) # Lower the symmetry if dmrg.ifs2proj: bsite = bsite.reduceQsymsToN() ksite = ksite.reduceQsymsToN() if isite == nsite-1: exphop = numpy.zeros(nop,dtype=dmrg_dtype) for iop in range(nop): if debug: print ' iop=',iop,' of nop=',nop # COP cop.load(fop,'site'+str(isite)+'/op'+str(iop)) # LOP lop.load(f0,'opers'+str(iop)) #--- kernel --- tmp = qtensor.tensordot(bsite,lop,axes=([0],[1]),ifc1=True) tmp = qtensor.tensordot(cop,tmp,axes=([0,2],[2,0])) tmp = qtensor.tensordot(tmp,ksite,axes=([1,3],[1,0])) #--- kernel --- tmp.dump(f1,'opers'+str(iop)) # Get the expectation value <psi|Hx|psi> at the boundary if isite == nsite-1: exphop[iop] = tmp.value f0.close() f1.close() # final isite tf = time.time() if debug: print ' isite =',os.path.split(f1name)[-1],\ ' nop =',nop,' t = %.2f s'%(tf-ti) elif status == 'R': mpo_dmrg_initQt.genBopsQt(dmrg,fname,nop,nsite,ifslc=False) rop = qtensor.qtensor() cop = qtensor.qtensor() for isite in range(nsite-1,-1,-1): if debug: print ' isite=',isite,' of nsite=',nsite ti = time.time() f0 = h5py.File(prefix+str(isite+1),"r") f1name = prefix+str(isite) f1 = h5py.File(f1name,"w") bsite = mpo_dmrg_io.loadSite(fbmps,isite,dmrg.ifQt) ksite = mpo_dmrg_io.loadSite(fkmps,isite,dmrg.ifQt) # Lower the symmetry if dmrg.ifs2proj: bsite = bsite.reduceQsymsToN() ksite = ksite.reduceQsymsToN() if isite == 0: exphop = numpy.zeros(nop,dtype=dmrg_dtype) for iop in range(nop): if debug: print ' iop=',iop,' of nop=',nop # COP cop.load(fop,'site'+str(isite)+'/op'+str(iop)) cop.qsyms[0] = -cop.qsyms[0] cop.qsyms[1] = -cop.qsyms[1] cop.status[0] = ~cop.status[0] cop.status[1] = ~cop.status[1] # ROP rop.load(f0,'opers'+str(iop)) #--- kernel --- tmp = qtensor.tensordot(bsite,rop,axes=([2],[1]),ifc1=True) tmp = qtensor.tensordot(cop,tmp,axes=([1,2],[2,1])) tmp = qtensor.tensordot(tmp,ksite,axes=([1,3],[1,2])) #--- kernel --- tmp.dump(f1,'opers'+str(iop)) # Get the expectation value <psi|Hx|psi> at the boundary if isite == 0: exphop[iop] = tmp.value f0.close() f1.close() # final isite tf = time.time() if debug: print ' isite =',os.path.split(f1name)[-1],\ ' nop =',nop,' t = %.2f s'%(tf-ti) # final sysutil_io.deleteDIR(path,1,debug) t1=time.time() if debug: print ' time for evalProps = %.2f s'%(t1-t0) return exphop
def left_sweep_projection(flmps0,flmps1,qtarget,thresh,Dcut,debug,\ ifBlockSingletEmbedding=False,\ ifBlockSymScreen=False,\ ifpermute=False): nsite,ne,ms,sval = qtarget print '\n[mpo_dmrg_samps.left_sweep_projection] (nsite,ne,ms,sval)=',(nsite,ne,ms,sval),\ 'ifBlockSE=',ifBlockSingletEmbedding flmps1['nsite'] = nsite # Reduced qnums - the basis states are ordered to maximize the efficiency. qnumsN = numpy.array([[0.,0.,1.],[1.,0.5,1.],[2.,0.,1.]]) qnums1 = [None]*(nsite+1) # Do not use singlet embedding for the first purification sweep, # otherwise, it will mess up with other spins. if ifBlockSingletEmbedding == False: qnums1[0] = numpy.array([[0.,0.,1.]]) wmat = numpy.identity(1) # Quantum numbers ne_eff = ne sval_eff = sval ms_eff = ms else: qnums1[0] = numpy.array([[2*sval,sval,1.]]) # (N,S)=(2*S,S) ndeg = int(2*sval+1) wmat = numpy.zeros((ndeg,1)) # Coupled to a noninteracting state |S(-M)> # to create a broken symetry state |S(-M)>*|SM> # which has the singlet component |00>/sqrt(2*S+1). im = int(-ms+sval) wmat[im] = 1.0 # Quantum numbers ne_eff = ne + 2.0*sval sval_eff = 0.0 ms_eff = 0.0 # Start conversion for the following sites for isite in range(nsite): site0 = mpo_dmrg_io.loadSite(flmps0,isite,False) # Expand |M> basis to |SM> basis wA0 = numpy.tensordot(wmat,site0,axes=([1],[0])) # Ll,lNr->LNr # Transform to combined basis if ifpermute and ifBlockSingletEmbedding and isite==0: wA0 = wA0.transpose(1,0,2) sgn = numpy.array([1.,(-1.0)**(int(2*sval)),(-1.0)**(int(2*sval)),1.]) wA0 = numpy.einsum('n,nvr->nvr',sgn,wA0) qnumsL,wA1 = util_tensor.spinCouple(wA0,qnumsN,qnums1[isite]) else: qnumsL,wA1 = util_tensor.spinCouple(wA0,qnums1[isite],qnumsN) # Quasi-RDM (reduced): each dim is {[(SaSb)S,na*nb]} classes,qrdmL = util_tensor.quasiRDM(wA1,qnumsL) if debug: print 'trace=',numpy.trace(qrdmL),qrdmL.shape,len(classes) # Decimation: [(SaSb)S,na*nb]=>[S,r] for all possible Sa,Sb # Therefore, rotL is a matrix with dimensions [(SaSb)S,na*nb]*[S,r]! if ifBlockSingletEmbedding and isite==0: dwts,qred,rotL,sigs = mpo_dmrg_qparser.rdm_blkdiag(qrdmL,classes,-1.e-10,-1,debug) else: dwts,qred,rotL,sigs = mpo_dmrg_qparser.rdm_blkdiag(qrdmL,classes,thresh,Dcut,debug) # # | (0) w[i-1] = <l[i-1]|m[i-1]>: expansion coeff of |l(NM)> to |m(NSM)> # |-----|------| (1) update formula for w[i] # | L | (2) L[i] - isometry # | | | (3) L[i]*t[CG] - new site A[i] # | t[CG] | # | / \ | | # | w---A[i]---A[i+1]--- # |------------| # # Screen qnumsl,rotL = util_spinsym.symScreen(isite,nsite,rotL,qred,sigs,ne_eff,sval_eff,debug,\ ifBlockSymScreen=ifBlockSymScreen) if ifpermute and ifBlockSingletEmbedding and isite==0: # A[lnr] site1 = util_tensor.expandRotL(rotL,qnumsN,qnums1[isite],qnumsL,qnumsl) # srotR = rotL^\dagger * wA srotR = numpy.tensordot(site1,wA0,axes=([0,1],[0,1])) # LNR,LNr->Rr site1 = site1.transpose(1,0,2) else: # A[lnr] site1 = util_tensor.expandRotL(rotL,qnums1[isite],qnumsN,qnumsL,qnumsl) # srotR = rotL^\dagger * wA srotR = numpy.tensordot(site1,wA0,axes=([0,1],[0,1])) # LNR,LNr->Rr # Print: print ' ---> isite=',isite,'site0=',site0.shape,'rotL=',rotL.shape,\ 'site1=',site1.shape #,'srotR=',srotR.shape print # Check left canonical form if debug: tmp = numpy.tensordot(site1,site1,axes=([0,1],[0,1])) # lna,lnb->ab d = tmp.shape[0] diff = numpy.linalg.norm(tmp-numpy.identity(d)) print ' diff=',diff assert diff < 1.e-10 # Save if isite < nsite-1: wmat = srotR.copy() qnums1[isite+1] = qnumsl.copy() flmps1.create_dataset('rotL'+str(isite),data=rotL ,compression='lzf') flmps1.create_dataset('site'+str(isite),data=site1,compression='lzf') else: # Special treatment for the last site by invoking symmetry selection! info,qnumsl,rotL,site1,wt = util_tensor.lastSite(rotL,site1,srotR,qnumsl,ne_eff,ms_eff,sval_eff) if info: qnums1[nsite] = qnumsl.copy() flmps1.create_dataset('rotL'+str(isite),data=rotL ,compression='lzf') flmps1.create_dataset('site'+str(isite),data=site1,compression='lzf') # Dump on flmps1 if info: print '\nfinalize left_sweep ...' qnumsm = [None]*(nsite+1) for isite in range(nsite): qnumsm[isite] = util_spinsym.expandSM(qnums1[isite]) # Left case qnumsm[nsite] = numpy.array([[ne,ms]]) # DUMP qnums for isite in range(nsite+1): flmps1['qnum'+str(isite)] = qnumsm[isite] flmps1['qnumNS'+str(isite)] = qnums1[isite] # Check for isite in range(nsite+1): dimr = util_spinsym.dim_red(qnums1[isite]) print ' ibond/dim(NS)/dim(NSM)=',isite,dimr,len(qnumsm[isite]) return wt
def right_sweep_projection(flmps0,flmps1,qtarget,thresh,Dcut,debug,\ ifBlockSingletEmbedding=False): nsite,ne,ms,sval = qtarget print '\n[mpo_dmrg_samps.right_sweep_projection] (nsite,,ne,ms,sval)=',(nsite,ne,ms,sval) flmps1['nsite'] = nsite # Reduced qnums - the basis states are ordered to maximize the efficiency. qnumsN = numpy.array([[0.,0.,1.],[1.,0.5,1.],[2.,0.,1.]]) qnums1 = [None]*(nsite+1) # This works perfectly. qnums1[nsite] = numpy.array([[0.,0.,1.]]) wmat = numpy.identity(1) for isite in range(nsite-1,-1,-1): site0 = mpo_dmrg_io.loadSite(flmps0,isite,False) # *** Change the direction for combinations site0 = numpy.einsum('lnr->rnl',site0) # Expand |M> basis to |SM> basis wA0 = numpy.tensordot(wmat,site0,axes=([1],[0])) # Ll,lNr->LNr # Transform to combined basis qnumsL,wA1 = util_tensor.spinCouple(wA0,qnums1[isite+1],qnumsN) # Quasi-RDM (reduced): each dim is {[(SaSb)S,na*nb]} classes,qrdmL = util_tensor.quasiRDM(wA1,qnumsL) if debug: print 'trace=',numpy.trace(qrdmL),qrdmL.shape,len(classes) # Decimation: [(SaSb)S,na*nb]=>[S,r] for all possible Sa,Sb # Therefore, rotL is a matrix with dimensions [(SaSb)S,na*nb]*[S,r]! dwts,qred,rotL,sigs = mpo_dmrg_qparser.rdm_blkdiag(qrdmL,classes,thresh,Dcut,debug) # # | (0) w[i-1] = <l[i-1]|m[i-1]>: expansion coeff of |l(NM)> to |m(NSM)> # |-----|------| (1) update formula for w[i] # | L | (2) L[i] - isometry # | | | (3) L[i]*t[CG] - new site A[i] # | t[CG] | # | / \ | | # | w---A[i]---A[i+1]--- # |------------| # # Screen qnumsl,rotL = util_spinsym.symScreen(isite,nsite,rotL,qred,sigs,ne,sval,debug,status='R') # A[lnr] site1 = util_tensor.expandRotL(rotL,qnums1[isite+1],qnumsN,qnumsL,qnumsl) # srotR = rotL^\dagger * wA srotR = numpy.tensordot(site1,wA0,axes=([0,1],[0,1])) # LNR,LNr->Rr # Print: print ' ---> isite=',isite,'site0=',site0.shape,'rotL=',rotL.shape,\ 'site1=',site1.shape #,'srotR=',srotR.shape print # Check left canonical form if debug: tmp = numpy.tensordot(site1,site1,axes=([0,1],[0,1])) # lna,lnb->ab d = tmp.shape[0] diff = numpy.linalg.norm(tmp-numpy.identity(d)) print ' diff=',diff assert diff < 1.e-10 # Save if isite > 0: wmat = srotR.copy() qnums1[isite] = qnumsl.copy() flmps1.create_dataset('rotL'+str(isite),data=rotL ,compression='lzf') # *** Change the direction back [Right canonical form] site1 = numpy.einsum('rnl->lnr',site1) flmps1.create_dataset('site'+str(isite),data=site1,compression='lzf') else: # In case of singlet embedding, the in bond of the first site # is expanded in the subroutine expandRotL in the left sweep for nonsinglet state. if ifBlockSingletEmbedding and abs(sval)>1.e-10: im = int(ms+sval) srotR = srotR[:,im].reshape(srotR.shape[0],1) # Special treatment for the last site by invoking symmetry selection! info,qnumsl,rotL,site1,wt = util_tensor.lastSite(rotL,site1,srotR,qnumsl,ne,ms,sval) if info: qnums1[0] = qnumsl.copy() flmps1.create_dataset('rotL'+str(isite),data=rotL ,compression='lzf') site1 = numpy.einsum('rnl->lnr',site1) flmps1.create_dataset('site'+str(isite),data=site1,compression='lzf') # Dump on flmps1 if info: print '\nfinalize right_sweep ...' qnumsm = [None]*(nsite+1) for isite in range(1,nsite+1): qnumsm[isite] = util_spinsym.expandSM(qnums1[isite]) # Right case qnumsm[0] = numpy.array([[ne,ms]]) # DUMP qnums for isite in range(nsite+1): flmps1['qnum'+str(isite)] = qnumsm[isite] flmps1['qnumNS'+str(isite)] = qnums1[isite] # Check for isite in range(nsite+1): dimr = util_spinsym.dim_red(qnums1[isite]) print ' ibond/dim(NS)/dim(NSM)=',isite,dimr,len(qnumsm[isite]) return wt