{ "cells": [ { "cell_type": "code", "execution_count": 6, "metadata": { "scrolled": false }, "outputs": [ { "name": "stdout", "output_type": "stream", "text": [ "H mult 20 E estimates: -4.8736033 -4.870925\n", "H mult 40 E estimates: -4.906122 -4.936535\n", "H mult 60 E estimates: -4.953764 -4.958339\n", "H mult 80 E estimates: -4.954471 -4.969118\n", "H mult 100 E estimates: -4.9802227 -4.975509\n", "H mult 120 E estimates: -4.978794 -4.979719\n", "H mult 140 E estimates: -4.9785876 -4.9826875\n", "H mult 160 E estimates: -4.979694 -4.9848804\n", "H mult 180 E estimates: -4.9799604 -4.986557\n", "H mult 200 E estimates: -4.9803667 -4.987875\n", "H mult 220 E estimates: -4.9814525 -4.9889355\n", "H mult 240 E estimates: -4.981823 -4.989806\n", "H mult 260 E estimates: -4.9822025 -4.990535\n", "H mult 280 E estimates: -4.990591 -4.9911537\n", "H mult 300 E estimates: -4.9908123 -4.9916887\n", "H mult 320 E estimates: -4.991027 -4.9921565\n", "H mult 340 E estimates: -4.9912343 -4.992571\n", "H mult 360 E estimates: -4.991433 -4.9929423\n", "H mult 380 E estimates: -4.9916234 -4.993278\n", "H mult 400 E estimates: -4.991859 -4.9935837\n", "H mult 420 E estimates: -4.9920354 -4.9938645\n", "H mult 440 E estimates: -4.9922 -4.9941235\n", "H mult 460 E estimates: -4.9923544 -4.994364\n", "H mult 480 E estimates: -4.992499 -4.994588\n", "H mult 500 E estimates: -4.992635 -4.9947977\n", "H mult 520 E estimates: -4.9927626 -4.9949946\n", "H mult 540 E estimates: -4.9949656 -4.9951797\n", "H mult 560 E estimates: -4.9952087 -4.995354\n", "H mult 580 E estimates: -4.99534 -4.995519\n", "H mult 600 E estimates: -4.9954615 -4.995675\n", "H mult 620 E estimates: -4.9956956 -4.995823\n", "H mult 640 E estimates: -4.9957986 -4.9959626\n", "H mult 660 E estimates: -4.996016 -4.9960957\n", "H mult 680 E estimates: -4.996103 -4.9962215\n", "H mult 700 E estimates: -4.996185 -4.996341\n", "H mult 720 E estimates: -4.9963756 -4.9964547\n", "H mult 740 E estimates: -4.996445 -4.996563\n", "H mult 760 E estimates: -4.9966187 -4.9966655\n", "H mult 780 E estimates: -4.9967337 -4.996763\n", "H mult 800 E estimates: -4.9968333 -4.9968557\n", "H mult 820 E estimates: -4.9969277 -4.996944\n", "H mult 840 E estimates: -4.9969716 -4.997028\n", "H mult 860 E estimates: -4.997099 -4.9971075\n", "H mult 880 E estimates: -4.997136 -4.9971833\n", "H mult 900 E estimates: -4.9972486 -4.997256\n", "H mult 920 E estimates: -4.9972806 -4.9973245\n", "H mult 940 E estimates: -4.9973793 -4.99739\n", "H mult 960 E estimates: -4.9974065 -4.997452\n", "H mult 980 E estimates: -4.997482 -4.997511\n", "H mult 1000 E estimates: -4.997507 -4.997567\n", "H mult 1020 E estimates: -4.997586 -4.9976206\n", "H mult 1040 E estimates: -4.997607 -4.9976716\n", "H mult 1060 E estimates: -4.9976287 -4.9977202\n", "H mult 1080 E estimates: -4.997694 -4.9977665\n", "H mult 1100 E estimates: -4.997712 -4.9978104\n", "H mult 1120 E estimates: -4.9977303 -4.9978523\n", "H mult 1140 E estimates: -4.9977484 -4.9978924\n", "H mult 1160 E estimates: -4.9978004 -4.9979305\n", "H mult 1180 E estimates: -4.997816 -4.997967\n", "H mult 1200 E estimates: -4.997832 -4.998001\n", "H mult 1220 E estimates: -4.9978476 -4.998034\n", "H mult 1240 E estimates: -4.997863 -4.9980655\n", "H mult 1260 E estimates: -4.9979153 -4.9980955\n", "H mult 1280 E estimates: -4.997929 -4.998124\n", "H mult 1300 E estimates: -4.997943 -4.998152\n", "H mult 1320 E estimates: -4.9979563 -4.998178\n", "H mult 1340 E estimates: -4.9979897 -4.998203\n", "H mult 1360 E estimates: -4.998002 -4.9982266\n", "H mult 1380 E estimates: -4.9980145 -4.9982495\n", "H mult 1400 E estimates: -4.9980264 -4.9982715\n", "H mult 1420 E estimates: -4.9980383 -4.9982924\n", "H mult 1440 E estimates: -4.99805 -4.9983125\n", "H mult 1460 E estimates: -4.9980617 -4.9983315\n", "H mult 1480 E estimates: -4.9981017 -4.99835\n", "H mult 1500 E estimates: -4.998112 -4.998368\n", "H mult 1520 E estimates: -4.9981227 -4.9983845\n", "H mult 1540 E estimates: -4.9981327 -4.9984007\n", "H mult 1560 E estimates: -4.9981427 -4.9984164\n", "H mult 1580 E estimates: -4.9981527 -4.9984317\n", "H mult 1600 E estimates: -4.9981627 -4.998446\n", "H mult 1620 E estimates: -4.9981723 -4.99846\n", "H mult 1640 E estimates: -4.998182 -4.998473\n", "H mult 1660 E estimates: -4.998192 -4.998486\n", "H mult 1680 E estimates: -4.9982095 -4.9984984\n", "H mult 1700 E estimates: -4.9982185 -4.9985104\n", "H mult 1720 E estimates: -4.9982276 -4.998522\n", "H mult 1740 E estimates: -4.998248 -4.998533\n", "H mult 1760 E estimates: -4.998256 -4.9985437\n", "H mult 1780 E estimates: -4.998265 -4.998554\n", "H mult 1800 E estimates: -4.998279 -4.9985642\n", "H mult 1820 E estimates: -4.998287 -4.998574\n", "H mult 1840 E estimates: -4.998295 -4.9985833\n", "H mult 1860 E estimates: -4.9983025 -4.9985924\n", "H mult 1880 E estimates: -4.99831 -4.998601\n", "H mult 1900 E estimates: -4.9983177 -4.9986095\n", "H mult 1920 E estimates: -4.9983253 -4.9986176\n", "H mult 1940 E estimates: -4.998333 -4.9986258\n", "H mult 1960 E estimates: -4.99834 -4.9986334\n", "H mult 1980 E estimates: -4.9983478 -4.998641\n", "H mult 2000 E estimates: -4.998355 -4.998648\n", "H mult 2020 E estimates: -4.998362 -4.9986553\n", "H mult 2040 E estimates: -4.998369 -4.998662\n", "H mult 2060 E estimates: -4.9983764 -4.9986687\n", "H mult 2080 E estimates: -4.998383 -4.9986753\n", "H mult 2100 E estimates: -4.9983997 -4.9986815\n", "H mult 2120 E estimates: -4.998406 -4.9986877\n", "H mult 2140 E estimates: -4.998412 -4.998694\n", "H mult 2160 E estimates: -4.998419 -4.9986997\n", "H mult 2180 E estimates: -4.998425 -4.9987054\n", "H mult 2200 E estimates: -4.998431 -4.998711\n", "H mult 2220 E estimates: -4.9984374 -4.9987164\n", "H mult 2240 E estimates: -4.9984436 -4.9987216\n", "H mult 2260 E estimates: -4.99845 -4.998727\n", "H mult 2280 E estimates: -4.998456 -4.998732\n", "H mult 2300 E estimates: -4.9984617 -4.998737\n", "H mult 2320 E estimates: -4.998468 -4.9987416\n", "H mult 2340 E estimates: -4.9984736 -4.9987464\n", "H mult 2360 E estimates: -4.99848 -4.998751\n", "H mult 2380 E estimates: -4.9984856 -4.9987555\n", "H mult 2400 E estimates: -4.9984913 -4.99876\n", "H mult 2420 E estimates: -4.998497 -4.9987645\n", "H mult 2440 E estimates: -4.9985027 -4.998769\n", "H mult 2460 E estimates: -4.9985085 -4.9987726\n", "H mult 2480 E estimates: -4.998514 -4.998777\n", "H mult 2500 E estimates: -4.99852 -4.9987807\n", "H mult 2520 E estimates: -4.998525 -4.9987845\n", "H mult 2540 E estimates: -4.998531 -4.9987884\n", "H mult 2560 E estimates: -4.998536 -4.998792\n", "H mult 2580 E estimates: -4.998542 -4.998796\n", "H mult 2600 E estimates: -4.998547 -4.9987993\n", "H mult 2620 E estimates: -4.9985523 -4.998803\n", "H mult 2640 E estimates: -4.9985576 -4.9988065\n", "H mult 2660 E estimates: -4.998563 -4.99881\n", "H mult 2680 E estimates: -4.998568 -4.998813\n", "H mult 2700 E estimates: -4.9985733 -4.9988165\n", "H mult 2720 E estimates: -4.9985785 -4.99882\n", "H mult 2740 E estimates: -4.998584 -4.998823\n", "H mult 2760 E estimates: -4.9985886 -4.998826\n", "H mult 2780 E estimates: -4.998598 -4.998829\n", "H mult 2800 E estimates: -4.9986033 -4.998832\n", "H mult 2820 E estimates: -4.998608 -4.998835\n", "H mult 2840 E estimates: -4.998613 -4.998838\n", "H mult 2860 E estimates: -4.998618 -4.998841\n", "H mult 2880 E estimates: -4.998623 -4.9988437\n", "H mult 2900 E estimates: -4.9986277 -4.9988465\n", "H mult 2920 E estimates: -4.9986324 -4.9988494\n", "H mult 2940 E estimates: -4.998637 -4.998852\n", "H mult 2960 E estimates: -4.9986415 -4.9988546\n", "H mult 2980 E estimates: -4.9986463 -4.998857\n", "H mult 3000 E estimates: -4.998651 -4.99886\n", "H mult 3020 E estimates: -4.9986553 -4.9988623\n", "H mult 3040 E estimates: -4.99866 -4.9988647\n", "H mult 3060 E estimates: -4.9986644 -4.998867\n", "H mult 3080 E estimates: -4.9986687 -4.9988694\n", "H mult 3100 E estimates: -4.9986734 -4.998872\n", "H mult 3120 E estimates: -4.9986777 -4.998874\n", "H mult 3140 E estimates: -4.998682 -4.9988766\n", "H mult 3160 E estimates: -4.9986863 -4.998879\n", "H mult 3180 E estimates: -4.9986906 -4.9988813\n", "H mult 3200 E estimates: -4.998695 -4.9988832\n", "H mult 3220 E estimates: -4.998699 -4.9988856\n", "H mult 3240 E estimates: -4.998703 -4.9988875\n", "H mult 3260 E estimates: -4.9987073 -4.99889\n", "H mult 3280 E estimates: -4.9987116 -4.998892\n", "H mult 3300 E estimates: -4.9987154 -4.998894\n", "H mult 3320 E estimates: -4.998719 -4.998896\n", "H mult 3340 E estimates: -4.9987235 -4.998898\n", "H mult 3360 E estimates: -4.9987273 -4.9989004\n", "H mult 3380 E estimates: -4.998731 -4.9989023\n", "H mult 3400 E estimates: -4.998735 -4.998904\n", "H mult 3420 E estimates: -4.9987392 -4.998906\n", "H mult 3440 E estimates: -4.9987483 -4.998908\n", "H mult 3460 E estimates: -4.998752 -4.99891\n", "H mult 3480 E estimates: -4.998756 -4.998912\n", "H mult 3500 E estimates: -4.9987597 -4.998914\n", "H mult 3520 E estimates: -4.998763 -4.998915\n", "H mult 3540 E estimates: -4.998767 -4.998917\n" ] }, { "name": "stdout", "output_type": "stream", "text": [ "H mult 3560 E estimates: -4.9987707 -4.998919\n", "H mult 3580 E estimates: -4.998774 -4.998921\n", "H mult 3600 E estimates: -4.998778 -4.9989223\n", "H mult 3620 E estimates: -4.998781 -4.9989243\n", "H mult 3640 E estimates: -4.9987845 -4.9989257\n", "H mult 3660 E estimates: -4.9987884 -4.9989276\n", "H mult 3680 E estimates: -4.9987917 -4.998929\n", "H mult 3700 E estimates: -4.998795 -4.998931\n", "H mult 3720 E estimates: -4.9987984 -4.9989324\n", "H mult 3740 E estimates: -4.9988017 -4.9989343\n", "H mult 3760 E estimates: -4.998805 -4.9989357\n", "H mult 3780 E estimates: -4.9988084 -4.998937\n", "H mult 3800 E estimates: -4.9988117 -4.998939\n", "H mult 3820 E estimates: -4.998815 -4.9989405\n", "H mult 3840 E estimates: -4.9988184 -4.998942\n", "H mult 3860 E estimates: -4.9988213 -4.9989433\n", "H mult 3880 E estimates: -4.9988246 -4.9989448\n", "H mult 3900 E estimates: -4.998834 -4.9989467\n", "H mult 3920 E estimates: -4.998837 -4.998948\n", "H mult 3940 E estimates: -4.9988403 -4.9989495\n", "H mult 3960 E estimates: -4.998843 -4.998951\n", "H mult 3980 E estimates: -4.998846 -4.9989524\n", "H mult 4000 E estimates: -4.998849 -4.998954\n", "H mult 4020 E estimates: -4.9988523 -4.9989552\n", "H mult 4040 E estimates: -4.998855 -4.9989567\n", "H mult 4060 E estimates: -4.998858 -4.998958\n", "H mult 4080 E estimates: -4.998861 -4.9989595\n", "H mult 4100 E estimates: -4.9988637 -4.9989605\n", "H mult 4120 E estimates: -4.9988666 -4.998962\n", "H mult 4140 E estimates: -4.998869 -4.9989634\n", "H mult 4160 E estimates: -4.998872 -4.998965\n", "H mult 4180 E estimates: -4.9988747 -4.998966\n", "H mult 4200 E estimates: -4.9988775 -4.998967\n", "H mult 4220 E estimates: -4.99888 -4.9989686\n", "H mult 4240 E estimates: -4.998883 -4.99897\n", "H mult 4260 E estimates: -4.998885 -4.9989715\n", "H mult 4280 E estimates: -4.998888 -4.9989724\n", "H mult 4300 E estimates: -4.9988966 -4.998974\n", "H mult 4320 E estimates: -4.9988995 -4.9989753\n", "H mult 4340 E estimates: -4.998902 -4.998976\n", "H mult 4360 E estimates: -4.998904 -4.9989777\n", "H mult 4380 E estimates: -4.998907 -4.9989786\n", "H mult 4400 E estimates: -4.9989095 -4.99898\n", "H mult 4420 E estimates: -4.998912 -4.998981\n", "H mult 4440 E estimates: -4.9989142 -4.9989824\n", "H mult 4460 E estimates: -4.9989166 -4.9989834\n", "H mult 4480 E estimates: -4.998919 -4.998985\n", "H mult 4500 E estimates: -4.9989214 -4.998986\n", "H mult 4520 E estimates: -4.998924 -4.998987\n", "H mult 4540 E estimates: -4.998926 -4.998988\n", "H mult 4560 E estimates: -4.9989285 -4.9989896\n", "H mult 4580 E estimates: -4.9989305 -4.9989905\n", "H mult 4600 E estimates: -4.998933 -4.998992\n", "H mult 4620 E estimates: -4.998935 -4.998993\n", "H mult 4640 E estimates: -4.9989376 -4.998994\n", "H mult 4660 E estimates: -4.9989395 -4.9989953\n", "H mult 4680 E estimates: -4.998942 -4.9989963\n", "H mult 4700 E estimates: -4.9989505 -4.998997\n", "H mult 4720 E estimates: -4.9989524 -4.9989986\n", "H mult 4740 E estimates: -4.998955 -4.9989996\n", "H mult 4760 E estimates: -4.9989567 -4.9990005\n", "H mult 4780 E estimates: -4.998959 -4.999002\n", "H mult 4800 E estimates: -4.998961 -4.999003\n", "H mult 4820 E estimates: -4.998963 -4.999004\n", "H mult 4840 E estimates: -4.998965 -4.999005\n", "H mult 4860 E estimates: -4.998967 -4.9990063\n", "H mult 4880 E estimates: -4.998969 -4.999007\n", "H mult 4900 E estimates: -4.998971 -4.999008\n", "H mult 4920 E estimates: -4.998973 -4.999009\n", "H mult 4940 E estimates: -4.998975 -4.99901\n", "H mult 4960 E estimates: -4.9989767 -4.9990115\n", "H mult 4980 E estimates: -4.9989786 -4.9990125\n", "H mult 5000 E estimates: -4.9989805 -4.9990134\n", "\n", "Lanczos iteration 10 E0 = -4.8897123\n", "Lanczos iteration 20 E0 = -4.972412\n", "Lanczos iteration 30 E0 = -4.987833\n", "Lanczos iteration 40 E0 = -4.9926105\n", "Lanczos iteration 50 E0 = -4.9952383\n", "Lanczos iteration 60 E0 = -4.996963\n", "Lanczos iteration 70 E0 = -4.9979033\n", "Lanczos iteration 80 E0 = -4.9983196\n", "Lanczos iteration 90 E0 = -4.998532\n", "Lanczos iteration 100 E0 = -4.998662\n", "Lanczos iteration 110 E0 = -4.9987807\n", "Lanczos iteration 120 E0 = -4.9988737\n", "Lanczos iteration 130 E0 = -4.9989295\n", "Lanczos iteration 140 E0 = -4.998964\n", "Lanczos iteration 150 E0 = -4.9989924\n", "Lanczos iteration 160 E0 = -4.9990287\n", "Lanczos iteration 170 E0 = -4.9990716\n", "Lanczos iteration 180 E0 = -4.999118\n", "Lanczos iteration 190 E0 = -4.9991536\n", "Lanczos iteration 200 E0 = -4.9991817\n", "\n", "Graphing final psi^2 (times ll^2) from gsprojection\n" ] }, { "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", "\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", "\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" ] }, "execution_count": 6, "metadata": {}, "output_type": "execute_result" } ], "source": [ "# Demonstration of the H-power (ground state projection) and\n", "# Lanczos methods, using the 2D Anderson Localization model\n", "# - a particle on an L*L lattice with random potentials.\n", "\n", " using Plots\n", " using LinearAlgebra\n", "\n", "# Constructs random potential in a vector of size nn=ll*ll\n", "# - v0 is the overall scale of uniformly distributed values\n", "# - v0 should be positive to ensure ground state has the\n", "# lowest |E|. A negative offset (-1) also helps here.\n", "\n", " function randompotential(nn::Int,v0::Float64)\n", " rpot=Vector{Float64}(undef,nn)\n", " for i=1:nn\n", " rpot[i]=v0*(rand()-1.)-1.\n", " end\n", " return rpot\n", " end\n", "\n", "# Creates a random initial state vector\n", "\n", " function randomstate(nn::Int)\n", " rvec=Vector{Float64}(undef,nn)\n", " for i=1:nn\n", " rvec[i]=rand()-0.5\n", " end\n", " rvec.=rvec./dot(rvec,rvec)^0.5\n", " return rvec\n", " end\n", "\n", "# For a 2D lattice of size ll*ll, operates with the Hamiltonian (including the \n", "# random potential energy stored in the vecor vv) on the vector f1, resulting in f2.\n", "# If nor==1, the state f2 is normalized and two energy estimates are returned.\n", "# - r is the ratio of two large wave function coefficients\n", "# - e is the expectation value \n", "\n", "function hoperation(ll::Int,vv,nor,f1,f2)\n", " f2.=vv.*f1\n", " for y0=0:ll-1\n", " for x0=0:ll-1\n", " j0=1+x0+y0*ll\n", " x1::Int=mod(x0-1,ll)\n", " x2::Int=mod(x0+1,ll)\n", " y1::Int=mod(y0-1,ll)\n", " y2::Int=mod(y0+1,ll)\n", " f2[1+x1+y0*ll] -= f1[j0]\n", " f2[1+x2+y0*ll] -= f1[j0]\n", " f2[1+x0+y1*ll] -= f1[j0]\n", " f2[1+x0+y2*ll] -= f1[j0]\n", " end\n", " end\n", " if nor==1\n", " mv,mi=findmax(abs.(f1)) # Julia base function, value and index of max\n", " r::Float64=f2[mi]/f1[mi] # This ratio should converge to the max (abs) eigenvalue\n", " e::Float64=dot(f1,f2) # This is the expectation value \n", " f2.=f2./dot(f2,f2)^0.5 # Now normalize |f2> = H|f1> \n", " return r,e \n", " else\n", " return nothing\n", " end \n", " end\n", "\n", "# Makes a 2D color-contour plot of the wave function psi\n", "# Returns the plot as p, which can be displayed later\n", "\n", "function plotstate(ll::Int,psi)\n", " x=zeros(Int,ll)\n", " y=zeros(Int,ll)\n", " z=zeros(Float64,ll,ll)\n", " for j=1:ll\n", " x[j]=j\n", " y[j]=j\n", " for i=1:ll\n", " z[i,j]=(ll^2)*psi[i+(j-1)*ll]^2\n", " end\n", " end\n", " p=contour(x,y,z,fill=true)\n", " return p\n", " end\n", "\n", "# Performs nit Lanczos iterations, starting from f0.\n", "# The method with normalized states is used.\n", "# Note that the elements in aa[:] and nn[:] start at 1,\n", "# so are shifted by one step relative to the notation used\n", "# in the lectures. After each 10 iterations, the tridiagonal\n", "# H matrix is constructed and its eigenvalues are calculated\n", "# and shown, to check the convergence.\n", "\n", " function lanczos(nit::Int,ll::Int,vv,f0)\n", " f1=copy(f0)\n", " f2=copy(f1)\n", " aa=zeros(Float64,nit)\n", " nn=zeros(Float64,nit+1)\n", " nn[1]=1\n", " hoperation(ll,vv,0,f0,f1)\n", " aa[1]=dot(f0,f1)\n", " f1.=f1.-aa[1].*f0\n", " nn[2]=dot(f1,f1)^0.5\n", " f1.=f1./nn[2]\n", " for m=2:nit\n", " hoperation(ll,vv,0,f1,f2)\n", " aa[m]=dot(f1,f2)\n", " f2.=f2.-aa[m].*f1-nn[m]*f0\n", " nn[m+1]=dot(f2,f2)^0.5\n", " f2.=f2./nn[m+1]\n", " f0.=f1\n", " f1.=f2\n", " if mod(m,10)==0\n", " mat=zeros(Float64,m,m)\n", " for i=1:m\n", " mat[i,i]=aa[i]\n", " if i>1\n", " mat[i-1,i]=nn[i]\n", " end\n", " if i0)\n", " ll=200 # length of the periodic ll*ll lattice\n", " nn=ll*ll # The number of lattice sites\n", " nit1=5000 # Number of iterations for ground state projection\n", " nit2=200 # Number of iterations in the Lanczos calculation\n", "\n", " vv=randompotential(nn,v0)\n", " psi=randomstate(nn)\n", " f0=copy(psi)\n", "\n", " p=gsprojection(nit1,ll,vv,f0)\n", " println()\n", "\n", " f0.=psi\n", " lanczos(nit2,ll,vv,f0)\n", "\n", " println()\n", " println(\"Graphing final psi^2 (times ll^2) from gsprojection\")\n", " plot(p)\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] }, { "cell_type": "code", "execution_count": null, "metadata": {}, "outputs": [], "source": [] } ], "metadata": { "kernelspec": { "display_name": "Julia 1.6.1", "language": "julia", "name": "julia-1.6" }, "language_info": { "file_extension": ".jl", "mimetype": "application/julia", "name": "julia", "version": "1.6.1" } }, "nbformat": 4, "nbformat_minor": 4 }