{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "# Tensor Network Renormalization Group for the 2D Ising model" ] }, { "cell_type": "code", "execution_count": 1, "metadata": {}, "outputs": [], "source": [ "using LinearAlgebra" ] }, { "cell_type": "code", "execution_count": 2, "metadata": {}, "outputs": [ { "name": "stderr", "output_type": "stream", "text": [ "┌ Info: Recompiling stale cache file /home/ssteinha/.julia/compiled/v1.0/Plots/ld3vC.ji for Plots [91a5bcdd-55d7-5caf-9e0b-520d859cae80]\n", "└ @ Base loading.jl:1187\n" ] } ], "source": [ "using Plots" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Define the basic functions for the Ising model" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Define Kronecker Delta" ] }, { "cell_type": "code", "execution_count": 3, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "kronecker_delta (generic function with 1 method)" ] }, "execution_count": 3, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function kronecker_delta(i::Int,j::Int)\n", " \n", " res = 0.\n", " \n", " if i==j\n", " \n", " res = 1.\n", " \n", " end\n", " \n", " res\n", " \n", "end" ] }, { "cell_type": "code", "execution_count": 4, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "1.0" ] }, "execution_count": 4, "metadata": {}, "output_type": "execute_result" } ], "source": [ "kronecker_delta(1,1)" ] }, { "cell_type": "code", "execution_count": 5, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "0.0" ] }, "execution_count": 5, "metadata": {}, "output_type": "execute_result" } ], "source": [ "kronecker_delta(1,2)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Define the weight functions of Ising model" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "$Z_2$ has two representations, $k=0$ and $k=1$. The weight on an edge is $\\cosh(\\beta)$ for $k=0$, $\\sinh(\\beta)$ for $k=1$." ] }, { "cell_type": "code", "execution_count": 6, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "weight (generic function with 1 method)" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function weight(beta::Float64,k::Int)\n", " \n", " res = 1/2 * (exp(beta) + (-1)^k * exp(-beta))\n", " \n", " res\n", " \n", "end" ] }, { "cell_type": "code", "execution_count": 7, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "2.352409615243247" ] }, "execution_count": 7, "metadata": {}, "output_type": "execute_result" } ], "source": [ "weight(1.5,0)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Define the initial tensor" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The tensor is 4-valent, each edge carrying a representation of $Z_2$, and a square root of an edge weight (for $\\beta$ and $k$). On the vertex, we have a Kronecker Delta of the sum of the 4 representations meeting at the vertex modulo 2; if this sum is 1, it vanishes." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Note that Julia starts counting array indices with 1, not 0. You need to account for that." ] }, { "cell_type": "code", "execution_count": 8, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "init_tensor (generic function with 1 method)" ] }, "execution_count": 8, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function init_tensor(beta::Float64)\n", " \n", " tensor = zeros(2,2,2,2)\n", " \n", " for i in 1:2, j in 1:2, k in 1:2, l in 1:2\n", " \n", " tensor[i,j,k,l] = sqrt(weight(beta,i-1) * weight(beta,j-1) * weight(beta,k-1) * weight(beta,l-1)) *\n", " kronecker_delta(mod(i+j+l+k,2),0)\n", " \n", " end\n", " \n", " return tensor\n", " \n", "end" ] }, { "cell_type": "code", "execution_count": 9, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "2×2×2×2 Array{Float64,4}:\n", "[:, :, 1, 1] =\n", " 2.3811 0.0 \n", " 0.0 1.81343\n", "\n", "[:, :, 2, 1] =\n", " 0.0 1.81343\n", " 1.81343 0.0 \n", "\n", "[:, :, 1, 2] =\n", " 0.0 1.81343\n", " 1.81343 0.0 \n", "\n", "[:, :, 2, 2] =\n", " 1.81343 0.0 \n", " 0.0 1.3811" ] }, "execution_count": 9, "metadata": {}, "output_type": "execute_result" } ], "source": [ "tensor = init_tensor(1.0)" ] }, { "cell_type": "code", "execution_count": 10, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "2" ] }, "execution_count": 10, "metadata": {}, "output_type": "execute_result" } ], "source": [ "length(tensor[:,1,1,1])" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Splitting tensors for the singular value decomposition" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "This will be similar to the introductory Julia File. We pair the indices of the tensor in the way we want to split it, and define a matrix. On this matrix we will apply a singular value decomposition, truncate here at 2 singular values, and define new 3-valent tensors from it." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Translating back just means splitting the combined index again into two separate ones." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "The function should return two 3-valent tensors and a list of singular values." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Normalize the singular values with respect to the largest one (so the largest singular value is always one)." ] }, { "cell_type": "code", "execution_count": 11, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "compute_SVD1 (generic function with 1 method)" ] }, "execution_count": 11, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function compute_SVD1(tensor::Array{Float64,4})\n", " \n", " #length of tensor\n", " \n", " len = length(tensor[:,1,1,1])\n", " \n", " mat = zeros(len^2,len^2)\n", " \n", " index1 = 0\n", " \n", " for i in 1:len, l in 1:len\n", " \n", " index1 += 1\n", " \n", " index2 = 0\n", " \n", " for j in 1:len, k in 1:len\n", " \n", " index2 += 1\n", " \n", " mat[index1,index2] = tensor[i,j,k,l]\n", " \n", " end\n", " \n", " end\n", " \n", " #Compute the SVD\n", " \n", " println(\"Compute SVD\")\n", " \n", " decomp = svd(mat)\n", " \n", " U = decomp.U\n", " S = decomp.S\n", " V = decomp.V' # Gives us V^\\dagger\n", " \n", " println(\"Print SVs\")\n", " \n", " println(S / S[1]) #prints normalized singular values\n", " \n", " #Define the new S1 and S2\n", " \n", " S1 = zeros(len,len,len) # Here I make the choice to keep the same bond dimension\n", " S2 = zeros(len,len,len)\n", " \n", " index1 = 0\n", " \n", " for i in 1:len, l in 1:len\n", " \n", " index1 += 1\n", " \n", " for p in 1:len\n", " \n", " S1[i,l,p] = U[index1,p] * sqrt(S[p] / S[1])\n", " \n", " end\n", " \n", " end\n", " \n", " index2 = 0\n", " \n", " for j in 1:len, k in 1:len\n", " \n", " index2 += 1\n", " \n", " for q in 1:len\n", " \n", " S2[j,k,q] = sqrt(S[q] / S[1]) * V[q,index2]\n", " \n", " end\n", " \n", " end\n", " \n", " return (S1,S2,S / S[1])\n", " \n", "end " ] }, { "cell_type": "code", "execution_count": 12, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Compute SVD\n", "Print SVs\n", "[1.0, 0.964028, 4.69534e-17, 4.17334e-17]\n" ] }, { "data": { "text/plain": [ "([-0.795551 0.0; 0.0 -0.605887]\n", "\n", "[0.0 -0.694272; -0.694272 0.0], [-0.795551 0.0; 0.0 -0.605887]\n", "\n", "[-0.0 -0.694272; -0.694272 0.0], [1.0, 0.964028, 4.69534e-17, 4.17334e-17])" ] }, "execution_count": 12, "metadata": {}, "output_type": "execute_result" } ], "source": [ "(S1test,S2test) = compute_SVD1(tensor)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We need to the same again, but this time for the other splitting." ] }, { "cell_type": "code", "execution_count": 13, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "compute_SVD2 (generic function with 1 method)" ] }, "execution_count": 13, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function compute_SVD2(tensor::Array{Float64,4})\n", " \n", " #length of tensor\n", " \n", " len = length(tensor[:,1,1,1])\n", " \n", " mat = zeros(len^2,len^2)\n", " \n", " index1 = 0\n", " \n", " for i in 1:len, j in 1:len\n", " \n", " index1 += 1\n", " \n", " index2 = 0\n", " \n", " for k in 1:len, l in 1:len\n", " \n", " index2 += 1\n", " \n", " mat[index1,index2] = tensor[i,j,k,l]\n", " \n", " end\n", " \n", " end\n", " \n", " #Compute the SVD\n", " \n", " println(\"Compute SVD\")\n", " \n", " decomp = svd(mat)\n", " \n", " U = decomp.U\n", " S = decomp.S\n", " V = decomp.V' # Gives us V^\\dagger\n", " \n", " println(\"Print SVs\")\n", " \n", " println(S / S[1]) #prints normalized singular values\n", " \n", " #Define the new S1 and S2\n", " \n", " S3 = zeros(len,len,len) # Here I make the choice to keep the same bond dimension\n", " S4 = zeros(len,len,len)\n", " \n", " index1 = 0\n", " \n", " for i in 1:len, j in 1:len\n", " \n", " index1 += 1\n", " \n", " for p in 1:len\n", " \n", " S3[i,j,p] = U[index1,p] * sqrt(S[p] / S[1])\n", " \n", " end\n", " \n", " end\n", " \n", " index2 = 0\n", " \n", " for k in 1:len, l in 1:len\n", " \n", " index2 += 1\n", " \n", " for q in 1:len\n", " \n", " S4[k,l,q] = sqrt(S[q] / S[1]) * V[q,index2]\n", " \n", " end\n", " \n", " end\n", " \n", " return (S3,S4,S / S[1])\n", " \n", "end " ] }, { "cell_type": "code", "execution_count": 14, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Compute SVD\n", "Print SVs\n", "[1.0, 0.964028, 1.4086e-16, 0.0]\n" ] }, { "data": { "text/plain": [ "([-0.795551 0.0; 0.0 -0.605887]\n", "\n", "[0.0 -0.694272; -0.694272 0.0], [-0.795551 0.0; 0.0 -0.605887]\n", "\n", "[-0.0 -0.694272; -0.694272 0.0], [1.0, 0.964028, 1.4086e-16, 0.0])" ] }, "execution_count": 14, "metadata": {}, "output_type": "execute_result" } ], "source": [ "(S3test,S4test) = compute_SVD2(tensor)" ] }, { "cell_type": "code", "execution_count": 15, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "2×2×2 Array{Float64,3}:\n", "[:, :, 1] =\n", " -0.795551 0.0 \n", " 0.0 -0.605887\n", "\n", "[:, :, 2] =\n", " 0.0 -0.694272\n", " -0.694272 0.0 " ] }, "execution_count": 15, "metadata": {}, "output_type": "execute_result" } ], "source": [ "S3test" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "To get the new effective tensor, we now have to contract the fine degrees of freedom." ] }, { "cell_type": "code", "execution_count": 16, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "contract (generic function with 1 method)" ] }, "execution_count": 16, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function contract(S1::Array{Float64,3},S2::Array{Float64,3},S3::Array{Float64,3},S4::Array{Float64,3})\n", " \n", " len_coarse = length(S1[1,1,:])\n", " \n", " len_fine = length(S1[:,1,1])\n", " \n", " new_ten = zeros(len_coarse,len_coarse,len_coarse,len_coarse)\n", " \n", " for p in 1:len_coarse, q in 1:len_coarse, r in 1:len_coarse, s in 1:len_coarse\n", " \n", " for i in 1:len_fine, j in 1:len_fine, k in 1:len_fine, l in 1:len_fine\n", " \n", " new_ten[p,q,r,s] += S1[i,j,p] * S2[l,k,r] * S3[k,j,q] * S4[i,l,s]\n", " \n", " end\n", " \n", " end\n", " \n", " return new_ten\n", " \n", "end" ] }, { "cell_type": "code", "execution_count": 17, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "2×2×2×2 Array{Float64,4}:\n", "[:, :, 1, 1] =\n", " 0.535325 0.0 \n", " 0.0 0.482014\n", "\n", "[:, :, 2, 1] =\n", " 0.0 0.482014\n", " 0.464675 0.0 \n", "\n", "[:, :, 1, 2] =\n", " 0.0 0.464675\n", " 0.482014 0.0 \n", "\n", "[:, :, 2, 2] =\n", " 0.482014 0.0 \n", " 0.0 0.464675" ] }, "execution_count": 17, "metadata": {}, "output_type": "execute_result" } ], "source": [ "contract(S1test,S2test,S3test,S4test)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "To compute the flow of renormalized tensors, we have to iterate these procedures a few times. We write this into a main function." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We have to specify the number of iterations of the algorithm, then define the initial tensor. We can pass $\\beta$ as a parameter of the main function." ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Then we do a for loop, in which we do the first splitting, second splitting and then contract the new three valent tensor. Then we repeat.\n", "\n", "For each iteration we should store the singular values computed in the singular value decomposition. We will afterwards plot them to distinguish the phases of the Ising model and finding the phase transition." ] }, { "cell_type": "code", "execution_count": 18, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "main (generic function with 1 method)" ] }, "execution_count": 18, "metadata": {}, "output_type": "execute_result" } ], "source": [ "function main(beta::Float64)\n", " \n", " iter = 30 #number of iterations\n", " \n", " svlist = zeros(iter,4)\n", " \n", " tensor = init_tensor(beta)\n", " \n", " for i in 1:iter\n", " \n", " println(\"Iteration \", i)\n", " \n", " println(\"Compute first SVD\")\n", " \n", " (S1,S2,sv1) = compute_SVD1(tensor);\n", " \n", " println(\"Compute second SVD\")\n", " \n", " (S3,S4,sv2) = compute_SVD2(tensor);\n", " \n", " svlist[i,1] = i\n", " svlist[i,2] = sv1[2]\n", " svlist[i,3] = sv1[3]\n", " svlist[i,4] = sv1[4]\n", " \n", " println(\"Compute new tensor\")\n", " \n", " tensor = contract(S1,S2,S3,S4)\n", " \n", " println()\n", " \n", " end\n", " \n", " return svlist\n", " \n", "end" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "If it is running, try the algorithm for different initial values of $\\beta$ and look for a phase transition. The actual critical inverse temperature of the 2D square lattice Ising model is $\\beta_c \\approx 0.4406...$. The plot the singular values for different betas." ] }, { "cell_type": "code", "execution_count": 19, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Iteration 1\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.53705, 3.31113e-17, 2.24788e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.53705, 6.62226e-17, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 2\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.436412, 0.131445, 0.0573641]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.436412, 0.131445, 0.0573641]\n", "Compute new tensor\n", "\n", "Iteration 3\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.397664, 0.0986671, 0.0392363]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.397664, 0.0986671, 0.0392363]\n", "Compute new tensor\n", "\n", "Iteration 4\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.329952, 0.105168, 0.0347005]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.329952, 0.105168, 0.0347005]\n", "Compute new tensor\n", "\n", "Iteration 5\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.257491, 0.0974826, 0.0251009]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.257491, 0.0974826, 0.0251009]\n", "Compute new tensor\n", "\n", "Iteration 6\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.184526, 0.0860175, 0.0158725]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.184526, 0.0860175, 0.0158725]\n", "Compute new tensor\n", "\n", "Iteration 7\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.121123, 0.0688877, 0.0083439]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.121123, 0.0688877, 0.0083439]\n", "Compute new tensor\n", "\n", "Iteration 8\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0732355, 0.0496831, 0.00363856]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0732355, 0.0496831, 0.00363856]\n", "Compute new tensor\n", "\n", "Iteration 9\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0413942, 0.0323024, 0.00133713]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0413942, 0.0323024, 0.00133713]\n", "Compute new tensor\n", "\n", "Iteration 10\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0222755, 0.0192138, 0.000427997]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0222755, 0.0192138, 0.000427997]\n", "Compute new tensor\n", "\n", "Iteration 11\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0116086, 0.0106833, 0.000124019]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0116086, 0.0106833, 0.000124019]\n", "Compute new tensor\n", "\n", "Iteration 12\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.00593499, 0.00567612, 3.36877e-5]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.00593499, 0.00567612, 3.36877e-5]\n", "Compute new tensor\n", "\n", "Iteration 13\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.00300213, 0.00293321, 8.80586e-6]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.00300213, 0.00293321, 8.80586e-6]\n", "Compute new tensor\n", "\n", "Iteration 14\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.00151, 0.00149218, 2.25318e-6]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.00151, 0.00149218, 2.25318e-6]\n", "Compute new tensor\n", "\n", "Iteration 15\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000757267, 0.000752734, 5.70021e-7]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000757267, 0.000752734, 5.70021e-7]\n", "Compute new tensor\n", "\n", "Iteration 16\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000379205, 0.000378054, 1.4336e-7]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000380584, 0.000376684, 1.4336e-7]\n", "Compute new tensor\n", "\n", "Iteration 17\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000279528, 3.76216e-5, 1.05163e-8]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000251024, 4.18936e-5, 1.05163e-8]\n", "Compute new tensor\n", "\n", "Iteration 18\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.000264893, 1.07738e-10, 2.8539e-14]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.8571e-7, 1.5332e-7, 2.8473e-14]\n", "Compute new tensor\n", "\n", "Iteration 19\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 2.42715e-6, 2.61246e-10, 6.34083e-16]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 6.58044e-6, 9.63587e-11, 6.34083e-16]\n", "Compute new tensor\n", "\n", "Iteration 20\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 4.8989e-14, 1.80498e-14, 8.8424e-28]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 3.99646e-6, 2.06996e-22, 8.27253e-28]\n", "Compute new tensor\n", "\n", "Iteration 21\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 4.4238e-10, 1.43908e-21, 6.36622e-31]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 9.07017e-12, 7.01885e-20, 6.36622e-31]\n", "Compute new tensor\n", "\n", "Iteration 22\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 6.3344e-11, 1.75162e-46, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 3.95584e-27, 9.94559e-29, 3.93346e-55]\n", "Compute new tensor\n", "\n", "Iteration 23\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 3.33452e-34, 3.32983e-35, 1.11034e-68]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 5.00578e-19, 2.17134e-84, 6.50547e-92]\n", "Compute new tensor\n", "\n", "Iteration 24\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 2.23092e-45, 1.38935e-70, 4.06071e-115]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.29197e-26, 6.78195e-123, -0.0]\n", "Compute new tensor\n", "\n", "Iteration 25\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 5.36869e-36, 0.0, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 4.14179e-74, 0.0, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 26\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.17017e-161, 1.75345e-192, 5.18065e-318]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 4.7155e-55, 4.27065e-300, -0.0]\n", "Compute new tensor\n", "\n", "Iteration 27\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 2.34903e-108, 0.0, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 7.67845e-239, 0.0, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 28\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.34302e-173, 0.0, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0, 0.0, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 29\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0, 0.0, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0, 0.0, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 30\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0, 0.0, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.0, 0.0, 0.0]\n", "Compute new tensor\n", "\n" ] } ], "source": [ "res0_3 = main(0.3);" ] }, { "cell_type": "code", "execution_count": 20, "metadata": {}, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "Iteration 1\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.885352, 0.0, -0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.885352, 0.0, -0.0]\n", "Compute new tensor\n", "\n", "Iteration 2\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.873306, 0.0531058, 0.0463776]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.873306, 0.0531058, 0.0463776]\n", "Compute new tensor\n", "\n", "Iteration 3\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.959152, 0.013882, 0.013315]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.959152, 0.013882, 0.013315]\n", "Compute new tensor\n", "\n", "Iteration 4\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.984828, 0.00685986, 0.00675578]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.984828, 0.00685986, 0.00675578]\n", "Compute new tensor\n", "\n", "Iteration 5\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.998294, 0.000782592, 0.000781257]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.998294, 0.000782592, 0.000781257]\n", "Compute new tensor\n", "\n", "Iteration 6\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.999856, 7.11969e-5, 7.11866e-5]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.999856, 7.11969e-5, 7.11866e-5]\n", "Compute new tensor\n", "\n", "Iteration 7\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.999998, 8.55842e-7, 8.5584e-7]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 0.999998, 8.55842e-7, 8.5584e-7]\n", "Compute new tensor\n", "\n", "Iteration 8\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 5.31568e-9, 5.31568e-9]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 5.31568e-9, 5.31568e-9]\n", "Compute new tensor\n", "\n", "Iteration 9\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.50925e-13, 7.50755e-13]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.50833e-13, 7.50792e-13]\n", "Compute new tensor\n", "\n", "Iteration 10\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 2.48636e-25, -0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.76642e-17, 3.64438e-17]\n", "Compute new tensor\n", "\n", "Iteration 11\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 3.92523e-17, 0.0]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.61841e-16, 4.22777e-33]\n", "Compute new tensor\n", "\n", "Iteration 12\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 9.64879e-17, 3.81728e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.98189e-17, 3.43078e-17]\n", "Compute new tensor\n", "\n", "Iteration 13\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 3.9376e-17, 1.33628e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 4.80741e-17, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 14\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.6356e-16, 4.90996e-18]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.57567e-32, 1.70494e-48]\n", "Compute new tensor\n", "\n", "Iteration 15\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.98139e-17, 5.80547e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.20921e-16, 1.08049e-17]\n", "Compute new tensor\n", "\n", "Iteration 16\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.50207e-17, 2.61026e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.65461e-16, 4.31663e-17]\n", "Compute new tensor\n", "\n", "Iteration 17\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.03644e-16, 4.42703e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.76288e-16, 5.92268e-17]\n", "Compute new tensor\n", "\n", "Iteration 18\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 6.03741e-17, 2.36823e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.15868e-17, 2.59492e-17]\n", "Compute new tensor\n", "\n", "Iteration 19\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.85046e-17, 7.85046e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 4.85826e-17, 3.08991e-32]\n", "Compute new tensor\n", "\n", "Iteration 20\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.07739e-16, 3.89693e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.80614e-16, 2.41289e-17]\n", "Compute new tensor\n", "\n", "Iteration 21\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.45289e-17, 3.64548e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 3.92523e-17, 5.44735e-34]\n", "Compute new tensor\n", "\n", "Iteration 22\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.83392e-16, 2.60528e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.85046e-17, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 23\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.89283e-17, 4.7826e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.17784e-16, 1.35831e-17]\n", "Compute new tensor\n", "\n", "Iteration 24\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.88106e-17, 2.52175e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.05278e-16, 5.28191e-17]\n", "Compute new tensor\n", "\n", "Iteration 25\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 4.70799e-17, 3.27261e-18]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.85046e-17, 0.0]\n", "Compute new tensor\n", "\n", "Iteration 26\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 5.55746e-17, 3.91611e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.15882e-17, 3.91219e-17]\n", "Compute new tensor\n", "\n", "Iteration 27\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.07025e-17, 2.34263e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 7.49006e-17, 5.50811e-18]\n", "Compute new tensor\n", "\n", "Iteration 28\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.94947e-17, 2.43618e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 1.31494e-16, 3.25866e-17]\n", "Compute new tensor\n", "\n", "Iteration 29\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 3.92523e-17, 8.88698e-34]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 2.16763e-31, -0.0]\n", "Compute new tensor\n", "\n", "Iteration 30\n", "Compute first SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 6.18019e-17, 2.3593e-17]\n", "Compute second SVD\n", "Compute SVD\n", "Print SVs\n", "[1.0, 1.0, 8.92162e-17, 1.51768e-17]\n", "Compute new tensor\n", "\n" ] } ], "source": [ "res0_7 = main(0.7);" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Plot the singular values, actually the second largest one." ] }, { "cell_type": "code", "execution_count": 21, "metadata": {}, "outputs": [ { "data": { "text/plain": [ "2-element Array{Array{Float64,1},1}:\n", " [0.53705, 0.436412, 0.397664, 0.329952, 0.257491, 0.184526, 0.121123, 0.0732355, 0.0413942, 0.0222755 … 4.4238e-10, 6.3344e-11, 3.33452e-34, 2.23092e-45, 5.36869e-36, 1.17017e-161, 2.34903e-108, 1.34302e-173, 0.0, 0.0]\n", " [0.885352, 0.873306, 0.959152, 0.984828, 0.998294, 0.999856, 0.999998, 1.0, 1.0, 1.0 … 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0] " ] }, "execution_count": 21, "metadata": {}, "output_type": "execute_result" } ], "source": [ "y = [res0_3[:,2],res0_7[:,2]]" ] }, { "cell_type": "code", "execution_count": 22, "metadata": {}, "outputs": [ { "data": { "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "5\n", "\n", "\n", "10\n", "\n", "\n", "15\n", "\n", "\n", "20\n", "\n", "\n", "25\n", "\n", "\n", "30\n", "\n", "\n", "0.00\n", "\n", "\n", "0.25\n", "\n", "\n", "0.50\n", "\n", "\n", "0.75\n", "\n", "\n", "1.00\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "y1\n", "\n", "\n", "\n", "y2\n", "\n", "\n" ] }, "execution_count": 22, "metadata": {}, "output_type": "execute_result" } ], "source": [ "plot(res0_3[:,1],y)" ] }, { "cell_type": "code", "execution_count": 24, "metadata": {}, "outputs": [ { "data": { "image/svg+xml": [ "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", " \n", " \n", " \n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "5\n", "\n", "\n", "10\n", "\n", "\n", "15\n", "\n", "\n", "20\n", "\n", "\n", "25\n", "\n", "\n", "30\n", "\n", "\n", "0.00\n", "\n", "\n", "0.25\n", "\n", "\n", "0.50\n", "\n", "\n", "0.75\n", "\n", "\n", "1.00\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "\n", "y1\n", "\n", "\n", "\n", "y2\n", "\n", "\n", "\n", "y3\n", "\n", "\n", "\n", "y4\n", "\n", "\n", "\n", "y5\n", "\n", "\n", "\n", "y6\n", "\n", "\n" ] }, "execution_count": 24, "metadata": {}, "output_type": "execute_result" } ], "source": [ "plot(res0_3[:,1],[res0_3[:,2],res0_3[:,3],res0_3[:,4],res0_7[:,2],res0_7[:,3],res0_7[:,4]])" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Julia 1.0.1", "language": "julia", "name": "julia-1.0" }, "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.0.1" } }, "nbformat": 4, "nbformat_minor": 2 }