2021-01-03 18:45:58 +01:00
|
|
|
#+TITLE: Quantum Monte Carlo
|
|
|
|
#+AUTHOR: Anthony Scemama, Claudia Filippi
|
2021-01-26 10:02:38 +01:00
|
|
|
#+LANGUAGE: en
|
2021-01-28 15:04:57 +01:00
|
|
|
#+INFOJS_OPT: toc:t mouse:underline path:org-info.js
|
2021-01-03 18:45:58 +01:00
|
|
|
#+STARTUP: latexpreview
|
2021-01-21 12:49:12 +01:00
|
|
|
#+LATEX_CLASS: report
|
|
|
|
#+LATEX_HEADER_EXTRA: \usepackage{minted}
|
2021-01-11 20:54:40 +01:00
|
|
|
#+HTML_HEAD: <link rel="stylesheet" title="Standard" href="worg.css" type="text/css" />
|
2021-01-26 10:02:38 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+OPTIONS: H:4 num:t toc:t \n:nil @:t ::t |:t ^:t -:t f:t *:t <:t
|
2021-01-21 12:49:12 +01:00
|
|
|
#+OPTIONS: TeX:t LaTeX:t skip:nil d:nil todo:t pri:nil tags:not-in-toc
|
2021-01-21 23:24:21 +01:00
|
|
|
# EXCLUDE_TAGS: Python solution
|
|
|
|
# EXCLUDE_TAGS: Fortran solution
|
2021-01-21 12:49:12 +01:00
|
|
|
|
|
|
|
#+BEGIN_SRC elisp :output none :exports none
|
|
|
|
(setq org-latex-listings 'minted
|
|
|
|
org-latex-packages-alist '(("" "minted"))
|
|
|
|
org-latex-pdf-process
|
|
|
|
'("pdflatex -shell-escape -interaction nonstopmode -output-directory %o %f"
|
|
|
|
"pdflatex -shell-escape -interaction nonstopmode -output-directory %o %f"
|
|
|
|
"pdflatex -shell-escape -interaction nonstopmode -output-directory %o %f"))
|
|
|
|
(setq org-latex-minted-options '(("breaklines" "true")
|
|
|
|
("breakanywhere" "true")))
|
|
|
|
(setq org-latex-minted-options
|
|
|
|
'(("frame" "lines")
|
|
|
|
("fontsize" "\\scriptsize")
|
|
|
|
("linenos" "")))
|
|
|
|
(org-beamer-export-to-pdf)
|
|
|
|
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+RESULTS:
|
|
|
|
: /home/scemama/TREX/qmc-lttc/QMC.pdf
|
2021-01-03 18:45:58 +01:00
|
|
|
|
|
|
|
* Introduction
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
This web site is the QMC tutorial of the LTTC winter school
|
|
|
|
[[https://www.irsamc.ups-tlse.fr/lttc/Luchon][Tutorials in Theoretical Chemistry]].
|
|
|
|
|
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-28 15:04:57 +01:00
|
|
|
gives the exact energy of the hydrogen atom and of the H_2 molecule.
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
Code examples will be given in Python and Fortran. You can use
|
|
|
|
whatever language you prefer to write the program.
|
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-03 18:45:58 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
All the quantities are expressed in /atomic units/ (energies,
|
|
|
|
coordinates, etc).
|
|
|
|
|
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-21 23:24:21 +01:00
|
|
|
We will first verify that, for a given value of $a$, $\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-21 23:24:21 +01:00
|
|
|
To do that, we will check if 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-21 23:24:21 +01:00
|
|
|
is constant.
|
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
|
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
$$ \langle f \rangle_p = \int_{-\infty}^\infty p(x)\, f(x)\,dx. $$
|
2021-01-11 18:41:36 +01:00
|
|
|
|
|
|
|
Recall that a probability density function $p(x)$ is non-negative
|
|
|
|
and integrates to one:
|
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
$$ \int_{-\infty}^\infty p(x)\,dx = 1. $$
|
2021-01-11 18:41:36 +01:00
|
|
|
|
|
|
|
|
|
|
|
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
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
Write all the functions of this section in a single file :
|
|
|
|
~hydrogen.py~ if you use Python, or ~hydrogen.f90~ is you use
|
|
|
|
Fortran.
|
|
|
|
|
2021-01-26 00:22:37 +01:00
|
|
|
#+begin_note
|
|
|
|
- When computing a square root in $\mathbb{R}$, *always* make sure
|
|
|
|
that the argument of the square root is non-negative.
|
|
|
|
- When you divide, *always* make sure that you will not divide by zero
|
|
|
|
|
|
|
|
If a /floating-point exception/ can occur, you should make a test
|
|
|
|
to catch the error.
|
|
|
|
#+end_note
|
|
|
|
|
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-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-21 23:24:21 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle none
|
|
|
|
import numpy as np
|
|
|
|
|
|
|
|
def potential(r):
|
|
|
|
# TODO
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
|
|
|
#+BEGIN_SRC f90 :tangle none
|
|
|
|
double precision function potential(r)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: r(3)
|
|
|
|
! TODO
|
|
|
|
end function potential
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
2021-01-11 20:54:40 +01:00
|
|
|
#+BEGIN_SRC python :results none
|
2021-01-03 18:45:58 +01:00
|
|
|
import numpy as np
|
|
|
|
|
|
|
|
def potential(r):
|
2021-01-26 00:22:37 +01:00
|
|
|
distance = np.sqrt(np.dot(r,r))
|
|
|
|
assert (distance > 0)
|
|
|
|
return -1. / distance
|
2021-01-11 20:54:40 +01:00
|
|
|
#+END_SRC
|
2021-01-07 10:01:55 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-11 20:54:40 +01:00
|
|
|
#+BEGIN_SRC f90
|
2021-01-03 18:45:58 +01:00
|
|
|
double precision function potential(r)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: r(3)
|
2021-01-26 00:22:37 +01:00
|
|
|
double precision :: distance
|
|
|
|
distance = dsqrt( r(1)*r(1) + r(2)*r(2) + r(3)*r(3) )
|
|
|
|
if (distance > 0.d0) then
|
|
|
|
potential = -1.d0 / dsqrt( r(1)*r(1) + r(2)*r(2) + r(3)*r(3) )
|
|
|
|
else
|
|
|
|
stop 'potential at r=0.d0 diverges'
|
|
|
|
end if
|
2021-01-03 18:45:58 +01:00
|
|
|
end function potential
|
2021-01-11 20:54:40 +01:00
|
|
|
#+END_SRC
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*** Exercise 2
|
2021-01-12 01:01:52 +01:00
|
|
|
#+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-11 20:54:40 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
def psi(a, r):
|
|
|
|
# TODO
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
double precision function psi(a, r)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, r(3)
|
|
|
|
! TODO
|
|
|
|
end function psi
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results none
|
|
|
|
def psi(a, r):
|
|
|
|
return np.exp(-a*np.sqrt(np.dot(r,r)))
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
*Fortran*
|
2021-01-11 20:54:40 +01:00
|
|
|
#+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-27 13:04:32 +01:00
|
|
|
*** Exercise 3
|
2021-01-12 01:01:52 +01:00
|
|
|
#+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-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
def kinetic(a,r):
|
|
|
|
# TODO
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
double precision function kinetic(a,r)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, r(3)
|
|
|
|
! TODO
|
|
|
|
end function kinetic
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results none
|
|
|
|
def kinetic(a,r):
|
|
|
|
distance = np.sqrt(np.dot(r,r))
|
|
|
|
assert (distance > 0.)
|
|
|
|
return -0.5 * (a**2 - (2.*a)/distance)
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
*Fortran*
|
2021-01-11 20:54:40 +01:00
|
|
|
#+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)
|
2021-01-26 00:22:37 +01:00
|
|
|
double precision :: distance
|
|
|
|
distance = dsqrt( r(1)*r(1) + r(2)*r(2) + r(3)*r(3) )
|
|
|
|
if (distance > 0.d0) then
|
|
|
|
kinetic = -0.5d0 * (a*a - (2.d0*a) / distance)
|
|
|
|
else
|
|
|
|
stop 'kinetic energy diverges at r=0'
|
|
|
|
end if
|
2021-01-03 18:45:58 +01:00
|
|
|
end function kinetic
|
2021-01-11 20:54:40 +01:00
|
|
|
#+END_SRC
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*** Exercise 4
|
2021-01-12 01:01:52 +01:00
|
|
|
#+begin_exercise
|
2021-01-21 23:24:21 +01:00
|
|
|
Write a function which computes the local energy at $\mathbf{r}$,
|
|
|
|
using the previously defined functions.
|
|
|
|
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-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
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
def e_loc(a,r):
|
|
|
|
#TODO
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
double precision function e_loc(a,r)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, r(3)
|
|
|
|
! TODO
|
|
|
|
end function e_loc
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results none
|
|
|
|
def e_loc(a,r):
|
|
|
|
return kinetic(a,r) + potential(r)
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
*Fortran*
|
2021-01-11 20:54:40 +01:00
|
|
|
#+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-27 13:04:32 +01:00
|
|
|
*** Exercise 5
|
|
|
|
|
|
|
|
#+begin_exercise
|
|
|
|
Find the theoretical value of $a$ for which $\Psi$ is an eigenfunction of $\hat{H}$.
|
|
|
|
#+end_exercise
|
|
|
|
|
|
|
|
**** Solution :solution:
|
|
|
|
|
|
|
|
\begin{eqnarray*}
|
|
|
|
E &=& \frac{\hat{H} \Psi}{\Psi} = - \frac{1}{2} \frac{\Delta \Psi}{\Psi} -
|
|
|
|
\frac{1}{|\mathbf{r}|} \\
|
|
|
|
&=& -\frac{1}{2}\left(a^2 - \frac{2a}{\mathbf{|r|}} \right) -
|
|
|
|
\frac{1}{|\mathbf{r}|} \\
|
|
|
|
&=&
|
|
|
|
-\frac{1}{2} a^2 + \frac{a-1}{\mathbf{|r|}}
|
|
|
|
\end{eqnarray*}
|
|
|
|
|
|
|
|
$a=1$ cancels the $1/|r|$ term, and makes the energy constant,
|
|
|
|
equal to -0.5 atomic units.
|
|
|
|
|
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-12 01:01:52 +01:00
|
|
|
|
2021-01-26 00:22:37 +01:00
|
|
|
#+begin_note
|
|
|
|
The potential and the kinetic energy both diverge at $r=0$, so we
|
|
|
|
choose a grid which does not contain the origin.
|
|
|
|
#+end_note
|
|
|
|
|
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
|
2021-01-21 23:24:21 +01:00
|
|
|
local energy along the $x$ axis. In Python, you can use matplotlib
|
|
|
|
for example. In Fortran, it is convenient to write in a text file
|
|
|
|
the values of $x$ and $E_L(\mathbf{r})$ for each point, and use
|
|
|
|
Gnuplot to plot the files.
|
2021-01-12 01:01:52 +01:00
|
|
|
#+end_exercise
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle 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)
|
2021-01-21 23:24:21 +01:00
|
|
|
plt.figure(figsize=(10,5))
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
# TODO
|
|
|
|
|
|
|
|
plt.tight_layout()
|
|
|
|
plt.legend()
|
|
|
|
plt.savefig("plot_py.png")
|
|
|
|
#+end_src
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+begin_src f90 :tangle none
|
2021-01-21 23:24:21 +01:00
|
|
|
program plot
|
|
|
|
implicit none
|
|
|
|
double precision, external :: e_loc
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
double precision :: x(50), dx
|
|
|
|
integer :: i, j
|
|
|
|
|
|
|
|
dx = 10.d0/(size(x)-1)
|
|
|
|
do i=1,size(x)
|
|
|
|
x(i) = -5.d0 + (i-1)*dx
|
|
|
|
end do
|
|
|
|
|
|
|
|
! TODO
|
|
|
|
|
|
|
|
end program plot
|
2021-01-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-21 23:24:21 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
To compile and run:
|
2021-01-21 23:24:21 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
#+begin_src sh :exports both
|
2021-01-21 23:24:21 +01:00
|
|
|
gfortran hydrogen.f90 plot_hydrogen.f90 -o plot_hydrogen
|
|
|
|
./plot_hydrogen > data
|
2021-01-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-21 23:24:21 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
To plot the data using Gnuplot:
|
2021-01-21 23:24:21 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src gnuplot :file plot.png :exports code
|
2021-01-21 23:24:21 +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-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-21 23:24:21 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results none
|
|
|
|
import numpy as np
|
|
|
|
import matplotlib.pyplot as plt
|
|
|
|
|
|
|
|
from hydrogen import e_loc
|
|
|
|
|
|
|
|
x=np.linspace(-5,5)
|
|
|
|
plt.figure(figsize=(10,5))
|
|
|
|
|
|
|
|
for a in [0.1, 0.2, 0.5, 1., 1.5, 2.]:
|
|
|
|
y=np.array([ e_loc(a, np.array([t,0.,0.]) ) for t in x])
|
|
|
|
plt.plot(x,y,label=f"a={a}")
|
|
|
|
|
|
|
|
plt.tight_layout()
|
|
|
|
plt.legend()
|
|
|
|
plt.savefig("plot_py.png")
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+RESULTS:
|
|
|
|
|
|
|
|
[[./plot_py.png]]
|
|
|
|
|
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+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-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :exports none
|
2021-01-03 18:45:58 +01:00
|
|
|
gfortran hydrogen.f90 plot_hydrogen.f90 -o plot_hydrogen
|
|
|
|
./plot_hydrogen > data
|
2021-01-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src gnuplot :file plot.png :exports results
|
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-25 23:52:53 +01:00
|
|
|
#+end_src
|
|
|
|
#+RESULTS:
|
|
|
|
[[file:plot.png]]
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-25 23:52:53 +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
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle none
|
|
|
|
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.]:
|
|
|
|
# TODO
|
|
|
|
print(f"a = {a} \t E = {E}")
|
|
|
|
|
|
|
|
#+end_src
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
|
|
|
#+begin_src f90
|
|
|
|
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
|
|
|
|
|
|
|
|
do j=1,size(a)
|
|
|
|
! TODO
|
|
|
|
print *, 'a = ', a(j), ' E = ', energy
|
|
|
|
end do
|
|
|
|
|
|
|
|
end program energy_hydrogen
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
To compile the Fortran and run it:
|
|
|
|
|
|
|
|
#+begin_src sh :results output :exports code
|
|
|
|
|
|
|
|
gfortran hydrogen.f90 energy_hydrogen.f90 -o energy_hydrogen
|
|
|
|
./energy_hydrogen
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results none :exports both
|
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.
|
2021-01-25 23:52:53 +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
|
|
|
|
E += w * e_loc(a,r)
|
|
|
|
norm += w
|
2021-01-11 18:41:36 +01:00
|
|
|
E = E / norm
|
|
|
|
print(f"a = {a} \t E = {E}")
|
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-25 23:52:53 +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-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+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-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :results output :exports results
|
2021-01-03 18:45:58 +01:00
|
|
|
gfortran hydrogen.f90 energy_hydrogen.f90 -o energy_hydrogen
|
|
|
|
./energy_hydrogen
|
2021-01-25 23:52:53 +01:00
|
|
|
#+end_src
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-25 23:52:53 +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
|
|
|
|
|
|
|
|
** 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
|
2021-01-25 23:52:53 +01:00
|
|
|
energy associated with $\Psi$ around its 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-19 10:45:36 +01:00
|
|
|
which can be simplified as
|
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
$$ \sigma^2(E_L) = \langle E_L^2 \rangle_{\Psi^2} - \langle E_L \rangle_{\Psi^2}^2.$$
|
2021-01-07 11:07:18 +01:00
|
|
|
|
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-19 10:45:36 +01:00
|
|
|
*** Exercise (optional)
|
|
|
|
#+begin_exercise
|
|
|
|
Prove that :
|
2021-01-25 23:52:53 +01:00
|
|
|
$$\left( \langle E - \langle E \rangle_{\Psi^2} \rangle_{\Psi^2} \right)^2 = \langle E^2 \rangle_{\Psi^2} - \langle E \rangle_{\Psi^2}^2 $$
|
2021-01-19 10:45:36 +01:00
|
|
|
#+end_exercise
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
|
|
|
|
$\bar{E} = \langle E \rangle$ is a constant, so $\langle \bar{E}
|
|
|
|
\rangle = \bar{E}$ .
|
|
|
|
|
|
|
|
\begin{eqnarray*}
|
|
|
|
\langle E - \bar{E} \rangle^2 & = &
|
|
|
|
\langle E^2 - 2 E \bar{E} + \bar{E}^2 \rangle \\
|
|
|
|
&=& \langle E^2 \rangle - 2 \langle E \bar{E} \rangle + \langle \bar{E}^2 \rangle \\
|
|
|
|
&=& \langle E^2 \rangle - 2 \langle E \rangle \bar{E} + \bar{E}^2 \\
|
|
|
|
&=& \langle E^2 \rangle - 2 \bar{E}^2 + \bar{E}^2 \\
|
|
|
|
&=& \langle E^2 \rangle - \bar{E}^2 \\
|
|
|
|
&=& \langle E^2 \rangle - \langle E \rangle^2 \\
|
|
|
|
\end{eqnarray*}
|
2021-01-12 01:01:52 +01:00
|
|
|
*** Exercise
|
|
|
|
#+begin_exercise
|
2021-01-27 15:21:44 +01:00
|
|
|
Add the calculation of the variance to the previous code, and
|
|
|
|
compute a numerical estimate of the variance of the local energy 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$.
|
2021-01-12 01:01:52 +01:00
|
|
|
#+end_exercise
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+begin_src python :results none :tangle none
|
2021-01-27 15:21:44 +01:00
|
|
|
import numpy as np from hydrogen import e_loc, psi
|
2021-01-25 23:52:53 +01:00
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
interval = np.linspace(-5,5,num=50) delta =
|
|
|
|
(interval[1]-interval[0])**3
|
2021-01-25 23:52:53 +01:00
|
|
|
|
|
|
|
r = np.array([0.,0.,0.])
|
|
|
|
|
|
|
|
for a in [0.1, 0.2, 0.5, 0.9, 1., 1.5, 2.]:
|
|
|
|
# TODO
|
2021-01-27 15:21:44 +01:00
|
|
|
print(f"a = {a} \t E = {E:10.8f} \t \sigma^2 = {s2:10.8f}")
|
|
|
|
#+end_src
|
2021-01-25 23:52:53 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
|
|
|
#+begin_src f90 :tangle none
|
2021-01-27 15:21:44 +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 double precision :: e, energy2 integer :: i, k, l, j
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
a = (/ 0.1d0, 0.2d0, 0.5d0, 1.d0, 1.5d0, 2.d0 /)
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
dx = 10.d0/(size(x)-1) do i=1,size(x) x(i) = -5.d0 + (i-1)*dx end do
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
delta = dx**3
|
|
|
|
|
|
|
|
r(:) = 0.d0
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
do j=1,size(a) ! TODO print *, 'a = ', a(j), ' E = ', energy, ' s2 =
|
|
|
|
', s2 end do
|
2021-01-27 13:04:32 +01:00
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
end program variance_hydrogen #+end_src
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
To compile and run:
|
|
|
|
|
|
|
|
#+begin_src sh :results output :exports both
|
|
|
|
gfortran hydrogen.f90 variance_hydrogen.f90 -o variance_hydrogen
|
2021-01-27 15:21:44 +01:00
|
|
|
./variance_hydrogen #+end_src
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+begin_src python :results none :exports both
|
2021-01-27 15:21:44 +01:00
|
|
|
import numpy as np from hydrogen import e_loc, psi
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
interval = np.linspace(-5,5,num=50) delta =
|
|
|
|
(interval[1]-interval[0])**3
|
2021-01-03 18:45:58 +01:00
|
|
|
|
|
|
|
r = np.array([0.,0.,0.])
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
for a in [0.1, 0.2, 0.5, 0.9, 1., 1.5, 2.]: E = 0. E2 = 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 E2 += w * El*El norm += w E = E / norm E2 = E2 /
|
|
|
|
norm s2 = E2 - E*E print(f"a = {a} \t E = {E:10.8f} \t \sigma^2 =
|
|
|
|
{s2:10.8f}") #+end_src
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
#+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-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-25 23:52:53 +01:00
|
|
|
#+begin_src f90
|
2021-01-27 15:21:44 +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 double precision :: e, energy2 integer :: i, k, l, j
|
2021-01-03 18:45:58 +01:00
|
|
|
|
|
|
|
a = (/ 0.1d0, 0.2d0, 0.5d0, 1.d0, 1.5d0, 2.d0 /)
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
dx = 10.d0/(size(x)-1) do i=1,size(x) x(i) = -5.d0 + (i-1)*dx end do
|
2021-01-03 18:45:58 +01:00
|
|
|
|
|
|
|
delta = dx**3
|
|
|
|
|
|
|
|
r(:) = 0.d0
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
do j=1,size(a) energy = 0.d0 energy2 = 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 e = e_loc(a(j), r)
|
|
|
|
energy = energy + w * e energy2 = energy2 + w * e * e norm =
|
|
|
|
norm + w end do end do end do energy = energy / norm energy2 =
|
|
|
|
energy2 / norm s2 = energy2 - energy*energy print *, 'a = ',
|
|
|
|
a(j), ' E = ', energy, ' s2 = ', s2 end do
|
2021-01-07 10:01:55 +01:00
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
end program variance_hydrogen #+end_src
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :results output :exports results
|
2021-01-03 18:45:58 +01:00
|
|
|
gfortran hydrogen.f90 variance_hydrogen.f90 -o variance_hydrogen
|
2021-01-27 15:21:44 +01:00
|
|
|
./variance_hydrogen #+end_src
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: a = 0.10000000000000001 E = -0.24518438948809140 s2 = 2.6965218719722767E-002
|
|
|
|
: a = 0.20000000000000001 E = -0.26966057967803236 s2 = 3.7197072370201284E-002
|
|
|
|
: a = 0.50000000000000000 E = -0.38563576125173815 s2 = 5.3185967578480653E-002
|
|
|
|
: a = 1.0000000000000000 E = -0.50000000000000000 s2 = 0.0000000000000000
|
|
|
|
: a = 1.5000000000000000 E = -0.39242967082602065 s2 = 0.31449670909172917
|
|
|
|
: a = 2.0000000000000000 E = -8.0869806678448772E-002 s2 = 1.8068814270846534
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-26 00:22:37 +01:00
|
|
|
* Variational Monte Carlo
|
2021-01-03 18:45:58 +01:00
|
|
|
|
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-26 00:22:37 +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
|
2021-01-26 00:22:37 +01:00
|
|
|
distribution according to the [[https://en.wikipedia.org/wiki/Central_limit_theorem][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-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-26 00:22:37 +01:00
|
|
|
#+BEGIN_SRC python :results none :tangle none
|
|
|
|
from math import sqrt
|
|
|
|
def ave_error(arr):
|
|
|
|
#TODO
|
|
|
|
return (average, error)
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
|
|
|
#+BEGIN_SRC f90
|
|
|
|
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
|
|
|
|
! TODO
|
|
|
|
end subroutine ave_error
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results none :exports code
|
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)
|
2021-01-26 00:22:37 +01:00
|
|
|
assert(M>0)
|
|
|
|
if M == 1:
|
|
|
|
return (arr[0], 0.)
|
|
|
|
else:
|
|
|
|
average = sum(arr)/M
|
|
|
|
variance = 1./(M-1) * sum( [ (x - average)**2 for x in arr ] )
|
|
|
|
error = sqrt(variance/M)
|
|
|
|
return (average, error)
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
|
|
|
#+BEGIN_SRC f90 :exports both
|
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
|
2021-01-26 00:22:37 +01:00
|
|
|
if (n < 1) then
|
|
|
|
stop 'n<1 in ave_error'
|
|
|
|
else if (n == 1) then
|
2021-01-07 10:01:55 +01:00
|
|
|
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-26 00:22:37 +01:00
|
|
|
#+END_SRC
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-26 10:02:38 +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-26 10:02:38 +01:00
|
|
|
At every Monte Carlo iteration:
|
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=
|
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
One Monte Carlo run will consist of $N$ Monte Carlo iterations. Once all the
|
|
|
|
iterations have been computed, the run returns the average energy
|
|
|
|
$\bar{E}_k$ over the $N$ iterations of the run.
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
To compute the statistical error, perform $M$ independent 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.
|
2021-01-12 01:01:52 +01:00
|
|
|
|
|
|
|
*** 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-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-26 10:02:38 +01:00
|
|
|
#+begin_note
|
|
|
|
To draw a uniform random number in Python, you can use
|
|
|
|
the [[https://numpy.org/doc/stable/reference/random/generated/numpy.random.uniform.html][~random.uniform~]] function of Numpy.
|
|
|
|
#+end_note
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+BEGIN_SRC python :tangle none :exports code
|
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
|
|
|
def MonteCarlo(a, nmax):
|
|
|
|
# TODO
|
2021-01-11 18:41:36 +01:00
|
|
|
|
|
|
|
a = 0.9
|
|
|
|
nmax = 100000
|
2021-01-27 13:04:32 +01:00
|
|
|
#TODO
|
2021-01-11 18:41:36 +01:00
|
|
|
print(f"E = {E} +/- {deltaE}")
|
2021-01-26 10:02:38 +01:00
|
|
|
#+END_SRC
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-26 10:02:38 +01:00
|
|
|
#+begin_note
|
|
|
|
To draw a uniform random number in Fortran, you can use
|
|
|
|
the [[https://gcc.gnu.org/onlinedocs/gfortran/RANDOM_005fNUMBER.html][~RANDOM_NUMBER~]] subroutine.
|
|
|
|
#+end_note
|
|
|
|
|
|
|
|
#+begin_note
|
|
|
|
When running Monte Carlo calculations, the number of steps is
|
|
|
|
usually very large. We expect =nmax= to be possibly larger than 2
|
|
|
|
billion, so we use 8-byte integers (=integer*8=) to represent it, as
|
|
|
|
well as the index of the current step.
|
|
|
|
#+end_note
|
2021-01-12 10:55:00 +01:00
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle none
|
|
|
|
subroutine uniform_montecarlo(a,nmax,energy)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a
|
|
|
|
integer*8 , intent(in) :: nmax
|
|
|
|
double precision, intent(out) :: energy
|
|
|
|
|
|
|
|
integer*8 :: istep
|
|
|
|
|
|
|
|
double precision :: norm, r(3), w
|
|
|
|
|
|
|
|
double precision, external :: e_loc, psi
|
|
|
|
|
|
|
|
! TODO
|
|
|
|
end subroutine uniform_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9
|
|
|
|
integer*8 , parameter :: nmax = 100000
|
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
|
|
|
double precision :: X(nruns)
|
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
!TODO
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
|
|
|
end program qmc
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :results output :exports code
|
2021-01-26 10:02:38 +01:00
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 qmc_uniform.f90 -o qmc_uniform
|
|
|
|
./qmc_uniform
|
|
|
|
#+end_src
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results output :exports both
|
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
|
|
|
def MonteCarlo(a, nmax):
|
|
|
|
energy = 0.
|
|
|
|
normalization = 0.
|
|
|
|
for istep in range(nmax):
|
|
|
|
r = np.random.uniform(-5., 5., (3))
|
|
|
|
w = psi(a,r)
|
|
|
|
w = w*w
|
|
|
|
normalization += w
|
|
|
|
energy += w * e_loc(a,r)
|
|
|
|
return energy/normalization
|
|
|
|
|
|
|
|
a = 0.9
|
|
|
|
nmax = 100000
|
|
|
|
X = [MonteCarlo(a,nmax) for i in range(30)]
|
|
|
|
E, deltaE = ave_error(X)
|
|
|
|
print(f"E = {E} +/- {deltaE}")
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.4956255109300764 +/- 0.0007082875482711226
|
|
|
|
|
|
|
|
*Fortran*
|
|
|
|
#+begin_note
|
|
|
|
When running Monte Carlo calculations, the number of steps is
|
|
|
|
usually very large. We expect =nmax= to be possibly larger than 2
|
|
|
|
billion, so we use 8-byte integers (=integer*8=) to represent it, as
|
|
|
|
well as the index of the current step.
|
|
|
|
#+end_note
|
2021-01-26 10:02:38 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+BEGIN_SRC f90 :exports code
|
2021-01-07 10:01:55 +01:00
|
|
|
subroutine uniform_montecarlo(a,nmax,energy)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a
|
2021-01-12 10:55:00 +01:00
|
|
|
integer*8 , intent(in) :: nmax
|
2021-01-07 10:01:55 +01:00
|
|
|
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
|
2021-01-12 10:55:00 +01:00
|
|
|
integer*8 , parameter :: nmax = 100000
|
2021-01-07 10:01:55 +01:00
|
|
|
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-27 13:04:32 +01:00
|
|
|
#+END_SRC
|
2021-01-07 10:01:55 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :results output :exports results
|
2021-01-07 10:01:55 +01:00
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 qmc_uniform.f90 -o qmc_uniform
|
|
|
|
./qmc_uniform
|
2021-01-27 13:04:32 +01:00
|
|
|
#+end_src
|
2021-01-07 10:01:55 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49588321986667677 +/- 7.1758863546737969E-004
|
2021-01-07 10:01:55 +01:00
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
** Metropolis sampling with $\Psi^2$
|
2021-01-07 11:07:18 +01:00
|
|
|
:PROPERTIES:
|
2021-01-20 18:31:49 +01:00
|
|
|
:header-args:python: :tangle qmc_metropolis.py
|
|
|
|
:header-args:f90: :tangle qmc_metropolis.f90
|
2021-01-07 11:07:18 +01:00
|
|
|
:END:
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
We will now use the square of the wave function to sample random
|
|
|
|
points distributed with the probability density
|
2021-01-07 11:07:18 +01:00
|
|
|
\[
|
2021-01-20 18:31:49 +01:00
|
|
|
P(\mathbf{r}) = \left[\Psi(\mathbf{r})\right]^2
|
2021-01-07 11:07:18 +01:00
|
|
|
\]
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
The expression of the average energy is now simplified as the average of
|
2021-01-20 18:31:49 +01:00
|
|
|
the local energies, since the weights are taken care of by the
|
|
|
|
sampling :
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
$$
|
|
|
|
E \approx \frac{1}{M}\sum_{i=1}^M E_L(\mathbf{r}_i)
|
|
|
|
$$
|
|
|
|
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
To sample a chosen probability density, an efficient method is the
|
|
|
|
[[https://en.wikipedia.org/wiki/Metropolis%E2%80%93Hastings_algorithm][Metropolis-Hastings sampling algorithm]]. Starting from a random
|
|
|
|
initial position $\mathbf{r}_0$, we will realize a random walk as follows:
|
2021-01-07 11:07:18 +01:00
|
|
|
|
|
|
|
$$
|
2021-01-28 15:04:57 +01:00
|
|
|
\mathbf{r}_{n+1} = \mathbf{r}_{n} + \delta t\, \mathbf{u}
|
2021-01-07 11:07:18 +01:00
|
|
|
$$
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
where $\delta t$ is a fixed constant (the so-called /time-step/), and
|
2021-01-20 18:31:49 +01:00
|
|
|
$\mathbf{u}$ is a uniform random number in a 3-dimensional box
|
|
|
|
$(-1,-1,-1) \le \mathbf{u} \le (1,1,1)$. We will then add the
|
2021-01-26 10:02:38 +01:00
|
|
|
accept/reject step that guarantees that the distribution of the
|
2021-01-20 18:31:49 +01:00
|
|
|
$\mathbf{r}_n$ is $\Psi^2$:
|
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
1) Compute $\Psi$ at a new position $\mathbf{r'} = \mathbf{r}_n +
|
2021-01-28 15:04:57 +01:00
|
|
|
\delta t\, \mathbf{u}$
|
2021-01-26 10:02:38 +01:00
|
|
|
2) Compute the ratio $R = \frac{\left[\Psi(\mathbf{r'})\right]^2}{\left[\Psi(\mathbf{r}_{n})\right]^2}$
|
|
|
|
3) Draw a uniform random number $v \in [0,1]$
|
|
|
|
4) if $v \le R$, accept the move : set $\mathbf{r}_{n+1} = \mathbf{r'}$
|
|
|
|
5) else, reject the move : set $\mathbf{r}_{n+1} = \mathbf{r}_n$
|
|
|
|
6) evaluate the local energy at $\mathbf{r}_{n+1}$
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
#+begin_note
|
|
|
|
A common error is to remove the rejected samples from the
|
|
|
|
calculation of the average. *Don't do it!*
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
All samples should be kept, from both accepted and rejected moves.
|
|
|
|
#+end_note
|
|
|
|
|
|
|
|
If the time step is infinitely small, the ratio will be very close
|
|
|
|
to one and all the steps will be accepted. But the trajectory will
|
|
|
|
be infinitely too short to have statistical significance.
|
|
|
|
|
|
|
|
On the other hand, as the time step increases, the number of
|
|
|
|
accepted steps will decrease because the ratios might become
|
|
|
|
small. If the number of accepted steps is close to zero, then the
|
|
|
|
space is not well sampled either.
|
|
|
|
|
|
|
|
The time step should be adjusted so that it is as large as
|
|
|
|
possible, keeping the number of accepted steps not too small. To
|
2021-01-26 10:02:38 +01:00
|
|
|
achieve that, we define the acceptance rate as the number of
|
2021-01-20 18:31:49 +01:00
|
|
|
accepted steps over the total number of steps. Adjusting the time
|
|
|
|
step such that the acceptance rate is close to 0.5 is a good compromise.
|
|
|
|
|
|
|
|
|
|
|
|
*** Exercise
|
|
|
|
|
2021-01-12 01:01:52 +01:00
|
|
|
#+begin_exercise
|
2021-01-26 10:02:38 +01:00
|
|
|
Modify the program of the previous section to compute the energy,
|
|
|
|
sampled with $\Psi^2$.
|
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
Compute also the acceptance rate, so that you can adapt the time
|
|
|
|
step in order to have an acceptance rate close to 0.5 .
|
2021-01-26 10:02:38 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
Can you observe a reduction in the statistical error?
|
2021-01-12 01:01:52 +01:00
|
|
|
#+end_exercise
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-26 10:02:38 +01:00
|
|
|
#+BEGIN_SRC python :results output :tangle none
|
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-28 15:04:57 +01:00
|
|
|
def MonteCarlo(a,nmax,dt):
|
2021-01-26 10:02:38 +01:00
|
|
|
# TODO
|
|
|
|
return energy/nmax, N_accep/nmax
|
|
|
|
|
|
|
|
# Run simulation
|
|
|
|
a = 0.9
|
|
|
|
nmax = 100000
|
2021-01-28 15:04:57 +01:00
|
|
|
dt = 1.3
|
|
|
|
X0 = [ MonteCarlo(a,nmax,dt) for i in range(30)]
|
2021-01-26 10:02:38 +01:00
|
|
|
|
|
|
|
# Energy
|
|
|
|
X = [ x for (x, _) in X0 ]
|
|
|
|
E, deltaE = ave_error(X)
|
|
|
|
print(f"E = {E} +/- {deltaE}")
|
|
|
|
|
|
|
|
# Acceptance rate
|
|
|
|
X = [ x for (_, x) in X0 ]
|
|
|
|
A, deltaA = ave_error(X)
|
|
|
|
print(f"A = {A} +/- {deltaA}")
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-26 10:02:38 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle none
|
2021-01-28 15:04:57 +01:00
|
|
|
subroutine metropolis_montecarlo(a,nmax,dt,energy,accep)
|
2021-01-26 10:02:38 +01:00
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a
|
|
|
|
integer*8 , intent(in) :: nmax
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, intent(in) :: dt
|
2021-01-26 10:02:38 +01:00
|
|
|
double precision, intent(out) :: energy
|
|
|
|
double precision, intent(out) :: accep
|
|
|
|
|
|
|
|
integer*8 :: istep
|
|
|
|
|
|
|
|
double precision :: r_old(3), r_new(3), psi_old, psi_new
|
|
|
|
double precision :: v, ratio
|
|
|
|
integer*8 :: n_accep
|
|
|
|
double precision, external :: e_loc, psi, gaussian
|
|
|
|
|
|
|
|
! TODO
|
|
|
|
end subroutine metropolis_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9d0
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, parameter :: dt = 1.3d0
|
2021-01-26 10:02:38 +01:00
|
|
|
integer*8 , parameter :: nmax = 100000
|
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
|
|
|
double precision :: X(nruns), Y(nruns)
|
|
|
|
double precision :: ave, err
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
do irun=1,nruns
|
2021-01-28 15:04:57 +01:00
|
|
|
call metropolis_montecarlo(a,nmax,dt,X(irun),Y(irun))
|
2021-01-26 10:02:38 +01:00
|
|
|
enddo
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
|
|
|
|
|
|
|
call ave_error(Y,nruns,ave,err)
|
|
|
|
print *, 'A = ', ave, '+/-', err
|
|
|
|
end program qmc
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+begin_src sh :results output :exports both
|
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 qmc_metropolis.f90 -o qmc_metropolis
|
|
|
|
./qmc_metropolis
|
|
|
|
#+end_src
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results output :exports both
|
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
def MonteCarlo(a,nmax,dt):
|
2021-01-27 13:04:32 +01:00
|
|
|
energy = 0.
|
|
|
|
N_accep = 0
|
2021-01-28 15:04:57 +01:00
|
|
|
r_old = np.random.uniform(-dt, dt, (3))
|
2021-01-27 13:04:32 +01:00
|
|
|
psi_old = psi(a,r_old)
|
|
|
|
for istep in range(nmax):
|
2021-01-28 15:04:57 +01:00
|
|
|
r_new = r_old + np.random.uniform(-dt,dt,(3))
|
2021-01-27 13:04:32 +01:00
|
|
|
psi_new = psi(a,r_new)
|
|
|
|
ratio = (psi_new / psi_old)**2
|
|
|
|
v = np.random.uniform()
|
|
|
|
if v <= ratio:
|
|
|
|
N_accep += 1
|
|
|
|
r_old = r_new
|
|
|
|
psi_old = psi_new
|
|
|
|
energy += e_loc(a,r_old)
|
|
|
|
return energy/nmax, N_accep/nmax
|
|
|
|
|
|
|
|
# Run simulation
|
|
|
|
a = 0.9
|
|
|
|
nmax = 100000
|
2021-01-28 15:04:57 +01:00
|
|
|
dt = 1.3
|
|
|
|
X0 = [ MonteCarlo(a,nmax,dt) for i in range(30)]
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
# Energy
|
|
|
|
X = [ x for (x, _) in X0 ]
|
|
|
|
E, deltaE = ave_error(X)
|
|
|
|
print(f"E = {E} +/- {deltaE}")
|
|
|
|
|
|
|
|
# Acceptance rate
|
|
|
|
X = [ x for (_, x) in X0 ]
|
|
|
|
A, deltaA = ave_error(X)
|
|
|
|
print(f"A = {A} +/- {deltaA}")
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.4950720838131573 +/- 0.00019089638602238043
|
|
|
|
: A = 0.5172960000000001 +/- 0.0003443446549306529
|
|
|
|
|
|
|
|
*Fortran*
|
|
|
|
#+BEGIN_SRC f90 :exports code
|
2021-01-28 15:04:57 +01:00
|
|
|
subroutine metropolis_montecarlo(a,nmax,dt,energy,accep)
|
2021-01-07 11:07:18 +01:00
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a
|
2021-01-12 10:55:00 +01:00
|
|
|
integer*8 , intent(in) :: nmax
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, intent(in) :: dt
|
2021-01-07 11:07:18 +01:00
|
|
|
double precision, intent(out) :: energy
|
2021-01-20 18:31:49 +01:00
|
|
|
double precision, intent(out) :: accep
|
2021-01-07 11:07:18 +01:00
|
|
|
|
|
|
|
integer*8 :: istep
|
|
|
|
|
2021-01-26 10:02:38 +01:00
|
|
|
double precision :: r_old(3), r_new(3), psi_old, psi_new
|
|
|
|
double precision :: v, ratio
|
|
|
|
integer*8 :: n_accep
|
2021-01-07 11:07:18 +01:00
|
|
|
double precision, external :: e_loc, psi, gaussian
|
|
|
|
|
|
|
|
energy = 0.d0
|
2021-01-26 10:02:38 +01:00
|
|
|
n_accep = 0_8
|
2021-01-20 18:31:49 +01:00
|
|
|
call random_number(r_old)
|
2021-01-28 15:04:57 +01:00
|
|
|
r_old(:) = dt * (2.d0*r_old(:) - 1.d0)
|
2021-01-20 18:31:49 +01:00
|
|
|
psi_old = psi(a,r_old)
|
2021-01-07 11:07:18 +01:00
|
|
|
do istep = 1,nmax
|
2021-01-20 18:31:49 +01:00
|
|
|
call random_number(r_new)
|
2021-01-28 15:04:57 +01:00
|
|
|
r_new(:) = r_old(:) + dt * (2.d0*r_new(:) - 1.d0)
|
2021-01-20 18:31:49 +01:00
|
|
|
psi_new = psi(a,r_new)
|
|
|
|
ratio = (psi_new / psi_old)**2
|
|
|
|
call random_number(v)
|
2021-01-26 10:02:38 +01:00
|
|
|
if (v <= ratio) then
|
2021-01-20 18:31:49 +01:00
|
|
|
r_old(:) = r_new(:)
|
|
|
|
psi_old = psi_new
|
2021-01-26 10:02:38 +01:00
|
|
|
n_accep = n_accep + 1_8
|
2021-01-20 18:31:49 +01:00
|
|
|
endif
|
|
|
|
energy = energy + e_loc(a,r_old)
|
2021-01-07 11:07:18 +01:00
|
|
|
end do
|
2021-01-26 10:02:38 +01:00
|
|
|
energy = energy / dble(nmax)
|
|
|
|
accep = dble(n_accep) / dble(nmax)
|
2021-01-20 18:31:49 +01:00
|
|
|
end subroutine metropolis_montecarlo
|
2021-01-07 11:07:18 +01:00
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
2021-01-20 18:31:49 +01:00
|
|
|
double precision, parameter :: a = 0.9d0
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, parameter :: dt = 1.3d0
|
2021-01-12 10:55:00 +01:00
|
|
|
integer*8 , parameter :: nmax = 100000
|
2021-01-07 11:07:18 +01:00
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
2021-01-20 18:31:49 +01:00
|
|
|
double precision :: X(nruns), Y(nruns)
|
2021-01-07 11:07:18 +01:00
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
do irun=1,nruns
|
2021-01-28 15:04:57 +01:00
|
|
|
call metropolis_montecarlo(a,nmax,dt,X(irun),Y(irun))
|
2021-01-07 11:07:18 +01:00
|
|
|
enddo
|
2021-01-26 10:02:38 +01:00
|
|
|
|
2021-01-07 11:07:18 +01:00
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
2021-01-26 10:02:38 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
call ave_error(Y,nruns,ave,err)
|
|
|
|
print *, 'A = ', ave, '+/-', err
|
2021-01-07 11:07:18 +01:00
|
|
|
end program qmc
|
2021-01-26 10:02:38 +01:00
|
|
|
#+END_SRC
|
2021-01-07 11:07:18 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :results output :exports results
|
2021-01-20 18:31:49 +01:00
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 qmc_metropolis.f90 -o qmc_metropolis
|
|
|
|
./qmc_metropolis
|
2021-01-26 10:02:38 +01:00
|
|
|
#+end_src
|
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49515370205041676 +/- 1.7660819245720729E-004
|
|
|
|
: A = 0.51713866666666664 +/- 3.7072551835783688E-004
|
2021-01-20 18:31:49 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
** Gaussian random number generator
|
2021-01-20 19:12:05 +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-11 18:41:36 +01:00
|
|
|
|
2021-01-20 19:12:05 +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-03 18:45:58 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
Below is a Fortran implementation returning a Gaussian-distributed
|
|
|
|
n-dimensional vector $\mathbf{z}$. This will be useful for the
|
|
|
|
following sections.
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
*Fortran*
|
|
|
|
#+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
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
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-03 18:45:58 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
In Python, you can use the [[https://numpy.org/doc/stable/reference/random/generated/numpy.random.normal.html][~random.normal~]] function of Numpy.
|
|
|
|
** Generalized Metropolis algorithm
|
2021-01-12 10:55:00 +01:00
|
|
|
:PROPERTIES:
|
2021-01-20 19:12:05 +01:00
|
|
|
:header-args:python: :tangle vmc_metropolis.py
|
|
|
|
:header-args:f90: :tangle vmc_metropolis.f90
|
2021-01-12 10:55:00 +01:00
|
|
|
:END:
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
One can use more efficient numerical schemes to move the electrons,
|
|
|
|
but the Metropolis accepation step has to be adapted accordingly:
|
|
|
|
the acceptance
|
2021-01-20 19:12:05 +01:00
|
|
|
probability $A$ is chosen so that it is consistent with the
|
|
|
|
probability of leaving $\mathbf{r}_n$ and the probability of
|
|
|
|
entering $\mathbf{r}_{n+1}$:
|
|
|
|
|
|
|
|
\[ A(\mathbf{r}_{n} \rightarrow \mathbf{r}_{n+1}) = \min \left( 1,
|
|
|
|
\frac{T(\mathbf{r}_{n+1} \rightarrow \mathbf{r}_{n}) P(\mathbf{r}_{n+1})}
|
|
|
|
{T(\mathbf{r}_{n} \rightarrow \mathbf{r}_{n+1}) P(\mathbf{r}_{n})}
|
|
|
|
\right)
|
|
|
|
\]
|
|
|
|
where $T(\mathbf{r}_n \rightarrow \mathbf{r}_{n+1})$ is the
|
|
|
|
probability of transition from $\mathbf{r}_n$ to
|
|
|
|
$\mathbf{r}_{n+1}$.
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
In the previous example, we were using uniform random
|
|
|
|
numbers. Hence, the transition probability was
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
\[
|
2021-01-26 13:11:57 +01:00
|
|
|
T(\mathbf{r}_{n} \rightarrow \mathbf{r}_{n+1}) =
|
2021-01-20 19:12:05 +01:00
|
|
|
\text{constant}
|
|
|
|
\]
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
so the expression of $A$ was simplified to the ratios of the squared
|
2021-01-20 19:12:05 +01:00
|
|
|
wave functions.
|
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
Now, if instead of drawing uniform random numbers we
|
|
|
|
choose to draw Gaussian random numbers with zero mean and variance
|
2021-01-28 15:04:57 +01:00
|
|
|
$\delta t$, the transition probability becomes:
|
2021-01-20 19:12:05 +01:00
|
|
|
|
|
|
|
\[
|
2021-01-26 13:11:57 +01:00
|
|
|
T(\mathbf{r}_{n} \rightarrow \mathbf{r}_{n+1}) =
|
2021-01-28 15:04:57 +01:00
|
|
|
\frac{1}{(2\pi\,\delta t)^{3/2}} \exp \left[ - \frac{\left(
|
|
|
|
\mathbf{r}_{n+1} - \mathbf{r}_{n} \right)^2}{2\delta t} \right]
|
2021-01-20 19:12:05 +01:00
|
|
|
\]
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
|
|
|
|
To sample even better the density, we can "push" the electrons
|
|
|
|
into in the regions of high probability, and "pull" them away from
|
|
|
|
the low-probability regions. This will mechanically increase the
|
|
|
|
acceptance ratios and improve the sampling.
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
To do this, we can add the drift vector
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
\[
|
2021-01-26 13:11:57 +01:00
|
|
|
\frac{\nabla [ \Psi^2 ]}{\Psi^2} = 2 \frac{\nabla \Psi}{\Psi}.
|
|
|
|
\]
|
2021-01-20 19:12:05 +01:00
|
|
|
|
|
|
|
The numerical scheme becomes a drifted diffusion:
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
\[
|
2021-01-28 15:04:57 +01:00
|
|
|
\mathbf{r}_{n+1} = \mathbf{r}_{n} + \delta t\, \frac{\nabla
|
2021-01-20 19:12:05 +01:00
|
|
|
\Psi(\mathbf{r})}{\Psi(\mathbf{r})} + \chi
|
|
|
|
\]
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
where $\chi$ is a Gaussian random variable with zero mean and
|
2021-01-28 15:04:57 +01:00
|
|
|
variance $\delta t$.
|
2021-01-20 19:12:05 +01:00
|
|
|
The transition probability becomes:
|
|
|
|
|
|
|
|
\[
|
2021-01-26 13:11:57 +01:00
|
|
|
T(\mathbf{r}_{n} \rightarrow \mathbf{r}_{n+1}) =
|
2021-01-28 15:04:57 +01:00
|
|
|
\frac{1}{(2\pi\,\delta t)^{3/2}} \exp \left[ - \frac{\left(
|
2021-01-20 19:12:05 +01:00
|
|
|
\mathbf{r}_{n+1} - \mathbf{r}_{n} - \frac{\nabla
|
2021-01-28 15:04:57 +01:00
|
|
|
\Psi(\mathbf{r}_n)}{\Psi(\mathbf{r}_n)} \right)^2}{2\,\delta t} \right]
|
2021-01-20 19:12:05 +01:00
|
|
|
\]
|
|
|
|
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
*** Exercise 1
|
2021-01-26 13:11:57 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
#+begin_exercise
|
|
|
|
Write a function to compute the drift vector $\frac{\nabla \Psi(\mathbf{r})}{\Psi(\mathbf{r})}$.
|
|
|
|
#+end_exercise
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-26 13:11:57 +01:00
|
|
|
#+BEGIN_SRC python :tangle hydrogen.py :tangle none
|
|
|
|
def drift(a,r):
|
|
|
|
# TODO
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-26 13:11:57 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle hydrogen.f90 :tangle none
|
|
|
|
subroutine drift(a,r,b)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, r(3)
|
|
|
|
double precision, intent(out) :: b(3)
|
|
|
|
! TODO
|
|
|
|
end subroutine drift
|
2021-01-20 19:12:05 +01:00
|
|
|
#+END_SRC
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :tangle hydrogen.py
|
|
|
|
def drift(a,r):
|
|
|
|
ar_inv = -a/np.sqrt(np.dot(r,r))
|
|
|
|
return r * ar_inv
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
*Fortran*
|
2021-01-20 19:12:05 +01:00
|
|
|
#+BEGIN_SRC f90 :tangle hydrogen.f90
|
2021-01-11 18:41:36 +01:00
|
|
|
subroutine drift(a,r,b)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, r(3)
|
|
|
|
double precision, intent(out) :: b(3)
|
|
|
|
double precision :: ar_inv
|
|
|
|
ar_inv = -a / dsqrt(r(1)*r(1) + r(2)*r(2) + r(3)*r(3))
|
|
|
|
b(:) = r(:) * ar_inv
|
|
|
|
end subroutine drift
|
2021-01-20 19:12:05 +01:00
|
|
|
#+END_SRC
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
*** Exercise 2
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-20 19:12:05 +01:00
|
|
|
#+begin_exercise
|
|
|
|
Modify the previous program to introduce the drifted diffusion scheme.
|
|
|
|
(This is a necessary step for the next section).
|
|
|
|
#+end_exercise
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Python*
|
2021-01-26 13:11:57 +01:00
|
|
|
#+BEGIN_SRC python :results output :tangle none
|
2021-01-12 11:23:13 +01:00
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
def MonteCarlo(a,nmax,dt):
|
2021-01-26 13:11:57 +01:00
|
|
|
# TODO
|
|
|
|
|
|
|
|
# Run simulation
|
|
|
|
a = 0.9
|
|
|
|
nmax = 100000
|
2021-01-28 15:04:57 +01:00
|
|
|
dt = 1.3
|
|
|
|
X0 = [ MonteCarlo(a,nmax,dt) for i in range(30)]
|
2021-01-26 13:11:57 +01:00
|
|
|
|
|
|
|
# Energy
|
|
|
|
X = [ x for (x, _) in X0 ]
|
|
|
|
E, deltaE = ave_error(X)
|
|
|
|
print(f"E = {E} +/- {deltaE}")
|
|
|
|
|
|
|
|
# Acceptance rate
|
|
|
|
X = [ x for (_, x) in X0 ]
|
|
|
|
A, deltaA = ave_error(X)
|
|
|
|
print(f"A = {A} +/- {deltaA}")
|
|
|
|
#+END_SRC
|
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
|
|
|
#+BEGIN_SRC f90 :tangle none
|
2021-01-28 15:04:57 +01:00
|
|
|
subroutine variational_montecarlo(a,dt,nmax,energy,accep_rate)
|
2021-01-27 13:04:32 +01:00
|
|
|
implicit none
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, intent(in) :: a, dt
|
2021-01-27 13:04:32 +01:00
|
|
|
integer*8 , intent(in) :: nmax
|
|
|
|
double precision, intent(out) :: energy, accep_rate
|
|
|
|
|
|
|
|
integer*8 :: istep
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision :: sq_dt, chi(3)
|
2021-01-27 13:04:32 +01:00
|
|
|
double precision :: psi_old, psi_new
|
|
|
|
double precision :: r_old(3), r_new(3)
|
|
|
|
double precision :: d_old(3), d_new(3)
|
|
|
|
double precision, external :: e_loc, psi
|
|
|
|
|
|
|
|
! TODO
|
|
|
|
|
|
|
|
end subroutine variational_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, parameter :: dt = 1.0
|
2021-01-27 13:04:32 +01:00
|
|
|
integer*8 , parameter :: nmax = 100000
|
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
|
|
|
double precision :: X(nruns), accep(nruns)
|
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
do irun=1,nruns
|
2021-01-28 15:04:57 +01:00
|
|
|
call variational_montecarlo(a,dt,nmax,X(irun),accep(irun))
|
2021-01-27 13:04:32 +01:00
|
|
|
enddo
|
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
|
|
|
call ave_error(accep,nruns,ave,err)
|
|
|
|
print *, 'A = ', ave, '+/-', err
|
|
|
|
end program qmc
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+begin_src sh :results output :exports code
|
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 vmc_metropolis.f90 -o vmc_metropolis
|
|
|
|
./vmc_metropolis
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
**** Solution :solution:
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results output :exports both
|
2021-01-26 13:11:57 +01:00
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
def MonteCarlo(a,nmax,dt):
|
2021-01-26 13:11:57 +01:00
|
|
|
energy = 0.
|
2021-01-12 11:23:13 +01:00
|
|
|
accep_rate = 0.
|
2021-01-28 15:04:57 +01:00
|
|
|
sq_dt = np.sqrt(dt)
|
2021-01-11 20:54:40 +01:00
|
|
|
r_old = np.random.normal(loc=0., scale=1.0, size=(3))
|
|
|
|
d_old = drift(a,r_old)
|
|
|
|
d2_old = np.dot(d_old,d_old)
|
|
|
|
psi_old = psi(a,r_old)
|
|
|
|
for istep in range(nmax):
|
2021-01-12 11:23:13 +01:00
|
|
|
chi = np.random.normal(loc=0., scale=1.0, size=(3))
|
2021-01-28 15:04:57 +01:00
|
|
|
r_new = r_old + dt * d_old + sq_dt * chi
|
2021-01-11 20:54:40 +01:00
|
|
|
d_new = drift(a,r_new)
|
|
|
|
d2_new = np.dot(d_new,d_new)
|
|
|
|
psi_new = psi(a,r_new)
|
|
|
|
# Metropolis
|
|
|
|
prod = np.dot((d_new + d_old), (r_new - r_old))
|
2021-01-28 15:04:57 +01:00
|
|
|
argexpo = 0.5 * (d2_new - d2_old)*dt + prod
|
2021-01-11 20:54:40 +01:00
|
|
|
q = psi_new / psi_old
|
|
|
|
q = np.exp(-argexpo) * q*q
|
|
|
|
if np.random.uniform() < q:
|
2021-01-12 11:23:13 +01:00
|
|
|
accep_rate += 1.
|
2021-01-11 20:54:40 +01:00
|
|
|
r_old = r_new
|
|
|
|
d_old = d_new
|
|
|
|
d2_old = d2_new
|
|
|
|
psi_old = psi_new
|
2021-01-26 13:11:57 +01:00
|
|
|
energy += e_loc(a,r_old)
|
|
|
|
return energy/nmax, accep_rate/nmax
|
2021-01-11 20:54:40 +01:00
|
|
|
|
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
# Run simulation
|
2021-01-12 11:23:13 +01:00
|
|
|
a = 0.9
|
2021-01-11 20:54:40 +01:00
|
|
|
nmax = 100000
|
2021-01-28 15:04:57 +01:00
|
|
|
dt = 1.3
|
|
|
|
X0 = [ MonteCarlo(a,nmax,dt) for i in range(30)]
|
2021-01-26 13:11:57 +01:00
|
|
|
|
|
|
|
# Energy
|
|
|
|
X = [ x for (x, _) in X0 ]
|
|
|
|
E, deltaE = ave_error(X)
|
|
|
|
print(f"E = {E} +/- {deltaE}")
|
|
|
|
|
|
|
|
# Acceptance rate
|
|
|
|
X = [ x for (_, x) in X0 ]
|
|
|
|
A, deltaA = ave_error(X)
|
|
|
|
print(f"A = {A} +/- {deltaA}")
|
|
|
|
#+END_SRC
|
2021-01-11 20:54:40 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.4951317910667116 +/- 0.00014045774335059988
|
|
|
|
: A = 0.7200673333333333 +/- 0.00045942791345632793
|
2021-01-11 20:54:40 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*Fortran*
|
2021-01-26 13:11:57 +01:00
|
|
|
#+BEGIN_SRC f90
|
2021-01-28 15:04:57 +01:00
|
|
|
subroutine variational_montecarlo(a,dt,nmax,energy,accep_rate)
|
2021-01-26 13:11:57 +01:00
|
|
|
implicit none
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, intent(in) :: a, dt
|
2021-01-26 13:11:57 +01:00
|
|
|
integer*8 , intent(in) :: nmax
|
|
|
|
double precision, intent(out) :: energy, accep_rate
|
|
|
|
|
|
|
|
integer*8 :: istep
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision :: sq_dt, chi(3), d2_old, prod, u
|
2021-01-12 11:23:13 +01:00
|
|
|
double precision :: psi_old, psi_new, d2_new, argexpo, q
|
|
|
|
double precision :: r_old(3), r_new(3)
|
|
|
|
double precision :: d_old(3), d_new(3)
|
|
|
|
double precision, external :: e_loc, psi
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
sq_dt = dsqrt(dt)
|
2021-01-20 19:12:05 +01:00
|
|
|
|
2021-01-12 11:23:13 +01:00
|
|
|
! Initialization
|
2021-01-11 18:41:36 +01:00
|
|
|
energy = 0.d0
|
2021-01-12 11:23:13 +01:00
|
|
|
accep_rate = 0.d0
|
|
|
|
call random_gauss(r_old,3)
|
|
|
|
call drift(a,r_old,d_old)
|
|
|
|
d2_old = d_old(1)*d_old(1) + d_old(2)*d_old(2) + d_old(3)*d_old(3)
|
|
|
|
psi_old = psi(a,r_old)
|
|
|
|
|
2021-01-11 18:41:36 +01:00
|
|
|
do istep = 1,nmax
|
2021-01-12 11:23:13 +01:00
|
|
|
call random_gauss(chi,3)
|
2021-01-28 15:04:57 +01:00
|
|
|
r_new(:) = r_old(:) + dt * d_old(:) + chi(:)*sq_dt
|
2021-01-12 11:23:13 +01:00
|
|
|
call drift(a,r_new,d_new)
|
|
|
|
d2_new = d_new(1)*d_new(1) + d_new(2)*d_new(2) + d_new(3)*d_new(3)
|
|
|
|
psi_new = psi(a,r_new)
|
|
|
|
! Metropolis
|
|
|
|
prod = (d_new(1) + d_old(1))*(r_new(1) - r_old(1)) + &
|
2021-01-26 13:11:57 +01:00
|
|
|
(d_new(2) + d_old(2))*(r_new(2) - r_old(2)) + &
|
|
|
|
(d_new(3) + d_old(3))*(r_new(3) - r_old(3))
|
2021-01-28 15:04:57 +01:00
|
|
|
argexpo = 0.5d0 * (d2_new - d2_old)*dt + prod
|
2021-01-12 11:23:13 +01:00
|
|
|
q = psi_new / psi_old
|
|
|
|
q = dexp(-argexpo) * q*q
|
|
|
|
call random_number(u)
|
|
|
|
if (u<q) then
|
|
|
|
accep_rate = accep_rate + 1.d0
|
|
|
|
r_old(:) = r_new(:)
|
|
|
|
d_old(:) = d_new(:)
|
|
|
|
d2_old = d2_new
|
|
|
|
psi_old = psi_new
|
|
|
|
end if
|
|
|
|
energy = energy + e_loc(a,r_old)
|
2021-01-11 18:41:36 +01:00
|
|
|
end do
|
2021-01-26 13:11:57 +01:00
|
|
|
energy = energy / dble(nmax)
|
|
|
|
accep_rate = dble(accep_rate) / dble(nmax)
|
2021-01-11 18:41:36 +01:00
|
|
|
end subroutine variational_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9
|
2021-01-28 15:04:57 +01:00
|
|
|
double precision, parameter :: dt = 1.0
|
2021-01-12 11:23:13 +01:00
|
|
|
integer*8 , parameter :: nmax = 100000
|
2021-01-11 18:41:36 +01:00
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
2021-01-12 11:23:13 +01:00
|
|
|
double precision :: X(nruns), accep(nruns)
|
2021-01-11 18:41:36 +01:00
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
do irun=1,nruns
|
2021-01-28 15:04:57 +01:00
|
|
|
call variational_montecarlo(a,dt,nmax,X(irun),accep(irun))
|
2021-01-11 18:41:36 +01:00
|
|
|
enddo
|
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
2021-01-12 11:23:13 +01:00
|
|
|
call ave_error(accep,nruns,ave,err)
|
|
|
|
print *, 'A = ', ave, '+/-', err
|
2021-01-11 18:41:36 +01:00
|
|
|
end program qmc
|
2021-01-26 13:11:57 +01:00
|
|
|
#+END_SRC
|
2021-01-11 18:41:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
#+begin_src sh :results output :exports results
|
2021-01-12 11:23:13 +01:00
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 vmc_metropolis.f90 -o vmc_metropolis
|
|
|
|
./vmc_metropolis
|
2021-01-26 13:11:57 +01:00
|
|
|
#+end_src
|
2021-01-12 10:55:00 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49495906384751226 +/- 1.5257646086088266E-004
|
|
|
|
: A = 0.78861366666666666 +/- 3.7855335138754813E-004
|
2021-01-12 01:01:52 +01:00
|
|
|
|
2021-01-27 00:56:36 +01:00
|
|
|
* Diffusion Monte Carlo
|
2021-01-12 10:55:00 +01:00
|
|
|
:PROPERTIES:
|
|
|
|
:header-args:python: :tangle dmc.py
|
|
|
|
:header-args:f90: :tangle dmc.f90
|
|
|
|
:END:
|
2021-01-13 17:59:48 +01:00
|
|
|
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
** Schrödinger equation in imaginary time
|
|
|
|
|
|
|
|
Consider the time-dependent Schrödinger equation:
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
\[
|
|
|
|
i\frac{\partial \Psi(\mathbf{r},t)}{\partial t} = \hat{H} \Psi(\mathbf{r},t)
|
|
|
|
\]
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
We can expand $\Psi(\mathbf{r},0)$, in the basis of the eigenstates
|
|
|
|
of the time-independent Hamiltonian:
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
\[
|
|
|
|
\Psi(\mathbf{r},0) = \sum_k a_k\, \Phi_k(\mathbf{r}).
|
|
|
|
\]
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
The solution of the Schrödinger equation at time $t$ is
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
\[
|
|
|
|
\Psi(\mathbf{r},t) = \sum_k a_k \exp \left( -i\, E_k\, t \right) \Phi_k(\mathbf{r}).
|
|
|
|
\]
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
Now, let's replace the time variable $t$ by an imaginary time variable
|
|
|
|
$\tau=i\,t$, we obtain
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
\[
|
|
|
|
-\frac{\partial \psi(\mathbf{r}, \tau)}{\partial \tau} = \hat{H} \psi(\mathbf{r}, \tau)
|
|
|
|
\]
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
where $\psi(\mathbf{r},\tau) = \Psi(\mathbf{r},-i\tau) = \Psi(\mathbf{r},t)$
|
|
|
|
and
|
|
|
|
\[
|
|
|
|
\psi(\mathbf{r},\tau) = \sum_k a_k \exp( -E_k\, \tau) \phi_k(\mathbf{r}).
|
|
|
|
\]
|
|
|
|
For large positive values of $\tau$, $\psi$ is dominated by the
|
|
|
|
$k=0$ term, namely the lowest eigenstate.
|
|
|
|
So we can expect that simulating the differetial equation in
|
|
|
|
imaginary time will converge to the exact ground state of the
|
|
|
|
system.
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
** Diffusion and branching
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
The [[https://en.wikipedia.org/wiki/Diffusion_equation][diffusion equation]] of particles is given by
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
\[
|
|
|
|
\frac{\partial \phi(\mathbf{r},t)}{\partial t} = D\, \Delta \phi(\mathbf{r},t).
|
|
|
|
\]
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
The [[https://en.wikipedia.org/wiki/Reaction_rate][rate of reaction]] $v$ is the speed at which a chemical reaction
|
|
|
|
takes place. In a solution, the rate is given as a function of the
|
|
|
|
concentration $[A]$ by
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
\[
|
|
|
|
v = \frac{d[A]}{dt},
|
|
|
|
\]
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
where the concentration $[A]$ is proportional to the number of particles.
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
These two equations allow us to interpret the Schrödinger equation
|
|
|
|
in imaginary time as the combination of:
|
|
|
|
- a diffusion equation with a diffusion coefficient $D=1/2$ for the
|
|
|
|
kinetic energy, and
|
|
|
|
- a rate equation for the potential.
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
The diffusion equation can be simulated by a Brownian motion:
|
2021-01-28 15:04:57 +01:00
|
|
|
\[ \mathbf{r}_{n+1} = \mathbf{r}_{n} + \sqrt{\delta t}\, \chi \]
|
2021-01-27 13:04:32 +01:00
|
|
|
where $\chi$ is a Gaussian random variable, and the rate equation
|
|
|
|
can be simulated by creating or destroying particles over time (a
|
|
|
|
so-called branching process).
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
/Diffusion Monte Carlo/ (DMC) consists in obtaining the ground state of a
|
|
|
|
system by simulating the Schrödinger equation in imaginary time, by
|
|
|
|
the combination of a diffusion process and a branching process.
|
|
|
|
|
|
|
|
** Importance sampling
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
In a molecular system, the potential is far from being constant,
|
|
|
|
and diverges at inter-particle coalescence points. Hence, when the
|
|
|
|
rate equation is simulated, it results in very large fluctuations
|
|
|
|
in the numbers of particles, making the calculations impossible in
|
|
|
|
practice.
|
|
|
|
Fortunately, if we multiply the Schrödinger equation by a chosen
|
|
|
|
/trial wave function/ $\Psi_T(\mathbf{r})$ (Hartree-Fock, Kohn-Sham
|
2021-01-27 15:21:44 +01:00
|
|
|
determinant, CI wave function, /etc/), one obtains (see appendix
|
|
|
|
for details)
|
2021-01-27 13:04:32 +01:00
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \psi(\mathbf{r},\tau)}{\partial \tau} \Psi_T(\mathbf{r}) =
|
|
|
|
\left[ -\frac{1}{2} \Delta \psi(\mathbf{r},\tau) + V(\mathbf{r}) \psi(\mathbf{r},\tau) \right] \Psi_T(\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
Defining $\Pi(\mathbf{r},\tau) = \psi(\mathbf{r},\tau)
|
2021-01-27 13:04:32 +01:00
|
|
|
\Psi_T(\mathbf{r})$,
|
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \Pi(\mathbf{r},\tau)}{\partial \tau}
|
|
|
|
= -\frac{1}{2} \Delta \Pi(\mathbf{r},\tau) +
|
|
|
|
\nabla \left[ \Pi(\mathbf{r},\tau) \frac{\nabla \Psi_T(\mathbf{r})}{\Psi_T(\mathbf{r})}
|
|
|
|
\right] + E_L(\mathbf{r}) \Pi(\mathbf{r},\tau)
|
|
|
|
\]
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
The new "kinetic energy" can be simulated by the drifted diffusion
|
|
|
|
scheme presented in the previous section (VMC).
|
2021-01-27 13:04:32 +01:00
|
|
|
The new "potential" is the local energy, which has smaller fluctuations
|
2021-01-28 15:04:57 +01:00
|
|
|
when $\Psi_T$ gets closer to the exact wave function. It can be simulated by
|
|
|
|
changing the number of particles according to $\exp\left[ -\delta t\,
|
|
|
|
\left(E_L(\mathbf{r}) - E_\text{ref}\right)\right]$
|
|
|
|
where $E_{\text{ref}}$ is a constant introduced so that the average
|
|
|
|
of this term is close to one, keeping the number of particles rather
|
|
|
|
constant.
|
2021-01-27 00:56:36 +01:00
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
This equation generates the /N/-electron density $\Pi$, which is the
|
|
|
|
product of the ground state with the trial wave function. It
|
|
|
|
introduces the constraint that $\Pi(\mathbf{r},\tau)=0$ where
|
|
|
|
$\Psi_T(\mathbf{r})=0$. If the wave function has the same sign
|
|
|
|
everywhere, as in the hydrogen atom or the H_2 molecule, this
|
|
|
|
constraint is harmless and we can obtain the exact energy. But for
|
|
|
|
systems where the wave function has nodes (systems with 3 or more
|
|
|
|
electrons), this scheme introduces an error known as the /fixed node
|
|
|
|
error/.
|
2021-01-28 15:04:57 +01:00
|
|
|
|
2021-01-27 13:04:32 +01:00
|
|
|
*** Appendix : Details of the Derivation
|
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \psi(\mathbf{r},\tau)}{\partial \tau} \Psi_T(\mathbf{r}) =
|
|
|
|
\left[ -\frac{1}{2} \Delta \psi(\mathbf{r},\tau) + V(\mathbf{r}) \psi(\mathbf{r},\tau) \right] \Psi_T(\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \big[ \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big]}{\partial \tau}
|
|
|
|
= -\frac{1}{2} \Big( \Delta \big[
|
|
|
|
\psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big] -
|
|
|
|
\psi(\mathbf{r},\tau) \Delta \Psi_T(\mathbf{r}) - 2
|
|
|
|
\nabla \psi(\mathbf{r},\tau) \nabla \Psi_T(\mathbf{r}) \Big) + V(\mathbf{r}) \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \big[ \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big]}{\partial \tau}
|
|
|
|
= -\frac{1}{2} \Delta \big[\psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big] +
|
|
|
|
\frac{1}{2} \psi(\mathbf{r},\tau) \Delta \Psi_T(\mathbf{r}) +
|
|
|
|
\Psi_T(\mathbf{r})\nabla \psi(\mathbf{r},\tau) \frac{\nabla \Psi_T(\mathbf{r})}{\Psi_T(\mathbf{r})} + V(\mathbf{r}) \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \big[ \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big]}{\partial \tau}
|
|
|
|
= -\frac{1}{2} \Delta \big[\psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big] +
|
|
|
|
\psi(\mathbf{r},\tau) \Delta \Psi_T(\mathbf{r}) +
|
|
|
|
\Psi_T(\mathbf{r})\nabla \psi(\mathbf{r},\tau) \frac{\nabla \Psi_T(\mathbf{r})}{\Psi_T(\mathbf{r})} + E_L(\mathbf{r}) \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})
|
|
|
|
\]
|
|
|
|
\[
|
|
|
|
-\frac{\partial \big[ \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big]}{\partial \tau}
|
|
|
|
= -\frac{1}{2} \Delta \big[\psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) \big] +
|
|
|
|
\nabla \left[ \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})
|
|
|
|
\frac{\nabla \Psi_T(\mathbf{r})}{\Psi_T(\mathbf{r})}
|
|
|
|
\right] + E_L(\mathbf{r}) \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
|
|
|
Defining $\Pi(\mathbf{r},t) = \psi(\mathbf{r},\tau)
|
|
|
|
\Psi_T(\mathbf{r})$,
|
|
|
|
|
|
|
|
\[
|
|
|
|
-\frac{\partial \Pi(\mathbf{r},\tau)}{\partial \tau}
|
|
|
|
= -\frac{1}{2} \Delta \Pi(\mathbf{r},\tau) +
|
|
|
|
\nabla \left[ \Pi(\mathbf{r},\tau) \frac{\nabla \Psi_T(\mathbf{r})}{\Psi_T(\mathbf{r})}
|
|
|
|
\right] + E_L(\mathbf{r}) \Pi(\mathbf{r},\tau)
|
|
|
|
\]
|
2021-01-26 17:06:12 +01:00
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
|
|
|
|
** Fixed-node DMC energy
|
|
|
|
|
|
|
|
Now that we have a process to sample $\Pi(\mathbf{r},\tau) =
|
|
|
|
\psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})$, we can compute the exact
|
|
|
|
energy of the system, within the fixed-node constraint, as:
|
|
|
|
|
|
|
|
\[
|
|
|
|
E = \lim_{\tau \to \infty} \frac{\int \Pi(\mathbf{r},\tau) E_L(\mathbf{r}) d\mathbf{r}}
|
|
|
|
{\int \Pi(\mathbf{r},\tau) d\mathbf{r}} = \lim_{\tau \to
|
|
|
|
\infty} E(\tau).
|
|
|
|
\]
|
|
|
|
|
|
|
|
Indeed, this is equivalent to
|
|
|
|
|
|
|
|
\[
|
|
|
|
E(\tau) = \frac{\int \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) E_L(\mathbf{r}) d\mathbf{r}}
|
|
|
|
{\int \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) d\mathbf{r}}
|
|
|
|
= \frac{\int \psi(\mathbf{r},\tau) \hat{H} \Psi_T(\mathbf{r}) d\mathbf{r}}
|
|
|
|
{\int \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) d\mathbf{r}}
|
|
|
|
= \frac{\langle \psi_\tau | \hat{H} | \Psi_T \rangle}
|
|
|
|
{\langle \psi_\tau | \Psi_T \rangle}
|
|
|
|
\]
|
|
|
|
|
|
|
|
As $\hat{H}$ is Hermitian,
|
|
|
|
|
|
|
|
\[
|
|
|
|
E(\tau) = \frac{\langle \psi_\tau | \hat{H} | \Psi_T \rangle}
|
|
|
|
{\langle \psi_\tau | \Psi_T \rangle}
|
|
|
|
= \frac{\langle \Psi_T | \hat{H} | \psi_\tau \rangle}
|
|
|
|
{\langle \Psi_T | \psi_\tau \rangle}
|
|
|
|
= E[\psi_\tau] \frac{\langle \Psi_T | \psi_\tau \rangle}
|
|
|
|
{\langle \Psi_T | \psi_\tau \rangle}
|
|
|
|
= E[\psi_\tau]
|
|
|
|
\]
|
|
|
|
|
|
|
|
So computing the energy within DMC consists in generating the
|
|
|
|
density $\Pi$, and simply averaging the local energies computed
|
|
|
|
with the trial wave function.
|
|
|
|
|
2021-01-27 15:21:44 +01:00
|
|
|
** Pure Diffusion Monte Carlo
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
Instead of having a variable number of particles to simulate the
|
|
|
|
branching process, one can choose to sample $[\Psi_T(\mathbf{r})]^2$ instead of
|
|
|
|
$\psi(\mathbf{r},\tau) \Psi_T(\mathbf{r})$, and consider the term
|
|
|
|
$\exp \left( -\delta t\,( E_L(\mathbf{r}) - E_{\text{ref}} \right)$ as a
|
|
|
|
cumulative product of branching weights:
|
|
|
|
|
|
|
|
\[
|
|
|
|
W(\mathbf{r}_n, \tau)
|
|
|
|
= \exp \left( \int_0^\tau - (E_L(\mathbf{r}_t) - E_{\text{ref}}) dt \right)
|
|
|
|
\approx \prod_{i=1}^{n} \exp \left( -\delta t\, (E_L(\mathbf{r}_i) - E_{\text{ref}}) \right).
|
|
|
|
\]
|
|
|
|
|
|
|
|
The wave function becomes
|
|
|
|
|
|
|
|
\[
|
|
|
|
\psi(\mathbf{r},\tau) = \Psi_T(\mathbf{r}) W(\mathbf{r},\tau)
|
|
|
|
\]
|
|
|
|
|
|
|
|
and the expression of the fixed-node DMC energy is
|
|
|
|
|
|
|
|
\begin{eqnarray*}
|
|
|
|
E(\tau) & = & \frac{\int \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) E_L(\mathbf{r}) d\mathbf{r}}
|
|
|
|
{\int \psi(\mathbf{r},\tau) \Psi_T(\mathbf{r}) d\mathbf{r}} \\
|
|
|
|
& = & \frac{\int \left[ \Psi_T(\mathbf{r}) \right]^2 W(\mathbf{r},\tau) E_L(\mathbf{r}) d\mathbf{r}}
|
|
|
|
{\int \left[ \Psi_T(\mathbf{r}) \right]^2 W(\mathbf{r},\tau) d\mathbf{r}} \\
|
|
|
|
\end{eqnarray*}
|
|
|
|
|
|
|
|
This algorithm is less stable than the branching algorithm: it
|
|
|
|
requires to have a value of $E_\text{ref}$ which is close to the
|
|
|
|
fixed-node energy, and a good trial wave function. Its big
|
|
|
|
advantage is that it is very easy to program starting from a VMC code.
|
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
** TODO Hydrogen atom
|
2021-01-13 17:59:48 +01:00
|
|
|
|
2021-01-21 23:24:21 +01:00
|
|
|
*** Exercise
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-13 17:59:48 +01:00
|
|
|
#+begin_exercise
|
|
|
|
Modify the Metropolis VMC program to introduce the PDMC weight.
|
2021-01-28 15:04:57 +01:00
|
|
|
In the limit $\delta t \rightarrow 0$, you should recover the exact
|
2021-01-13 17:59:48 +01:00
|
|
|
energy of H for any value of $a$.
|
|
|
|
#+end_exercise
|
|
|
|
|
2021-01-27 00:56:36 +01:00
|
|
|
**** Python
|
|
|
|
#+BEGIN_SRC python :results output
|
2021-01-13 17:59:48 +01:00
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
def MonteCarlo(a,nmax,tau,Eref):
|
|
|
|
energy = 0.
|
|
|
|
normalization = 0.
|
|
|
|
accep_rate = 0
|
2021-01-13 17:59:48 +01:00
|
|
|
sq_tau = np.sqrt(tau)
|
|
|
|
r_old = np.random.normal(loc=0., scale=1.0, size=(3))
|
|
|
|
d_old = drift(a,r_old)
|
|
|
|
d2_old = np.dot(d_old,d_old)
|
|
|
|
psi_old = psi(a,r_old)
|
|
|
|
w = 1.0
|
|
|
|
for istep in range(nmax):
|
|
|
|
chi = np.random.normal(loc=0., scale=1.0, size=(3))
|
|
|
|
el = e_loc(a,r_old)
|
|
|
|
w *= np.exp(-tau*(el - Eref))
|
2021-01-26 13:11:57 +01:00
|
|
|
normalization += w
|
|
|
|
energy += w * el
|
2021-01-13 17:59:48 +01:00
|
|
|
|
|
|
|
r_new = r_old + tau * d_old + sq_tau * chi
|
|
|
|
d_new = drift(a,r_new)
|
|
|
|
d2_new = np.dot(d_new,d_new)
|
|
|
|
psi_new = psi(a,r_new)
|
|
|
|
# Metropolis
|
|
|
|
prod = np.dot((d_new + d_old), (r_new - r_old))
|
|
|
|
argexpo = 0.5 * (d2_new - d2_old)*tau + prod
|
|
|
|
q = psi_new / psi_old
|
|
|
|
q = np.exp(-argexpo) * q*q
|
|
|
|
# PDMC weight
|
|
|
|
if np.random.uniform() < q:
|
2021-01-26 13:11:57 +01:00
|
|
|
accep_rate += 1
|
2021-01-13 17:59:48 +01:00
|
|
|
r_old = r_new
|
|
|
|
d_old = d_new
|
|
|
|
d2_old = d2_new
|
|
|
|
psi_old = psi_new
|
2021-01-26 13:11:57 +01:00
|
|
|
return energy/normalization, accep_rate/nmax
|
2021-01-13 17:59:48 +01:00
|
|
|
|
|
|
|
|
|
|
|
a = 0.9
|
|
|
|
nmax = 10000
|
|
|
|
tau = .1
|
2021-01-26 13:11:57 +01:00
|
|
|
E_ref = -0.5
|
|
|
|
X = [MonteCarlo(a,nmax,tau,E_ref) for i in range(30)]
|
2021-01-13 17:59:48 +01:00
|
|
|
E, deltaE = ave_error([x[0] for x in X])
|
|
|
|
A, deltaA = ave_error([x[1] for x in X])
|
|
|
|
print(f"E = {E} +/- {deltaE}\nA = {A} +/- {deltaA}")
|
2021-01-27 00:56:36 +01:00
|
|
|
#+END_SRC
|
2021-01-13 17:59:48 +01:00
|
|
|
|
2021-01-27 00:56:36 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49654807434947584 +/- 0.0006868522447409156
|
|
|
|
: A = 0.9876193891840709 +/- 0.00041857361650995804
|
2021-01-13 17:59:48 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
**** Fortran
|
|
|
|
#+BEGIN_SRC f90
|
2021-01-13 17:59:48 +01:00
|
|
|
subroutine variational_montecarlo(a,tau,nmax,energy,accep_rate)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, tau
|
|
|
|
integer*8 , intent(in) :: nmax
|
|
|
|
double precision, intent(out) :: energy, accep_rate
|
|
|
|
|
|
|
|
integer*8 :: istep
|
|
|
|
double precision :: norm, sq_tau, chi(3), d2_old, prod, u
|
|
|
|
double precision :: psi_old, psi_new, d2_new, argexpo, q
|
|
|
|
double precision :: r_old(3), r_new(3)
|
|
|
|
double precision :: d_old(3), d_new(3)
|
|
|
|
double precision, external :: e_loc, psi
|
|
|
|
|
|
|
|
sq_tau = dsqrt(tau)
|
|
|
|
|
|
|
|
! Initialization
|
|
|
|
energy = 0.d0
|
|
|
|
norm = 0.d0
|
|
|
|
accep_rate = 0.d0
|
|
|
|
call random_gauss(r_old,3)
|
|
|
|
call drift(a,r_old,d_old)
|
|
|
|
d2_old = d_old(1)*d_old(1) + d_old(2)*d_old(2) + d_old(3)*d_old(3)
|
|
|
|
psi_old = psi(a,r_old)
|
|
|
|
|
|
|
|
do istep = 1,nmax
|
|
|
|
call random_gauss(chi,3)
|
|
|
|
r_new(:) = r_old(:) + tau * d_old(:) + chi(:)*sq_tau
|
|
|
|
call drift(a,r_new,d_new)
|
|
|
|
d2_new = d_new(1)*d_new(1) + d_new(2)*d_new(2) + d_new(3)*d_new(3)
|
|
|
|
psi_new = psi(a,r_new)
|
|
|
|
! Metropolis
|
|
|
|
prod = (d_new(1) + d_old(1))*(r_new(1) - r_old(1)) + &
|
|
|
|
(d_new(2) + d_old(2))*(r_new(2) - r_old(2)) + &
|
|
|
|
(d_new(3) + d_old(3))*(r_new(3) - r_old(3))
|
|
|
|
argexpo = 0.5d0 * (d2_new - d2_old)*tau + prod
|
|
|
|
q = psi_new / psi_old
|
|
|
|
q = dexp(-argexpo) * q*q
|
|
|
|
call random_number(u)
|
|
|
|
if (u<q) then
|
|
|
|
accep_rate = accep_rate + 1.d0
|
|
|
|
r_old(:) = r_new(:)
|
|
|
|
d_old(:) = d_new(:)
|
|
|
|
d2_old = d2_new
|
|
|
|
psi_old = psi_new
|
|
|
|
end if
|
|
|
|
norm = norm + 1.d0
|
|
|
|
energy = energy + e_loc(a,r_old)
|
|
|
|
end do
|
|
|
|
energy = energy / norm
|
|
|
|
accep_rate = accep_rate / norm
|
|
|
|
end subroutine variational_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9
|
|
|
|
double precision, parameter :: tau = 1.0
|
|
|
|
integer*8 , parameter :: nmax = 100000
|
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
|
|
|
double precision :: X(nruns), accep(nruns)
|
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
do irun=1,nruns
|
|
|
|
call variational_montecarlo(a,tau,nmax,X(irun),accep(irun))
|
|
|
|
enddo
|
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
|
|
|
call ave_error(accep,nruns,ave,err)
|
|
|
|
print *, 'A = ', ave, '+/-', err
|
|
|
|
end program qmc
|
2021-01-26 13:11:57 +01:00
|
|
|
#+END_SRC
|
2021-01-13 17:59:48 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+begin_src sh :results output :exports both
|
2021-01-13 17:59:48 +01:00
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 vmc_metropolis.f90 -o vmc_metropolis
|
|
|
|
./vmc_metropolis
|
2021-01-26 13:11:57 +01:00
|
|
|
#+end_src
|
2021-01-13 17:59:48 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49499990423528023 +/- 1.5958250761863871E-004
|
|
|
|
: A = 0.78861366666666655 +/- 3.5096729498002445E-004
|
2021-01-13 17:59:48 +01:00
|
|
|
|
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
** TODO H_2
|
2021-01-13 17:59:48 +01:00
|
|
|
|
|
|
|
We will now consider the H_2 molecule in a minimal basis composed of the
|
|
|
|
$1s$ orbitals of the hydrogen atoms:
|
|
|
|
|
|
|
|
$$
|
|
|
|
\Psi(\mathbf{r}_1, \mathbf{r}_2) =
|
|
|
|
\exp(-(\mathbf{r}_1 - \mathbf{R}_A)) +
|
|
|
|
$$
|
|
|
|
where $\mathbf{r}_1$ and $\mathbf{r}_2$ denote the electron
|
|
|
|
coordinates and $\mathbf{R}_A$ and $\mathbf{R}_B$ the coordinates of
|
|
|
|
the nuclei.
|
2021-01-03 18:45:58 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
|
2021-01-28 15:04:57 +01:00
|
|
|
* Old sections to be removed :noexport:
|
2021-01-20 18:31:49 +01:00
|
|
|
|
|
|
|
** Gaussian sampling :noexport:
|
|
|
|
:PROPERTIES:
|
|
|
|
:header-args:python: :tangle qmc_gaussian.py
|
|
|
|
:header-args:f90: :tangle qmc_gaussian.f90
|
|
|
|
:END:
|
|
|
|
|
|
|
|
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.
|
|
|
|
|
|
|
|
Instead of drawing uniform random numbers, we will draw Gaussian
|
|
|
|
random numbers centered on 0 and with a variance of 1.
|
2021-01-20 19:12:05 +01:00
|
|
|
|
2021-01-20 18:31:49 +01:00
|
|
|
Now the sampling probability can be inserted into the equation of the energy:
|
|
|
|
|
|
|
|
\[
|
|
|
|
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}}
|
|
|
|
\]
|
|
|
|
|
|
|
|
with the Gaussian probability
|
|
|
|
|
|
|
|
\[
|
|
|
|
P(\mathbf{r}) = \frac{1}{(2 \pi)^{3/2}}\exp\left( -\frac{\mathbf{r}^2}{2} \right).
|
|
|
|
\]
|
|
|
|
|
|
|
|
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}
|
|
|
|
$$
|
|
|
|
|
|
|
|
|
|
|
|
*** Exercise
|
|
|
|
|
|
|
|
#+begin_exercise
|
|
|
|
Modify the program of the previous section to sample with
|
|
|
|
Gaussian-distributed random numbers. Can you see an reduction in
|
|
|
|
the statistical error?
|
|
|
|
#+end_exercise
|
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
**** Python
|
|
|
|
#+BEGIN_SRC python :results output
|
2021-01-20 18:31:49 +01:00
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
|
|
|
norm_gauss = 1./(2.*np.pi)**(1.5)
|
|
|
|
def gaussian(r):
|
|
|
|
return norm_gauss * np.exp(-np.dot(r,r)*0.5)
|
|
|
|
|
|
|
|
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
|
|
|
|
|
|
|
|
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-26 13:11:57 +01:00
|
|
|
#+END_SRC
|
2021-01-20 18:31:49 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49511014287471955 +/- 0.00012402022172236656
|
2021-01-20 18:31:49 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
**** Fortran
|
|
|
|
#+BEGIN_SRC f90
|
2021-01-20 18:31:49 +01:00
|
|
|
double precision function gaussian(r)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: r(3)
|
|
|
|
double precision, parameter :: norm_gauss = 1.d0/(2.d0*dacos(-1.d0))**(1.5d0)
|
|
|
|
gaussian = norm_gauss * dexp( -0.5d0 * (r(1)*r(1) + r(2)*r(2) + r(3)*r(3) ))
|
|
|
|
end function gaussian
|
|
|
|
|
|
|
|
|
|
|
|
subroutine gaussian_montecarlo(a,nmax,energy)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a
|
|
|
|
integer*8 , intent(in) :: nmax
|
|
|
|
double precision, intent(out) :: energy
|
|
|
|
|
|
|
|
integer*8 :: istep
|
|
|
|
|
|
|
|
double precision :: norm, r(3), w
|
|
|
|
|
|
|
|
double precision, external :: e_loc, psi, gaussian
|
|
|
|
|
|
|
|
energy = 0.d0
|
|
|
|
norm = 0.d0
|
|
|
|
do istep = 1,nmax
|
|
|
|
call random_gauss(r,3)
|
|
|
|
w = psi(a,r)
|
|
|
|
w = w*w / gaussian(r)
|
|
|
|
norm = norm + w
|
|
|
|
energy = energy + w * e_loc(a,r)
|
|
|
|
end do
|
|
|
|
energy = energy / norm
|
|
|
|
end subroutine gaussian_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9
|
|
|
|
integer*8 , parameter :: nmax = 100000
|
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
|
|
|
double precision :: X(nruns)
|
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
do irun=1,nruns
|
|
|
|
call gaussian_montecarlo(a,nmax,X(irun))
|
|
|
|
enddo
|
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
|
|
|
end program qmc
|
2021-01-26 13:11:57 +01:00
|
|
|
#+END_SRC
|
2021-01-20 18:31:49 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+begin_src sh :results output :exports both
|
2021-01-20 18:31:49 +01:00
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 qmc_gaussian.f90 -o qmc_gaussian
|
|
|
|
./qmc_gaussian
|
2021-01-26 13:11:57 +01:00
|
|
|
#+end_src
|
2021-01-20 18:31:49 +01:00
|
|
|
|
2021-01-26 13:11:57 +01:00
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.49517104619091717 +/- 1.0685523607878961E-004
|
2021-01-20 19:12:05 +01:00
|
|
|
|
|
|
|
** Improved sampling with $\Psi^2$ :noexport:
|
|
|
|
|
|
|
|
*** Importance sampling
|
|
|
|
:PROPERTIES:
|
|
|
|
:header-args:python: :tangle vmc.py
|
|
|
|
:header-args:f90: :tangle vmc.f90
|
|
|
|
:END:
|
|
|
|
|
|
|
|
To generate the probability density $\Psi^2$, we consider a
|
|
|
|
diffusion process characterized by a time-dependent density
|
|
|
|
$[\Psi(\mathbf{r},t)]^2$, which obeys the Fokker-Planck equation:
|
|
|
|
|
|
|
|
\[
|
|
|
|
\frac{\partial \Psi^2}{\partial t} = \sum_i D
|
|
|
|
\frac{\partial}{\partial \mathbf{r}_i} \left(
|
|
|
|
\frac{\partial}{\partial \mathbf{r}_i} - F_i(\mathbf{r}) \right)
|
|
|
|
[\Psi(\mathbf{r},t)]^2.
|
|
|
|
\]
|
|
|
|
|
|
|
|
$D$ is the diffusion constant and $F_i$ is the i-th component of a
|
|
|
|
drift velocity caused by an external potential. For a stationary
|
|
|
|
density, \( \frac{\partial \Psi^2}{\partial t} = 0 \), so
|
|
|
|
|
|
|
|
\begin{eqnarray*}
|
|
|
|
0 & = & \sum_i D
|
|
|
|
\frac{\partial}{\partial \mathbf{r}_i} \left(
|
|
|
|
\frac{\partial}{\partial \mathbf{r}_i} - F_i(\mathbf{r}) \right)
|
|
|
|
[\Psi(\mathbf{r})]^2 \\
|
|
|
|
0 & = & \sum_i D
|
|
|
|
\frac{\partial}{\partial \mathbf{r}_i} \left(
|
|
|
|
\frac{\partial [\Psi(\mathbf{r})]^2}{\partial \mathbf{r}_i} -
|
|
|
|
F_i(\mathbf{r})\,[\Psi(\mathbf{r})]^2 \right) \\
|
|
|
|
0 & = &
|
|
|
|
\frac{\partial^2 \Psi^2}{\partial \mathbf{r}_i^2} -
|
|
|
|
\frac{\partial F_i }{\partial \mathbf{r}_i}[\Psi(\mathbf{r})]^2 -
|
|
|
|
\frac{\partial \Psi^2}{\partial \mathbf{r}_i} F_i(\mathbf{r})
|
|
|
|
\end{eqnarray*}
|
|
|
|
|
|
|
|
we search for a drift function which satisfies
|
|
|
|
|
|
|
|
\[
|
|
|
|
\frac{\partial^2 \Psi^2}{\partial \mathbf{r}_i^2} =
|
|
|
|
[\Psi(\mathbf{r})]^2 \frac{\partial F_i }{\partial \mathbf{r}_i} +
|
|
|
|
\frac{\partial \Psi^2}{\partial \mathbf{r}_i} F_i(\mathbf{r})
|
|
|
|
\]
|
|
|
|
|
|
|
|
to obtain a second derivative on the left, we need the drift to be
|
|
|
|
of the form
|
|
|
|
\[
|
|
|
|
F_i(\mathbf{r}) = g(\mathbf{r}) \frac{\partial \Psi^2}{\partial \mathbf{r}_i}
|
|
|
|
\]
|
|
|
|
|
|
|
|
\[
|
|
|
|
\frac{\partial^2 \Psi^2}{\partial \mathbf{r}_i^2} =
|
|
|
|
[\Psi(\mathbf{r})]^2 \frac{\partial
|
|
|
|
g(\mathbf{r})}{\partial \mathbf{r}_i}\frac{\partial \Psi^2}{\partial \mathbf{r}_i} +
|
|
|
|
[\Psi(\mathbf{r})]^2 g(\mathbf{r}) \frac{\partial^2
|
|
|
|
\Psi^2}{\partial \mathbf{r}_i^2} +
|
|
|
|
\frac{\partial \Psi^2}{\partial \mathbf{r}_i}
|
|
|
|
g(\mathbf{r}) \frac{\partial \Psi^2}{\partial \mathbf{r}_i}
|
|
|
|
\]
|
|
|
|
|
|
|
|
$g = 1 / \Psi^2$ satisfies this equation, so
|
|
|
|
|
|
|
|
\[
|
|
|
|
F(\mathbf{r}) = \frac{\nabla [\Psi(\mathbf{r})]^2}{[\Psi(\mathbf{r})]^2} = 2 \frac{\nabla
|
|
|
|
\Psi(\mathbf{r})}{\Psi(\mathbf{r})} = 2 \nabla \left( \log \Psi(\mathbf{r}) \right)
|
|
|
|
\]
|
|
|
|
|
|
|
|
In statistical mechanics, Fokker-Planck trajectories are generated
|
|
|
|
by a Langevin equation:
|
|
|
|
|
|
|
|
\[
|
|
|
|
\frac{\partial \mathbf{r}(t)}{\partial t} = 2D \frac{\nabla
|
|
|
|
\Psi(\mathbf{r}(t))}{\Psi} + \eta
|
|
|
|
\]
|
|
|
|
|
|
|
|
where $\eta$ is a normally-distributed fluctuating random force.
|
|
|
|
|
|
|
|
Discretizing this differential equation gives the following drifted
|
|
|
|
diffusion scheme:
|
|
|
|
|
|
|
|
\[
|
|
|
|
\mathbf{r}_{n+1} = \mathbf{r}_{n} + \tau\, 2D \frac{\nabla
|
|
|
|
\Psi(\mathbf{r})}{\Psi(\mathbf{r})} + \chi
|
|
|
|
\]
|
|
|
|
where $\chi$ is a Gaussian random variable with zero mean and
|
|
|
|
variance $\tau\,2D$.
|
|
|
|
|
|
|
|
**** Exercise 2
|
|
|
|
|
|
|
|
#+begin_exercise
|
|
|
|
Sample $\Psi^2$ approximately using the drifted diffusion scheme,
|
|
|
|
with a diffusion constant $D=1/2$. You can use a time step of
|
|
|
|
0.001 a.u.
|
|
|
|
#+end_exercise
|
|
|
|
|
|
|
|
*Python*
|
|
|
|
#+BEGIN_SRC python :results output
|
|
|
|
from hydrogen import *
|
|
|
|
from qmc_stats import *
|
|
|
|
|
|
|
|
def MonteCarlo(a,tau,nmax):
|
|
|
|
sq_tau = np.sqrt(tau)
|
|
|
|
|
|
|
|
# Initialization
|
|
|
|
E = 0.
|
|
|
|
N = 0.
|
|
|
|
r_old = np.random.normal(loc=0., scale=1.0, size=(3))
|
|
|
|
|
|
|
|
for istep in range(nmax):
|
|
|
|
d_old = drift(a,r_old)
|
|
|
|
chi = np.random.normal(loc=0., scale=1.0, size=(3))
|
|
|
|
r_new = r_old + tau * d_old + chi*sq_tau
|
|
|
|
N += 1.
|
|
|
|
E += e_loc(a,r_new)
|
|
|
|
r_old = r_new
|
|
|
|
return E/N
|
|
|
|
|
|
|
|
|
|
|
|
a = 0.9
|
|
|
|
nmax = 100000
|
|
|
|
tau = 0.2
|
|
|
|
X = [MonteCarlo(a,tau,nmax) for i in range(30)]
|
|
|
|
E, deltaE = ave_error(X)
|
|
|
|
print(f"E = {E} +/- {deltaE}")
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.4858534479298907 +/- 0.00010203236131158794
|
|
|
|
|
|
|
|
*Fortran*
|
|
|
|
#+BEGIN_SRC f90
|
|
|
|
subroutine variational_montecarlo(a,tau,nmax,energy)
|
|
|
|
implicit none
|
|
|
|
double precision, intent(in) :: a, tau
|
|
|
|
integer*8 , intent(in) :: nmax
|
|
|
|
double precision, intent(out) :: energy
|
|
|
|
|
|
|
|
integer*8 :: istep
|
|
|
|
double precision :: norm, r_old(3), r_new(3), d_old(3), sq_tau, chi(3)
|
|
|
|
double precision, external :: e_loc
|
|
|
|
|
|
|
|
sq_tau = dsqrt(tau)
|
|
|
|
|
|
|
|
! Initialization
|
|
|
|
energy = 0.d0
|
|
|
|
norm = 0.d0
|
|
|
|
call random_gauss(r_old,3)
|
|
|
|
|
|
|
|
do istep = 1,nmax
|
|
|
|
call drift(a,r_old,d_old)
|
|
|
|
call random_gauss(chi,3)
|
|
|
|
r_new(:) = r_old(:) + tau * d_old(:) + chi(:)*sq_tau
|
|
|
|
norm = norm + 1.d0
|
|
|
|
energy = energy + e_loc(a,r_new)
|
|
|
|
r_old(:) = r_new(:)
|
|
|
|
end do
|
|
|
|
energy = energy / norm
|
|
|
|
end subroutine variational_montecarlo
|
|
|
|
|
|
|
|
program qmc
|
|
|
|
implicit none
|
|
|
|
double precision, parameter :: a = 0.9
|
|
|
|
double precision, parameter :: tau = 0.2
|
|
|
|
integer*8 , parameter :: nmax = 100000
|
|
|
|
integer , parameter :: nruns = 30
|
|
|
|
|
|
|
|
integer :: irun
|
|
|
|
double precision :: X(nruns)
|
|
|
|
double precision :: ave, err
|
|
|
|
|
|
|
|
do irun=1,nruns
|
|
|
|
call variational_montecarlo(a,tau,nmax,X(irun))
|
|
|
|
enddo
|
|
|
|
call ave_error(X,nruns,ave,err)
|
|
|
|
print *, 'E = ', ave, '+/-', err
|
|
|
|
end program qmc
|
|
|
|
#+END_SRC
|
|
|
|
|
|
|
|
#+begin_src sh :results output :exports both
|
|
|
|
gfortran hydrogen.f90 qmc_stats.f90 vmc.f90 -o vmc
|
|
|
|
./vmc
|
|
|
|
#+end_src
|
|
|
|
|
|
|
|
#+RESULTS:
|
|
|
|
: E = -0.48584030499187431 +/- 1.0411743995438257E-004
|
|
|
|
|
2021-01-25 23:52:53 +01:00
|
|
|
|
2021-01-28 16:02:19 +01:00
|
|
|
* TODO [0/3] Last things to do
|
2021-01-25 23:52:53 +01:00
|
|
|
|
2021-01-28 16:02:19 +01:00
|
|
|
- [ ] Give some hints of how much time is required for each section
|
2021-01-25 23:52:53 +01:00
|
|
|
- [ ] Prepare 4 questions for the exam: multiple-choice questions
|
|
|
|
with 4 possible answers. Questions should be independent because
|
|
|
|
they will be asked in a random order.
|
|
|
|
- [ ] Propose a project for the students to continue the
|
|
|
|
programs. Idea: Modify the program to compute the exact energy of
|
|
|
|
the H$_2$ molecule at $R$=1.4010 bohr. Answer: 0.17406 a.u.
|