1
0
mirror of https://github.com/TREX-CoE/qmc-lttc.git synced 2024-12-01 18:18:41 +01:00
qmc-lttc/QMC.org

1344 lines
37 KiB
Org Mode
Raw Normal View History

2021-01-03 18:45:58 +01:00
#+TITLE: Quantum Monte Carlo
#+AUTHOR: Anthony Scemama, Claudia Filippi
2021-01-11 18:41:36 +01:00
# SETUPFILE: https://fniessen.github.io/org-html-themes/org/theme-readtheorg.setup
# SETUPFILE: https://fniessen.github.io/org-html-themes/org/theme-bigblow.setup
2021-01-03 18:45:58 +01:00
#+STARTUP: latexpreview
2021-01-11 20:54:40 +01:00
#+HTML_HEAD: <link rel="stylesheet" title="Standard" href="worg.css" type="text/css" />
2021-01-11 18:41:36 +01:00
2021-01-03 18:45:58 +01:00
* Introduction
2021-01-07 11:07:18 +01:00
We propose different exercises to understand quantum Monte Carlo (QMC)
methods. In the first section, we propose to compute the energy of a
hydrogen atom using numerical integration. The goal of this section is
to introduce the /local energy/.
Then we introduce the variational Monte Carlo (VMC) method which
computes a statistical estimate of the expectation value of the energy
associated with a given wave function.
Finally, we introduce the diffusion Monte Carlo (DMC) method which
2021-01-11 18:41:36 +01:00
gives the exact energy of the $H_2$ molecule.
2021-01-07 11:07:18 +01:00
Code examples will be given in Python and Fortran. Whatever language
can be chosen.
2021-01-11 18:41:36 +01:00
We consider the stationary solution of the Schrödinger equation, so
the wave functions considered here are real: for an $N$ electron
system where the electrons move in the 3-dimensional space,
$\Psi : \mathbb{R}^{3N} \rightarrow \mathbb{R}$. In addition, $\Psi$
is defined everywhere, continuous and infinitely differentiable.
2021-01-11 20:54:40 +01:00
*Note*
#+begin_important
In Fortran, when you use a double precision constant, don't forget
to put d0 as a suffix (for example 2.0d0), or it will be
interpreted as a single precision value
#+end_important
2021-01-03 18:45:58 +01:00
* Numerical evaluation of the energy
2021-01-07 11:07:18 +01:00
In this section we consider the Hydrogen atom with the following
wave function:
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
\Psi(\mathbf{r}) = \exp(-a |\mathbf{r}|)
$$
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
We will first verify that $\Psi$ is an eigenfunction of the Hamiltonian
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
\hat{H} = \hat{T} + \hat{V} = - \frac{1}{2} \Delta - \frac{1}{|\mathbf{r}|}
$$
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
when $a=1$, by checking that $\hat{H}\Psi(\mathbf{r}) = E\Psi(\mathbf{r})$ for
2021-01-11 18:41:36 +01:00
all $\mathbf{r}$. We will check that the local energy, defined as
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
E_L(\mathbf{r}) = \frac{\hat{H} \Psi(\mathbf{r})}{\Psi(\mathbf{r})},
$$
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
is constant. We will also see that when $a \ne 1$ the local energy
is not constant, so $\hat{H} \Psi \ne E \Psi$.
2021-01-03 18:45:58 +01:00
2021-01-11 18:41:36 +01:00
The probabilistic /expected value/ of an arbitrary function $f(x)$
with respect to a probability density function $p(x)$ is given by
$$ \langle f \rangle_p = \int_{-\infty}^\infty p(x)\, f(x)\,dx $$.
Recall that a probability density function $p(x)$ is non-negative
and integrates to one:
$$ \int_{-\infty}^\infty p(x)\,dx = 1 $$.
The electronic energy of a system is the expectation value of the
2021-01-11 20:54:40 +01:00
local energy $E(\mathbf{r})$ with respect to the 3N-dimensional
2021-01-11 18:41:36 +01:00
electron density given by the square of the wave function:
2021-01-11 20:54:40 +01:00
\begin{eqnarray*}
E & = & \frac{\langle \Psi| \hat{H} | \Psi\rangle}{\langle \Psi |\Psi \rangle}
= \frac{\int \Psi(\mathbf{r})\, \hat{H} \Psi(\mathbf{r})\, d\mathbf{r}}{\int \left[\Psi(\mathbf{r}) \right]^2 d\mathbf{r}} \\
& = & \frac{\int \left[\Psi(\mathbf{r})\right]^2\, \frac{\hat{H} \Psi(\mathbf{r})}{\Psi(\mathbf{r})}\,d\mathbf{r}}{\int \left[\Psi(\mathbf{r}) \right]^2 d\mathbf{r}}
= \frac{\int \left[\Psi(\mathbf{r})\right]^2\, E_L(\mathbf{r})\,d\mathbf{r}}{\int \left[\Psi(\mathbf{r}) \right]^2 d\mathbf{r}}
2021-01-11 18:41:36 +01:00
= \langle E_L \rangle_{\Psi^2}
2021-01-11 20:54:40 +01:00
\end{eqnarray*}
2021-01-11 18:41:36 +01:00
2021-01-03 18:45:58 +01:00
** Local energy
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle hydrogen.py
:header-args:f90: :tangle hydrogen.f90
:END:
2021-01-12 01:01:52 +01:00
*** Exercise 1
#+begin_exercise
Write a function which computes the potential at $\mathbf{r}$.
2021-01-07 11:07:18 +01:00
The function accepts a 3-dimensional vector =r= as input arguments
and returns the potential.
2021-01-12 01:01:52 +01:00
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
$\mathbf{r}=\left( \begin{array}{c} x \\ y\\ z\end{array} \right)$, so
2021-01-07 11:07:18 +01:00
$$
2021-01-11 20:54:40 +01:00
V(\mathbf{r}) = -\frac{1}{\sqrt{x^2 + y^2 + z^2}}
2021-01-07 11:07:18 +01:00
$$
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Python*
#+BEGIN_SRC python :results none
2021-01-03 18:45:58 +01:00
import numpy as np
def potential(r):
return -1. / np.sqrt(np.dot(r,r))
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-07 10:01:55 +01:00
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Fortran*
#+BEGIN_SRC f90
2021-01-03 18:45:58 +01:00
double precision function potential(r)
implicit none
double precision, intent(in) :: r(3)
potential = -1.d0 / dsqrt( r(1)*r(1) + r(2)*r(2) + r(3)*r(3) )
end function potential
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
*** Exercise 2
#+begin_exercise
Write a function which computes the wave function at $\mathbf{r}$.
2021-01-07 11:07:18 +01:00
The function accepts a scalar =a= and a 3-dimensional vector =r= as
input arguments, and returns a scalar.
2021-01-12 01:01:52 +01:00
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Python*
#+BEGIN_SRC python :results none
2021-01-03 18:45:58 +01:00
def psi(a, r):
return np.exp(-a*np.sqrt(np.dot(r,r)))
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Fortran*
#+BEGIN_SRC f90
2021-01-03 18:45:58 +01:00
double precision function psi(a, r)
implicit none
double precision, intent(in) :: a, r(3)
psi = dexp(-a * dsqrt( r(1)*r(1) + r(2)*r(2) + r(3)*r(3) ))
end function psi
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
*** Exercise 3
#+begin_exercise
Write a function which computes the local kinetic energy at $\mathbf{r}$.
2021-01-07 11:07:18 +01:00
The function accepts =a= and =r= as input arguments and returns the
local kinetic energy.
2021-01-12 01:01:52 +01:00
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
The local kinetic energy is defined as $$-\frac{1}{2}\frac{\Delta \Psi}{\Psi}.$$
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
We differentiate $\Psi$ with respect to $x$:
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
\[\Psi(\mathbf{r}) = \exp(-a\,|\mathbf{r}|) \]
\[\frac{\partial \Psi}{\partial x}
= \frac{\partial \Psi}{\partial |\mathbf{r}|} \frac{\partial |\mathbf{r}|}{\partial x}
= - \frac{a\,x}{|\mathbf{r}|} \Psi(\mathbf{r}) \]
2021-01-07 11:07:18 +01:00
and we differentiate a second time:
$$
\frac{\partial^2 \Psi}{\partial x^2} =
2021-01-11 20:54:40 +01:00
\left( \frac{a^2\,x^2}{|\mathbf{r}|^2} -
\frac{a(y^2+z^2)}{|\mathbf{r}|^{3}} \right) \Psi(\mathbf{r}).
2021-01-07 11:07:18 +01:00
$$
The Laplacian operator $\Delta = \frac{\partial^2}{\partial x^2} +
\frac{\partial^2}{\partial y^2} + \frac{\partial^2}{\partial z^2}$
applied to the wave function gives:
$$
2021-01-11 20:54:40 +01:00
\Delta \Psi (\mathbf{r}) = \left(a^2 - \frac{2a}{\mathbf{|r|}} \right) \Psi(\mathbf{r})
2021-01-07 11:07:18 +01:00
$$
So the local kinetic energy is
$$
2021-01-11 20:54:40 +01:00
-\frac{1}{2} \frac{\Delta \Psi}{\Psi} (\mathbf{r}) = -\frac{1}{2}\left(a^2 - \frac{2a}{\mathbf{|r|}} \right)
2021-01-07 11:07:18 +01:00
$$
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Python*
#+BEGIN_SRC python :results none
2021-01-03 18:45:58 +01:00
def kinetic(a,r):
return -0.5 * (a**2 - (2.*a)/np.sqrt(np.dot(r,r)))
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Fortran*
#+BEGIN_SRC f90
2021-01-03 18:45:58 +01:00
double precision function kinetic(a,r)
implicit none
double precision, intent(in) :: a, r(3)
kinetic = -0.5d0 * (a*a - (2.d0*a) / &
dsqrt( r(1)*r(1) + r(2)*r(2) + r(3)*r(3) ) )
end function kinetic
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
*** Exercise 4
#+begin_exercise
Write a function which computes the local energy at $\mathbf{r}$.
2021-01-07 11:07:18 +01:00
The function accepts =x,y,z= as input arguments and returns the
local energy.
2021-01-12 01:01:52 +01:00
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
2021-01-11 20:54:40 +01:00
E_L(\mathbf{r}) = -\frac{1}{2} \frac{\Delta \Psi}{\Psi} (\mathbf{r}) + V(\mathbf{r})
2021-01-07 11:07:18 +01:00
$$
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Python*
#+BEGIN_SRC python :results none
2021-01-03 18:45:58 +01:00
def e_loc(a,r):
return kinetic(a,r) + potential(r)
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Fortran*
#+BEGIN_SRC f90
2021-01-03 18:45:58 +01:00
double precision function e_loc(a,r)
implicit none
double precision, intent(in) :: a, r(3)
double precision, external :: kinetic, potential
e_loc = kinetic(a,r) + potential(r)
end function e_loc
2021-01-11 20:54:40 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-11 18:41:36 +01:00
** Plot of the local energy along the $x$ axis
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle plot_hydrogen.py
:header-args:f90: :tangle plot_hydrogen.f90
:END:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
*** Exercise
#+begin_exercise
For multiple values of $a$ (0.1, 0.2, 0.5, 1., 1.5, 2.), plot the
local energy along the $x$ axis.
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
*Python*
#+BEGIN_SRC python :results none
2021-01-03 18:45:58 +01:00
import numpy as np
import matplotlib.pyplot as plt
from hydrogen import e_loc
x=np.linspace(-5,5)
def make_array(a):
y=np.array([ e_loc(a, np.array([t,0.,0.]) ) for t in x])
return y
plt.figure(figsize=(10,5))
for a in [0.1, 0.2, 0.5, 1., 1.5, 2.]:
y = make_array(a)
plt.plot(x,y,label=f"a={a}")
plt.tight_layout()
plt.legend()
plt.savefig("plot_py.png")
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
2021-01-07 10:01:55 +01:00
2021-01-12 01:01:52 +01:00
[[./plot_py.png]]
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
2021-01-12 01:01:52 +01:00
*Fortran*
#+begin_src f90
2021-01-03 18:45:58 +01:00
program plot
implicit none
double precision, external :: e_loc
double precision :: x(50), energy, dx, r(3), a(6)
integer :: i, j
a = (/ 0.1d0, 0.2d0, 0.5d0, 1.d0, 1.5d0, 2.d0 /)
dx = 10.d0/(size(x)-1)
do i=1,size(x)
x(i) = -5.d0 + (i-1)*dx
end do
r(:) = 0.d0
do j=1,size(a)
print *, '# a=', a(j)
do i=1,size(x)
r(1) = x(i)
energy = e_loc( a(j), r )
print *, x(i), energy
end do
print *, ''
print *, ''
end do
end program plot
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
To compile and run:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+begin_src sh :exports both
2021-01-03 18:45:58 +01:00
gfortran hydrogen.f90 plot_hydrogen.f90 -o plot_hydrogen
./plot_hydrogen > data
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
To plot the data using gnuplot:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+begin_src gnuplot :file plot.png :exports both
2021-01-03 18:45:58 +01:00
set grid
set xrange [-5:5]
set yrange [-2:1]
plot './data' index 0 using 1:2 with lines title 'a=0.1', \
'./data' index 1 using 1:2 with lines title 'a=0.2', \
'./data' index 2 using 1:2 with lines title 'a=0.5', \
'./data' index 3 using 1:2 with lines title 'a=1.0', \
'./data' index 4 using 1:2 with lines title 'a=1.5', \
'./data' index 5 using 1:2 with lines title 'a=2.0'
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
[[file:plot.png]]
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
** Numerical estimation of the energy
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle energy_hydrogen.py
:header-args:f90: :tangle energy_hydrogen.f90
:END:
2021-01-11 20:54:40 +01:00
If the space is discretized in small volume elements $\mathbf{r}_i$
of size $\delta \mathbf{r}$, the expression of $\langle E_L \rangle_{\Psi^2}$
becomes a weighted average of the local energy, where the weights
are the values of the probability density at $\mathbf{r}_i$
multiplied by the volume element:
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
2021-01-11 18:41:36 +01:00
\langle E \rangle_{\Psi^2} \approx \frac{\sum_i w_i E_L(\mathbf{r}_i)}{\sum_i w_i}, \;\;
2021-01-07 11:07:18 +01:00
w_i = \left[\Psi(\mathbf{r}_i)\right]^2 \delta \mathbf{r}
$$
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
#+begin_note
The energy is biased because:
2021-01-11 18:41:36 +01:00
- The volume elements are not infinitely small (discretization error)
- The energy is evaluated only inside the box (incompleteness of the space)
2021-01-11 20:54:40 +01:00
#+end_note
2021-01-07 10:01:55 +01:00
2021-01-12 01:01:52 +01:00
*** Exercise
#+begin_exercise
Compute a numerical estimate of the energy in a grid of
$50\times50\times50$ points in the range $(-5,-5,-5) \le
\mathbf{r} \le (5,5,5)$.
#+end_exercise
*Python*
#+BEGIN_SRC python :results none
2021-01-03 18:45:58 +01:00
import numpy as np
from hydrogen import e_loc, psi
2021-01-11 18:41:36 +01:00
interval = np.linspace(-5,5,num=50)
delta = (interval[1]-interval[0])**3
2021-01-07 11:07:18 +01:00
2021-01-11 18:41:36 +01:00
r = np.array([0.,0.,0.])
2021-01-07 11:07:18 +01:00
2021-01-11 18:41:36 +01:00
for a in [0.1, 0.2, 0.5, 0.9, 1., 1.5, 2.]:
E = 0.
norm = 0.
for x in interval:
r[0] = x
for y in interval:
r[1] = y
for z in interval:
r[2] = z
w = psi(a,r)
w = w * w * delta
E += w * e_loc(a,r)
norm += w
E = E / norm
print(f"a = {a} \t E = {E}")
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-07 11:07:18 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
: a = 0.1 E = -0.24518438948809218
: a = 0.2 E = -0.26966057967803525
: a = 0.5 E = -0.3856357612517407
: a = 0.9 E = -0.49435709786716214
: a = 1.0 E = -0.5
: a = 1.5 E = -0.39242967082602226
: a = 2.0 E = -0.08086980667844901
2021-01-07 11:07:18 +01:00
2021-01-12 01:01:52 +01:00
*Fortran*
#+begin_src f90
2021-01-03 18:45:58 +01:00
program energy_hydrogen
implicit none
double precision, external :: e_loc, psi
double precision :: x(50), w, delta, energy, dx, r(3), a(6), norm
integer :: i, k, l, j
a = (/ 0.1d0, 0.2d0, 0.5d0, 1.d0, 1.5d0, 2.d0 /)
dx = 10.d0/(size(x)-1)
do i=1,size(x)
x(i) = -5.d0 + (i-1)*dx
end do
delta = dx**3
r(:) = 0.d0
do j=1,size(a)
energy = 0.d0
norm = 0.d0
do i=1,size(x)
r(1) = x(i)
do k=1,size(x)
r(2) = x(k)
do l=1,size(x)
r(3) = x(l)
w = psi(a(j),r)
w = w * w * delta
energy = energy + w * e_loc(a(j), r)
norm = norm + w
end do
end do
end do
energy = energy / norm
print *, 'a = ', a(j), ' E = ', energy
end do
end program energy_hydrogen
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
To compile the Fortran and run it:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+begin_src sh :results output :exports both
2021-01-03 18:45:58 +01:00
gfortran hydrogen.f90 energy_hydrogen.f90 -o energy_hydrogen
./energy_hydrogen
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
: a = 0.10000000000000001 E = -0.24518438948809140
: a = 0.20000000000000001 E = -0.26966057967803236
: a = 0.50000000000000000 E = -0.38563576125173815
: a = 1.0000000000000000 E = -0.50000000000000000
: a = 1.5000000000000000 E = -0.39242967082602065
: a = 2.0000000000000000 E = -8.0869806678448772E-002
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
** Variance of the local energy
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle variance_hydrogen.py
:header-args:f90: :tangle variance_hydrogen.f90
:END:
2021-01-11 18:41:36 +01:00
The variance of the local energy is a functional of $\Psi$
which measures the magnitude of the fluctuations of the local
energy associated with $\Psi$ around the average:
2021-01-07 11:07:18 +01:00
$$
\sigma^2(E_L) = \frac{\int \left[\Psi(\mathbf{r})\right]^2\, \left[
E_L(\mathbf{r}) - E \right]^2 \, d\mathbf{r}}{\int \left[\Psi(\mathbf{r}) \right]^2 d\mathbf{r}}
$$
2021-01-11 18:41:36 +01:00
If the local energy is constant (i.e. $\Psi$ is an eigenfunction of
$\hat{H}$) the variance is zero, so the variance of the local
energy can be used as a measure of the quality of a wave function.
2021-01-12 01:01:52 +01:00
*** Exercise
#+begin_exercise
2021-01-07 11:07:18 +01:00
Compute a numerical estimate of the variance of the local energy
2021-01-12 01:01:52 +01:00
in a grid of $50\times50\times50$ points in the range
$(-5,-5,-5)
\le \mathbf{r} \le (5,5,5)$ for different values of $a$.
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Python*
2021-01-11 18:41:36 +01:00
#+begin_src python :results none
2021-01-03 18:45:58 +01:00
import numpy as np
from hydrogen import e_loc, psi
interval = np.linspace(-5,5,num=50)
delta = (interval[1]-interval[0])**3
r = np.array([0.,0.,0.])
for a in [0.1, 0.2, 0.5, 0.9, 1., 1.5, 2.]:
E = 0.
norm = 0.
for x in interval:
r[0] = x
for y in interval:
r[1] = y
for z in interval:
r[2] = z
w = psi(a, r)
w = w * w * delta
El = e_loc(a, r)
E += w * El
norm += w
2021-01-11 20:54:40 +01:00
E = E / norm
s2 = 0.
2021-01-03 18:45:58 +01:00
for x in interval:
r[0] = x
for y in interval:
r[1] = y
for z in interval:
r[2] = z
w = psi(a, r)
w = w * w * delta
El = e_loc(a, r)
2021-01-07 10:01:55 +01:00
s2 += w * (El - E)**2
2021-01-11 20:54:40 +01:00
s2 = s2 / norm
print(f"a = {a} \t E = {E:10.8f} \t \sigma^2 = {s2:10.8f}")
2021-01-07 11:07:18 +01:00
#+end_src
#+RESULTS:
: a = 0.1 E = -0.24518439 \sigma^2 = 0.02696522
: a = 0.2 E = -0.26966058 \sigma^2 = 0.03719707
: a = 0.5 E = -0.38563576 \sigma^2 = 0.05318597
: a = 0.9 E = -0.49435710 \sigma^2 = 0.00577812
: a = 1.0 E = -0.50000000 \sigma^2 = 0.00000000
: a = 1.5 E = -0.39242967 \sigma^2 = 0.31449671
: a = 2.0 E = -0.08086981 \sigma^2 = 1.80688143
2021-01-11 20:54:40 +01:00
*Fortran*
2021-01-07 11:07:18 +01:00
#+begin_src f90
2021-01-03 18:45:58 +01:00
program variance_hydrogen
implicit none
double precision, external :: e_loc, psi
double precision :: x(50), w, delta, energy, dx, r(3), a(6), norm, s2
integer :: i, k, l, j
a = (/ 0.1d0, 0.2d0, 0.5d0, 1.d0, 1.5d0, 2.d0 /)
dx = 10.d0/(size(x)-1)
do i=1,size(x)
x(i) = -5.d0 + (i-1)*dx
end do
delta = dx**3
r(:) = 0.d0
do j=1,size(a)
energy = 0.d0
norm = 0.d0
do i=1,size(x)
r(1) = x(i)
do k=1,size(x)
r(2) = x(k)
do l=1,size(x)
r(3) = x(l)
w = psi(a(j),r)
w = w * w * delta
energy = energy + w * e_loc(a(j), r)
norm = norm + w
end do
end do
end do
energy = energy / norm
s2 = 0.d0
norm = 0.d0
do i=1,size(x)
r(1) = x(i)
do k=1,size(x)
r(2) = x(k)
do l=1,size(x)
r(3) = x(l)
w = psi(a(j),r)
w = w * w * delta
s2 = s2 + w * ( e_loc(a(j), r) - energy )**2
norm = norm + w
end do
end do
end do
s2 = s2 / norm
print *, 'a = ', a(j), ' E = ', energy, ' s2 = ', s2
end do
2021-01-07 10:01:55 +01:00
2021-01-03 18:45:58 +01:00
end program variance_hydrogen
2021-01-07 11:07:18 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
To compile and run:
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
#+begin_src sh :results output :exports both
2021-01-03 18:45:58 +01:00
gfortran hydrogen.f90 variance_hydrogen.f90 -o variance_hydrogen
./variance_hydrogen
2021-01-07 11:07:18 +01:00
#+end_src
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
#+RESULTS:
: a = 0.10000000000000001 E = -0.24518438948809140 s2 = 2.6965218719733813E-002
: a = 0.20000000000000001 E = -0.26966057967803236 s2 = 3.7197072370217653E-002
: a = 0.50000000000000000 E = -0.38563576125173815 s2 = 5.3185967578488862E-002
: a = 1.0000000000000000 E = -0.50000000000000000 s2 = 0.0000000000000000
: a = 1.5000000000000000 E = -0.39242967082602065 s2 = 0.31449670909180444
: a = 2.0000000000000000 E = -8.0869806678448772E-002 s2 = 1.8068814270851303
2021-01-03 18:45:58 +01:00
* Variational Monte Carlo
2021-01-11 18:41:36 +01:00
Numerical integration with deterministic methods is very efficient
2021-01-11 20:54:40 +01:00
in low dimensions. When the number of dimensions becomes large,
instead of computing the average energy as a numerical integration
on a grid, it is usually more efficient to do a Monte Carlo sampling.
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
Moreover, a Monte Carlo sampling will alow us to remove the bias due
to the discretization of space, and compute a statistical confidence
interval.
2021-01-07 10:01:55 +01:00
2021-01-03 18:45:58 +01:00
** Computation of the statistical error
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle qmc_stats.py
:header-args:f90: :tangle qmc_stats.f90
:END:
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
To compute the statistical error, you need to perform $M$
independent Monte Carlo calculations. You will obtain $M$ different
estimates of the energy, which are expected to have a Gaussian
distribution by the central limit theorem.
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
The estimate of the energy is
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
E = \frac{1}{M} \sum_{i=1}^M E_M
$$
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
The variance of the average energies can be computed as
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
\sigma^2 = \frac{1}{M-1} \sum_{i=1}^{M} (E_M - E)^2
$$
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
And the confidence interval is given by
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
$$
E \pm \delta E, \text{ where } \delta E = \frac{\sigma}{\sqrt{M}}
$$
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
*** Exercise
#+begin_exercise
2021-01-07 11:07:18 +01:00
Write a function returning the average and statistical error of an
input array.
2021-01-12 01:01:52 +01:00
#+end_exercise
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Python*
2021-01-11 18:41:36 +01:00
#+BEGIN_SRC python :results none
2021-01-07 10:01:55 +01:00
from math import sqrt
2021-01-03 18:45:58 +01:00
def ave_error(arr):
M = len(arr)
assert (M>1)
average = sum(arr)/M
variance = 1./(M-1) * sum( [ (x - average)**2 for x in arr ] )
return (average, sqrt(variance/M))
2021-01-07 11:07:18 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
2021-01-11 20:54:40 +01:00
*Fortran*
2021-01-07 11:07:18 +01:00
#+BEGIN_SRC f90
2021-01-07 10:01:55 +01:00
subroutine ave_error(x,n,ave,err)
implicit none
integer, intent(in) :: n
double precision, intent(in) :: x(n)
double precision, intent(out) :: ave, err
double precision :: variance
if (n == 1) then
ave = x(1)
err = 0.d0
else
ave = sum(x(:)) / dble(n)
variance = sum( (x(:) - ave)**2 ) / dble(n-1)
err = dsqrt(variance/dble(n))
endif
end subroutine ave_error
2021-01-07 11:07:18 +01:00
#+END_SRC
2021-01-03 18:45:58 +01:00
** Uniform sampling in the box
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle qmc_uniform.py
:header-args:f90: :tangle qmc_uniform.f90
:END:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
We will now do our first Monte Carlo calculation to compute the
energy of the hydrogen atom.
2021-01-07 11:07:18 +01:00
At every Monte Carlo step:
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
- Draw a random point $\mathbf{r}_i$ in the box $(-5,-5,-5) \le
2021-01-07 11:07:18 +01:00
(x,y,z) \le (5,5,5)$
2021-01-12 01:01:52 +01:00
- Compute $[\Psi(\mathbf{r}_i)]^2$ and accumulate the result in a
variable =normalization=
- Compute $[\Psi(\mathbf{r}_i)]^2 \times E_L(\mathbf{r}_i)$, and accumulate the
result in a variable =energy=
One Monte Carlo run will consist of $N$ Monte Carlo steps. Once all the
steps have been computed, the run returns the average energy
$\bar{E}_k$ over the $N$ steps of the run.
To compute the statistical error, perform $M$ runs. The final
estimate of the energy will be the average over the $\bar{E}_k$,
and the variance of the $\bar{E}_k$ will be used to compute the
statistical error.
*** Exercise
2021-01-03 18:45:58 +01:00
2021-01-12 01:01:52 +01:00
#+begin_exercise
Parameterize the wave function with $a=0.9$. Perform 30
independent Monte Carlo runs, each with 100 000 Monte Carlo
steps. Store the final energies of each run and use this array to
compute the average energy and the associated error bar.
#+end_exercise
2021-01-07 10:01:55 +01:00
2021-01-12 01:01:52 +01:00
*Python*
#+BEGIN_SRC python :results output
2021-01-07 10:01:55 +01:00
from hydrogen import *
from qmc_stats import *
2021-01-11 18:41:36 +01:00
def MonteCarlo(a, nmax):
E = 0.
N = 0.
for istep in range(nmax):
r = np.random.uniform(-5., 5., (3))
w = psi(a,r)
w = w*w
N += w
E += w * e_loc(a,r)
return E/N
a = 0.9
nmax = 100000
X = [MonteCarlo(a,nmax) for i in range(30)]
E, deltaE = ave_error(X)
print(f"E = {E} +/- {deltaE}")
2021-01-12 01:01:52 +01:00
#+END_SRC
2021-01-07 11:07:18 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
: E = -0.4956255109300764 +/- 0.0007082875482711226
2021-01-07 11:07:18 +01:00
2021-01-12 01:01:52 +01:00
*Fortran*
#+BEGIN_SRC f90
2021-01-07 10:01:55 +01:00
subroutine uniform_montecarlo(a,nmax,energy)
implicit none
double precision, intent(in) :: a
integer , intent(in) :: nmax
double precision, intent(out) :: energy
2021-01-07 11:07:18 +01:00
2021-01-07 10:01:55 +01:00
integer*8 :: istep
2021-01-07 11:07:18 +01:00
2021-01-07 10:01:55 +01:00
double precision :: norm, r(3), w
double precision, external :: e_loc, psi
energy = 0.d0
norm = 0.d0
do istep = 1,nmax
call random_number(r)
r(:) = -5.d0 + 10.d0*r(:)
w = psi(a,r)
w = w*w
norm = norm + w
energy = energy + w * e_loc(a,r)
end do
energy = energy / norm
end subroutine uniform_montecarlo
program qmc
implicit none
double precision, parameter :: a = 0.9
integer , parameter :: nmax = 100000
integer , parameter :: nruns = 30
integer :: irun
double precision :: X(nruns)
double precision :: ave, err
do irun=1,nruns
call uniform_montecarlo(a,nmax,X(irun))
enddo
call ave_error(X,nruns,ave,err)
print *, 'E = ', ave, '+/-', err
end program qmc
2021-01-12 01:01:52 +01:00
#+END_SRC
2021-01-07 10:01:55 +01:00
2021-01-12 01:01:52 +01:00
#+begin_src sh :results output :exports both
2021-01-07 10:01:55 +01:00
gfortran hydrogen.f90 qmc_stats.f90 qmc_uniform.f90 -o qmc_uniform
./qmc_uniform
2021-01-12 01:01:52 +01:00
#+end_src
2021-01-07 10:01:55 +01:00
2021-01-12 01:01:52 +01:00
#+RESULTS:
: E = -0.49588321986667677 +/- 7.1758863546737969E-004
2021-01-07 10:01:55 +01:00
2021-01-03 18:45:58 +01:00
** Gaussian sampling
2021-01-07 11:07:18 +01:00
:PROPERTIES:
:header-args:python: :tangle qmc_gaussian.py
:header-args:f90: :tangle qmc_gaussian.f90
:END:
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
We will now improve the sampling and allow to sample in the whole
3D space, correcting the bias related to the sampling in the box.
2021-01-03 18:45:58 +01:00
2021-01-07 11:07:18 +01:00
Instead of drawing uniform random numbers, we will draw Gaussian
random numbers centered on 0 and with a variance of 1.
2021-01-07 10:01:55 +01:00
2021-01-07 11:07:18 +01:00
To obtain Gaussian-distributed random numbers, you can apply the
[[https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform][Box Muller transform]] to uniform random numbers:
2021-01-07 10:01:55 +01:00
2021-01-07 11:07:18 +01:00
\begin{eqnarray*}
z_1 &=& \sqrt{-2 \ln u_1} \cos(2 \pi u_2) \\
z_2 &=& \sqrt{-2 \ln u_1} \sin(2 \pi u_2)
\end{eqnarray*}
2021-01-12 01:01:52 +01:00
Here is a Fortran implementation returning a Gaussian-distributed
n-dimensional vector $\mathbf{z}$;
2021-01-11 20:54:40 +01:00
*Fortran*
2021-01-07 11:07:18 +01:00
#+BEGIN_SRC f90 :tangle qmc_stats.f90
subroutine random_gauss(z,n)
implicit none
integer, intent(in) :: n
double precision, intent(out) :: z(n)
double precision :: u(n+1)
double precision, parameter :: two_pi = 2.d0*dacos(-1.d0)
integer :: i
call random_number(u)
if (iand(n,1) == 0) then
! n is even
do i=1,n,2
z(i) = dsqrt(-2.d0*dlog(u(i)))
z(i+1) = z(i) + dsin( two_pi*u(i+1) )
z(i) = z(i) + dcos( two_pi*u(i+1) )
end do
else
! n is odd
do i=1,n-1,2
z(i) = dsqrt(-2.d0*dlog(u(i)))
z(i+1) = z(i) + dsin( two_pi*u(i+1) )
z(i) = z(i) + dcos( two_pi*u(i+1) )
end do
z(n) = dsqrt(-2.d0*dlog(u(n)))
z(n) = z(n) + dcos( two_pi*u(n+1) )
end if
end subroutine random_gauss
#+END_SRC
2021-01-12 01:01:52 +01:00
Now the sampling probability can be inserted into the equation of the energy:
2021-01-07 11:07:18 +01:00
\[
E = \frac{\int P(\mathbf{r}) \frac{\left[\Psi(\mathbf{r})\right]^2}{P(\mathbf{r})}\, \frac{\hat{H} \Psi(\mathbf{r})}{\Psi(\mathbf{r})}\,d\mathbf{r}}{\int P(\mathbf{r}) \frac{\left[\Psi(\mathbf{r}) \right]^2}{P(\mathbf{r})} d\mathbf{r}}
\]
2021-01-12 01:01:52 +01:00
with the Gaussian probability
2021-01-07 11:07:18 +01:00
\[
2021-01-12 01:01:52 +01:00
P(\mathbf{r}) = \frac{1}{(2 \pi)^{3/2}}\exp\left( -\frac{\mathbf{r}^2}{2} \right).
2021-01-07 11:07:18 +01:00
\]
As the coordinates are drawn with probability $P(\mathbf{r})$, the
average energy can be computed as
$$
E \approx \frac{\sum_i w_i E_L(\mathbf{r}_i)}{\sum_i w_i}, \;\;
w_i = \frac{\left[\Psi(\mathbf{r}_i)\right]^2}{P(\mathbf{r}_i)} \delta \mathbf{r}
$$
2021-01-12 01:01:52 +01:00
*** Exercise
#+begin_exercise
Modify the exercise of the previous section to sample with
Gaussian-distributed random numbers. Can you see an reduction in
the statistical error?
#+end_exercise
*Python*
#+BEGIN_SRC python :results output
2021-01-07 11:07:18 +01:00
from hydrogen import *
from qmc_stats import *
2021-01-07 10:01:55 +01:00
2021-01-03 18:45:58 +01:00
norm_gauss = 1./(2.*np.pi)**(1.5)
def gaussian(r):
2021-01-07 11:07:18 +01:00
return norm_gauss * np.exp(-np.dot(r,r)*0.5)
2021-01-03 18:45:58 +01:00
def MonteCarlo(a,nmax):
E = 0.
N = 0.
for istep in range(nmax):
r = np.random.normal(loc=0., scale=1.0, size=(3))
w = psi(a,r)
w = w*w / gaussian(r)
N += w
E += w * e_loc(a,r)
return E/N