Quantum Monte Carlo
Top: Quantum Monte CarloHELP / toggle view | ||
Previous | Next | Table of Contents | 1 |
Quantum Monte Carlo
Table of Contents
Quantum Monte Carlo
Table of Contents
1. Introduction
This website contains the QMC tutorial of the 2023 LTTC winter school Tutorials in Theoretical Chemistry.
We propose different exercises to understand quantum Monte Carlo (QMC) methods. In the first section, we start with the computation of the energy of a hydrogen atom using numerical integration. The goal of this section is to familiarize yourself with the concept of 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, and apply this approach to the hydrogen atom. Finally, we present the diffusion Monte Carlo (DMC) method which we use here to estimate the exact energy of the hydrogen atom and of the H2 molecule, starting from an approximate wave function.
Code examples will be given in Python3, C and Fortran. You can use whatever language you prefer to write the programs.
We consider the stationary solution of the Schrödinger equation, so
the wave functions considered here are real: for an
All the quantities are expressed in atomic units (energies, coordinates, etc).
1.1. Energy and local energy
For a given system with Hamiltonian
where
The electronic energy of a system,
For few dimensions, one can easily compute
To this aim, recall that the probabilistic expected value of an
arbitrary function
where the probability density function
Similarly, we can view the energy of a system,
where the probability density is given by the square of the wave function:
If we can sample
Here is an animation of the local energy in the N2 molecule along a QMC trajectory:
2. Numerical evaluation of the energy of the hydrogen atom
In this section, we consider the hydrogen atom with the following wave function:
We will first verify that, for a particular value of
To do that, we will compute the local energy and check whether it is constant.
2.1. Local energy
You will now program all quantities needed to compute the local energy of the Hydrogen atom for the given wave function.
Write all the functions of this section in a single file :
hydrogen.py
if you use Python, hydrogen.F90
if you use
Fortran, or hydrogen.c
and hydrogen.h
if you use C.
- When computing a square root in
, 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.
2.1.1. Test-driven development
Test-driven development (TDD) is a software development process in which automated tests are written before the code they are testing. The tests are then run to ensure that the code is working as expected, and the developer continues to write code until all tests pass. The benefits of TDD include:
- Early detection of bugs: By writing tests before the code, developers can identify and fix bugs early on in the development process, before the code becomes too complex to easily debug.
- Improved code quality: Writing tests forces developers to think about the behavior and functionality of the code they are writing, leading to more robust and maintainable code.
- Better documentation: Tests act as a form of documentation for the code, making it clear what the code is supposed to do and how it should behave.
2.1.2. Exercise 1
Write a function called potential
that takes a single argument
r
, which is an array of 3 double precision numbers. The function
calculates the distance from the origin using the Euclidean
distance formula, and then calculates the potential using this
distance.
If the distance is 0, the program prints a warning that the
potential diverges.
If the distance is greater than 0, the potential is calculated as
- Python
#!/usr/bin/env python3 import numpy as np def potential(r): # TODO def test_potential(): expected_output = -1./15. for r in [( 2., 5., 14.), (5., 14., 2.), (-2., 5.,-14.), (5.,-14.,-2.), ( 0., 9.,-12.), (9.,-12., 0.)]: assert potential(r) == expected_output r = (0., 0., 0.) assert potential(r) == -float("inf") print("potential ok") if __name__ == "__main__": test_potential()
- Fortran
If the extension of your file is
.F90
and not.f90
, the C preprocessor will be called before compiling. This enables the possibility to have conditional compilation with#ifdef
statements, activated with the-D
compiler option.To compile your code and activate the test, use:
gfortran -DTEST_H hydrogen.F90 -o test_hydrogen
double precision function potential(r) implicit none double precision, intent(in) :: r(3) ! TODO end function potential subroutine test_potential implicit none double precision :: r(3) double precision :: expected_output double precision, external :: potential expected_output = -1.d0/15.d0 r(:) = (/ 2.d0, 5.d0, 14.d0 /) if (potential(r) /= expected_output) stop 'Failed' r(:) = (/ 5.d0, 14.d0, 2.d0 /) if (potential(r) /= expected_output) stop 'Failed' r(:) = (/ -2.d0, 5.d0, -14.d0 /) if (potential(r) /= expected_output) stop 'Failed' r(:) = (/ 5.d0, -14.d0, -2.d0 /) if (potential(r) /= expected_output) stop 'Failed' r(:) = (/ 0.d0, 9.d0, 12.d0 /) if (potential(r) /= expected_output) stop 'Failed' r(:) = (/ 9.d0, -12.d0, 0.d0 /) if (potential(r) /= expected_output) stop 'Failed' r(:) = 0.d0 expected_output = -huge(1.d0) if (potential(r) /= expected_output) stop 'Failed r=0' print *, 'potential ok' end subroutine test_potential #ifdef TEST_H program test_h call test_potential end program test_h #endif
- C
Compile your code with the
-DTEST_H
option. It will activate the creation of the main function that will test your functions. Don't forget to use-lm
to link with the math library.gcc -DTEST_H hydrogen.c -lm -o test_hydrogen
#include <stdio.h> // printf #include <math.h> // sqrt #include <stdlib.h> // exit #include <assert.h> // assert double potential(double r[3]) { // TODO } static void test_potential() { double r[3]; double expected_output; expected_output = -1./15.; r[0] = 2.; r[1] = 5.; r[2] = 14.; assert (potential(r) == expected_output); r[0] = 5.; r[1] = 14.; r[2] = 2.; assert (potential(r) == expected_output); r[0] = -2.; r[1] = 5.; r[2] = -14.; assert (potential(r) == expected_output); r[0] = 5.; r[1] = -14.; r[2] = -2.; assert (potential(r) == expected_output); r[0] = 0.; r[1] = 9.; r[2] = 12.; assert (potential(r) == expected_output); r[0] = 9.; r[1] = -12.; r[2] = 0.; assert (potential(r) == expected_output); expected_output = -HUGE_VAL; r[0] = 0.; r[1] = 0.; r[2] = 0.; assert (potential(r) == expected_output); printf("potential ok\n"); } #ifdef TEST_H int main() { test_potential(); return 0; } #endif
2.1.3. Exercise 2
Write a function called psi
that takes two arguments: a
which is a
double precision number, and r
which is an array of 3 double
precision numbers. The function calculates and returns the value of the wave
function at the given point
2.1.4. Exercise 3
Write a function called kinetic
that takes two arguments: a
which is a double precision number, and r
which is an array of 3
double precision numbers. The function calculates the local kinetic energy
at the given point
It first calculates the distance from the origin using the Euclidean distance formula. If the distance is greater than 0, the kinetic energy is calculated using the formula given below. If the distance is 0, the program prints a warning.
The local kinetic energy is defined as
We differentiate
and we differentiate a second time:
The Laplacian operator
Therefore, the local kinetic energy is
2.1.5. Exercise 4
Write a function called e_loc
that takes two arguments: a
which is a double precision number and r
which is an array of 3
double precision numbers. The function calculates the local energy at
the given point
It uses two functions kinetic
and potential
to calculate the
kinetic energy and potential energy respectively and add them to
get the local energy.
- Python
def e_loc(a,r): #TODO
- Fortran
When you call a function in Fortran, you need to declare its return type. You might by accident choose a function name which is the same as an internal function of Fortran. So it is recommended to always use the keyword
external
to make sure the function you are calling is yours.double precision function e_loc(a,r) implicit none double precision, intent(in) :: a, r(3) double precision, external :: kinetic double precision, external :: potential ! TODO end function e_loc
- C
double e_loc(double a, double r[3]) { // TODO }
2.1.6. Exercise 5
Find the theoretical value of
2.2. Plot of the local energy along the axis
The program you will write in this section will be written in
another file (plot_hydrogen.py
, plot_hydrogen.F90
or plot_hydrogen.c
for
example).
It will use the functions previously defined.
If you use C, don't forget to write the header file corresponding
to the functions defined in the previous section.
In Python, you should put at the beginning of the file
#!/usr/bin/env python3 from hydrogen import e_loc
to be able to use the e_loc
function of the hydrogen.py
file.
It is better to use #!/usr/bin/env python3
than
#!/usr/bin/python
because:
- you are sure you are not using Python2, which is incompatible with Python3 syntax,
- if you are on a machine where you can load different
environments (VirtualEnv, module, etc), you will use the
python3
provided by your environment, and not the system's one.
In Fortran, you will need to compile all the source files together:
gfortran hydrogen.F90 plot_hydrogen.F90 -o plot_hydrogen
Similarly, in C
gcc hydrogen.c plot_hydrogen.c -lm -o plot_hydrogen
In C, you need the -lm
argument to link with the math library
that contains functions like sqrt
and exp
.
2.2.1. Exercise
For multiple values of
In Python, you can use matplotlib for example.
In Fortran, it is convenient to write in a text file
the values of
The potential and the kinetic energy both diverge at
- Python
#!/usr/bin/env python3 import numpy as np import matplotlib.pyplot as plt from hydrogen import e_loc x=np.linspace(-5,5) plt.figure(figsize=(10,5)) # TODO plt.tight_layout() plt.legend() plt.savefig("plot_py.png")
- Fortran
program plot implicit none double precision, external :: e_loc 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
To compile and run:
gfortran hydrogen.F90 plot_hydrogen.F90 -o plot_hydrogen ./plot_hydrogen > data
- C
#include <stdio.h> #include <math.h> #include "hydrogen.h" #define NPOINTS 50 #define NEXPO 6 int main() { double x[NPOINTS], energy, dx, r[3]; double a[NEXPO] = { 0.1, 0.2, 0.5, 1.0, 1.5, 2.0 }; int i, j; dx = 10.0/(NPOINTS-1); for (i = 0; i < NPOINTS; i++) { x[i] = -5.0 + i*dx; } // TODO return 0; }
To compile and run:
gcc hydrogen.c plot_hydrogen.c -lm -o plot_hydrogen ./plot_hydrogen > data
Plotting for Fortran of C
Plot the data using Gnuplot:
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'
2.3. Numerical estimation of the energy
If the space is discretized in small volume elements
The energy is biased because:
- The volume elements are not infinitely small (discretization error)
- The energy is evaluated only inside the box (incompleteness of the space)
2.3.1. Exercise
Compute a numerical estimate of the energy using a grid of
- Python
#!/usr/bin/env python3 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}")
- Fortran
program energy_hydrogen implicit none double precision, external :: e_loc, psi double precision :: x(50), w, delta, energy, dx, r(3), a(6), normalization 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
To compile the Fortran code and run it:
gfortran hydrogen.F90 energy_hydrogen.F90 -o energy_hydrogen ./energy_hydrogen
- C
#include <stdio.h> #include <math.h> #include "hydrogen.h" #define NPOINTS 50 #define NEXPO 6 int main() { double x[NPOINTS], energy, dx, r[3], delta, normalization, w; double a[NEXPO] = { 0.1, 0.2, 0.5, 1.0, 1.5, 2.0 }; dx = 10.0/(NPOINTS-1); for (int i = 0; i < NPOINTS; i++) { x[i] = -5.0 + i*dx; } for (int j = 0; j < NEXPO; j++) { // TODO printf("a = %f E = %f\n", a[j], energy); } }
To compile the C code and run it:
gcc hydrogen.c energy_hydrogen.c -lm -o energy_hydrogen ./energy_hydrogen
Hints if you are stuck
The program starts by defining some variables and arrays, including an array
a
that contains 6 different values of the parametera
which will be used in thee_loc
andpsi
functions to calculate the local energy and wave function respectively.The program then calculates the value of
dx
, which is the step size in , and sets up an arrayx
that contains 50 equally spaced points between -5 and 5. The program sets all elements of ther
array to 0, and then enters a nested loop structure. The outer loop iterates over the values ofa
in thea
array, and the next three loops iterate over the values ofx
in thex
array for the three dimensions. For each value ofa
andx
, the program sets the first element of ther
array to the current value ofx
, calls thepsi
function to calculate the wave function, calls thee_loc
function to calculate the local energy, and then accumulates the energy and the normalization factor.At the end of the outer loop, the program calculates the final energy by dividing the accumulated energy by the accumulated normalization factor, and prints the value of
a
and the corresponding energy.
2.4. Variance of the local energy
The variance of the local energy is a functional of
If the local energy is constant (i.e.
2.4.2. Exercise
Add the calculation of the variance to the previous code, and
compute a numerical estimate of the variance of the local energy using
a grid of
- Python
#!/usr/bin/env python3 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:10.8f} \t \sigma^2 = {s2:10.8f}")
- Fortran
program variance_hydrogen implicit none double precision :: x(50), w, delta, energy, energy2 double precision :: dx, r(3), a(6), normalization, e_tmp, s2 integer :: i, k, l, j double precision, external :: e_loc, psi 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, ' s2 = ', s2 end do end program variance_hydrogen
To compile and run:
gfortran hydrogen.F90 variance_hydrogen.F90 -o variance_hydrogen ./variance_hydrogen
- C
#include <stdio.h> #include <math.h> #include "hydrogen.h" #define NPOINTS 50 #define NEXPO 6 int main() { double x[NPOINTS], energy, dx, r[3], delta, normalization, w; double a[NEXPO] = { 0.1, 0.2, 0.5, 1.0, 1.5, 2.0 }; double energy2, e_tmp, s2; dx = 10.0/(NPOINTS-1); for (int i = 0; i < NPOINTS; i++) { x[i] = -5.0 + i*dx; } for (int j = 0; j < NEXPO; j++) { // TODO printf("a = %f E = %f s2 = %f\n", a[j], energy, s2); } }
To compile and run:
gcc hydrogen.c variance_hydrogen.c -lm -o variance_hydrogen ./variance_hydrogen
3. Variational Monte Carlo
Numerical integration with deterministic methods is very efficient 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 use Monte Carlo sampling.
Moreover, Monte Carlo sampling will allow us to remove the bias due to the discretization of space, and compute a statistical confidence interval.
3.1. Computation of the statistical error
To compute the statistical error, you need to perform
The estimate of the energy is
The variance of the average energies can be computed as
And the confidence interval is given by
3.1.1. Exercise
Write a function returning the average and statistical error of an input array.
- Python
#!/usr/bin/env python3 from math import sqrt def ave_error(arr): #TODO return (average, error)
- Fortran
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
- C
#include <stdio.h> #include <math.h> #include <stddef.h> // for size_t void ave_error(double* x, size_t n, double *ave, double *err) { // TODO }
Hints if you are stuck
Write a subroutine called
ave_error
that calculates the average and error of a given array of real numbers. The subroutine takes in three arguments: an arrayx
of real numbers, an integern
representing the size of the array, and two output argumentsave
anderr
representing the average and error of the array, respectively.The subroutine starts by checking if the input integer
n
is less than 1. If it is, the subroutine stops and prints an error message. Ifn
is equal to 1, the subroutine sets the average to the first element of the array and the error to zero. Ifn
is greater than 1, the subroutine calculates the average of the array by dividing the sum of the elements by the number of elements in the array. Then it calculates the variance of the array by taking the sum of the square of the difference between each element and the average and dividing byn-1
. Finally, it calculates the error by taking the square root of the variance divided byn
.
3.2. Uniform sampling in the box
We will now perform our first Monte Carlo calculation to compute the energy of the hydrogen atom.
Consider again the expression of the energy
Clearly, the square of the wave function is a good choice of probability density to sample but we will start with something simpler and rewrite the energy as
Here, we will sample a uniform probability
and zero outside the cube.
One Monte Carlo run will consist of
- Draw a random point
in the box - Compute
and accumulate the result in a variablenormalization
- Compute
, and accumulate the result in a variableenergy
Once all the iterations have been computed, the run returns the average energy
To compute the statistical error, perform
3.2.1. Exercise
Parameterize the wave function with
- Python
To draw a uniform random number in Python, you can use the
random.uniform
function of Numpy.#!/usr/bin/env python3 from hydrogen import * from qmc_stats import * def MonteCarlo(a, nmax): # TODO a = 1.2 nmax = 100000 #TODO print(f"E = {E} +/- {deltaE}")
- Fortran
To draw a uniform random number in Fortran, you can use the
RANDOM_NUMBER
subroutine.When running Monte Carlo calculations, the number of steps is usually very large. We expect
nmax
to be possibly larger than 2 billion. You would need to use 8-byte integers (integer*8
) to represent it, as well as the index of the current step. This would imply modifying also theave_error
function.subroutine uniform_montecarlo(a,nmax,energy) implicit none double precision, intent(in) :: a integer , intent(in) :: nmax double precision, intent(out) :: energy integer :: istep double precision :: normalization, r(3), w double precision, external :: e_loc, psi ! TODO end subroutine uniform_montecarlo program qmc implicit none double precision, parameter :: a = 1.2d0 integer , 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
gfortran hydrogen.F90 qmc_stats.F90 qmc_uniform.F90 -o qmc_uniform ./qmc_uniform
- C
To draw a uniform random number in C, you can use:
drand48()
, which is defined in thestdlib.h
header. To initialize randomly the generator, usesrand48(time(NULL))
using thetime
function fromtime.h
.#include <stdlib.h> #include <math.h> #include <stdio.h> #include <stddef.h> // for size_t #include <time.h> #include "hydrogen.h" #include "qmc_stats.h" // for ave_error void uniform_montecarlo(double a, size_t nmax, double *energy) { // TODO } int main(void) { #define a 1.2 #define nmax 100000 #define nruns 30 srand48(time(NULL)); // TODO printf("E = %f +/- %f\n", ave, err); return 0; }
Hints if you are stuck
Write first a subroutine called
uniform_montecarlo
that calculates the energy of the Hydrogen atom using the Monte Carlo method with a uniform distribution. The subroutine takes in three arguments: a real numbera
, an integernmax
representing the number of Monte Carlo steps, and an output argumentenergy
representing the calculated energy.The subroutine starts by initializing the energy and normalization factor to 0 and defines some variables such as
istep
,normalization
,r
andw
. The subroutine also makes use of two external functions:e_loc
andpsi
which were defined in previous examples.The subroutine then enters a loop that iterates for
nmax
times. On each iteration, the subroutine generates three random numbers between 0 and 1, and then uses these random numbers to calculate a random point in 3D space between -5 and 5. The subroutine then calls thepsi
function to calculate the wave function at that point and thee_loc
function to calculate the local energy at that point. The subroutine then accumulates the energy and normalization factor using the generated point and the results of thepsi
ande_loc
functions.At the end of the loop, the subroutine calculates the final energy by dividing the accumulated energy by the accumulated normalization factor.
Then, write a Fortran program called
qmc
that uses theuniform_montecarlo
subroutine to estimate the energy of the Hydrogen atom using the Monte Carlo method. The program starts by defining some parameters:a
,nmax
, andnruns
.The program then defines a variable
irun
which is used as a counter in a loop, an arrayX
of lengthnruns
to store the energies calculated by theuniform_montecarlo
subroutine, and variablesave
anderr
to store the average and error of the energies, respectively.The program then enters a loop that iterates for
nruns
times. On each iteration, the program calls theuniform_montecarlo
subroutine to calculate the energy of the Hydrogen atom and stores the result in theX
array.After the loop, the program calls the
ave_error
subroutine to calculate the average and error of the energies stored in theX
array and assigns the results toave
anderr
variables respectively.Finally, the program prints the average and error of the energies.
3.3. Metropolis sampling with Ψ2
We will now use the square of the wave function to sample random
points distributed with the probability density
The expression of the average energy is now simplified as the average of the local energies, since the weights are taken care of by the sampling:
To sample a chosen probability density, an efficient method is the
Metropolis-Hastings sampling algorithm. Starting from a random
initial position
according to the following algorithm.
At every step, we propose a new move according to a transition probability
For simplicity, we will move the electron in a 3-dimensional box of side
where
After having moved the electron, we add the
accept/reject step that guarantees that the distribution of the
which, for our choice of transition probability, becomes
Explain why the transition probability cancels out in the
expression of
Also note that we do not need to compute the norm of the wave function!
The algorithm is summarized as follows:
- Evaluate the local energy at
and add it to an accumulator - Compute a new position
- Evaluate
at the new position - Compute the ratio
- Draw a uniform random number
- if
, accept the move : set - else, reject the move : set
A common error is to remove the rejected samples from the calculation of the average. Don't do it!
All samples should be kept, from both accepted and rejected moves.
3.3.1. Optimal step size
If the box is infinitely small, the ratio will be very close to one and all the steps will be accepted. However, the moves will be very correlated and you will explore the configurational space very slowly.
On the other hand, if you propose too large moves, 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 size of the move should be adjusted so that it is as large as possible, keeping the number of accepted steps not too small. To achieve that, we define the acceptance rate as the number of 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 for the current problem.
Below, we use the symbol
3.3.2. Exercise
Modify the program of the previous section to compute the energy,
sampled with
Compute also the acceptance rate, so that you can adapt the time step in order to have an acceptance rate close to 0.5.
Can you observe a reduction in the statistical error?
- Python
#!/usr/bin/env python3 from hydrogen import * from qmc_stats import * def MonteCarlo(a,nmax,dt): # TODO return energy/nmax, N_accep/nmax # Run simulation a = 1.2 nmax = 100000 dt = #TODO X0 = [ MonteCarlo(a,nmax,dt) for i in range(30)] # 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}")
- Fortran
subroutine metropolis_montecarlo(a,nmax,dt,energy,accep) implicit none double precision, intent(in) :: a integer*8 , intent(in) :: nmax double precision, intent(in) :: dt double precision, intent(out) :: energy double precision, intent(out) :: accep integer*8 :: istep integer*8 :: n_accep double precision :: r_old(3), r_new(3), psi_old, psi_new double precision :: v, ratio double precision, external :: e_loc, psi, gaussian ! TODO end subroutine metropolis_montecarlo program qmc implicit none double precision, parameter :: a = 1.2d0 double precision, parameter :: dt = ! TODO integer*8 , parameter :: nmax = 100000 integer , parameter :: nruns = 30 integer :: irun double precision :: X(nruns), Y(nruns) double precision :: ave, err do irun=1,nruns call metropolis_montecarlo(a,nmax,dt,X(irun),Y(irun)) enddo 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
- C
#include <stdio.h> #include <stdlib.h> #include <stddef.h> // for size_t #include <time.h> #include <math.h> #include "hydrogen.h" #include "qmc_stats.h" void metropolis_montecarlo(double a, size_t nmax, double dt, double *energy, double *accep) { // TODO } int main(void) { #define a 1.2 #define nmax 100000 #define dt //TODO #define nruns 30 double energy[nruns]; double accep[nruns]; double ave, err; srand48(time(NULL)); for (size_t irun = 0; irun < nruns; irun++) { metropolis_montecarlo(a, nmax, dt, energy, accep); } ave_error(energy, nruns, &ave, &err); printf("E = %f +/- %f\n", ave, err); ave_error(accep, nruns, &ave, &err); printf("A = %f +/- %f\n", ave, err); return 0; }
3.4. Generalized Metropolis algorithm
One can use more efficient numerical schemes to move the electrons by choosing a smarter expression for the transition probability.
The Metropolis acceptance step has to be adapted keeping in mind that
the detailed balance condition is satisfied. This means that the acceptance
probability
where
In the previous example, we were using uniform sampling in a box centered at the current position. Hence, the transition probability was symmetric
so the expression of
Now, if instead of drawing uniform random numbers, we
choose to draw Gaussian random numbers with zero mean and variance
Furthermore, to sample the density even better, we can "push" the electrons into in the regions of high probability, and "pull" them away from the low-probability regions. This will increase the acceptance ratios and improve the sampling.
To do this, we can use the gradient of the probability density
and add the so-called drift vector,
The corresponding move is proposed as
where
Here is an illustration of a trajectory:
Averaging all the trajectories converges to the density of the
hydrogen atom.
The algorithm of the previous exercise is only slightly modified as:
- Evaluate the local energy at
and accumulate it - Compute a new position
- Evaluate
and at the new position - Compute the ratio
- Draw a uniform random number
- if
, accept the move : set - else, reject the move : set
3.4.1. Gaussian random number generator
To obtain Gaussian-distributed random numbers, you can apply the Box Muller transform to uniform random numbers:
This gives Gaussian-distributed random numbers with zero mean and a variance of one.
Below is a Fortran and a C implementation returning a
Gaussian-distributed n-dimensional vector random.normal
function of Numpy.
If you use Fortran of C, copy/paste the random_gauss
function in
the qmc_stats.F90
or qmc_stats.c
file.
- Fortran
subroutine random_gauss(z,n) implicit none integer, intent(in) :: n double precision, intent(out) :: z(n) double precision :: u(n+1) double precision, parameter :: two_pi = 2.d0*dacos(-1.d0) integer :: i call random_number(u) if (iand(n,1) == 0) then ! n is even do i=1,n,2 z(i) = dsqrt(-2.d0*dlog(u(i))) z(i+1) = z(i) * dsin( two_pi*u(i+1) ) z(i) = z(i) * dcos( two_pi*u(i+1) ) end do else ! n is odd do i=1,n-1,2 z(i) = dsqrt(-2.d0*dlog(u(i))) z(i+1) = z(i) * dsin( two_pi*u(i+1) ) z(i) = z(i) * dcos( two_pi*u(i+1) ) end do z(n) = dsqrt(-2.d0*dlog(u(n))) z(n) = z(n) * dcos( two_pi*u(n+1) ) end if end subroutine random_gauss
- C
void random_gauss(double* z, size_t n);
#include <stdlib.h> void random_gauss(double* z, size_t n) { double u[n+1]; double two_pi = 2.0 * acos(-1.0); size_t i; //generate random numbers for (i = 0; i <= n; i++) { u[i] = drand48(); } if (n % 2 == 0) { // n is even for (i = 0; i < n; i += 2) { z[i] = sqrt(-2.0 * log(u[i])); z[i+1] = z[i] * sin(two_pi * u[i+1]); z[i] = z[i] * cos(two_pi * u[i+1]); } } else { // n is odd for (i = 0; i < n-1; i += 2) { z[i] = sqrt(-2.0 * log(u[i])); z[i+1] = z[i] * sin(two_pi * u[i+1]); z[i] = z[i] * cos(two_pi * u[i+1]); } z[n-1] = sqrt(-2.0 * log(u[n-1])); z[n-1] = z[n-1] * cos(two_pi * u[n]); } }
3.4.2. Exercise 1
Write a function to compute the drift vector
3.4.3. Exercise 2
Modify the previous program to introduce the drift-diffusion scheme. (This is a necessary step for the next section on diffusion Monte Carlo).
- Python
#!/usr/bin/env python3 from hydrogen import * from qmc_stats import * def MonteCarlo(a,nmax,dt): # TODO # Run simulation a = 1.2 nmax = 100000 dt = # TODO X0 = [ MonteCarlo(a,nmax,dt) for i in range(30)] # 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}")
- Fortran
subroutine variational_montecarlo(a,nmax,dt,energy,accep) implicit none double precision, intent(in) :: a, dt integer*8 , intent(in) :: nmax double precision, intent(out) :: energy, accep integer*8 :: istep integer*8 :: n_accep double precision :: sq_dt, chi(3) 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 = 1.2d0 double precision, parameter :: dt = ! TODO 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,nmax,dt,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
gfortran hydrogen.F90 qmc_stats.F90 vmc_metropolis.F90 -o vmc_metropolis ./vmc_metropolis
- C
#include <stdio.h> #include <stdlib.h> #include <stddef.h> // for size_t #include <math.h> #include <time.h> #include "hydrogen.h" #include "qmc_stats.h" void variational_montecarlo(double a, size_t nmax, double dt, double *energy, double *accep) { // TODO } int main(void) { #define a 1.2 #define nmax 100000 #define dt // TODO #define nruns 30 double X[nruns]; double Y[nruns]; double ave, err; srand48(time(NULL)); for (size_t irun = 0; irun < nruns; irun++) { variational_montecarlo(a, nmax, dt, &X[irun], &Y[irun]); } ave_error(X, nruns, &ave, &err); printf("E = %f +/- %f\n", ave, err); ave_error(Y, nruns, &ave, &err); printf("A = %f +/- %f\n", ave, err); return 0; }
4. Diffusion Monte Carlo
As we have seen, Variational Monte Carlo is a powerful method to compute integrals in large dimensions. It is often used in cases where the expression of the wave function is such that the integrals can't be evaluated (multi-centered Slater-type orbitals, correlation factors, etc).
Diffusion Monte Carlo is different. It goes beyond the computation of the integrals associated with an input wave function, and aims at finding a near-exact numerical solution to the Schrödinger equation.
4.1. Schrödinger equation in imaginary time
Consider the time-dependent Schrödinger equation:
where we introduced a shift in the energy,
We can expand a given starting wave function,
The solution of the Schrödinger equation at time
Now, if we replace the time variable
where
For large positive values of
4.2. Relation to diffusion
The diffusion equation of particles is given by
where
it can be identified as the combination of:
- a diffusion equation (Laplacian)
- an equation whose solution is an exponential (potential)
The diffusion equation can be simulated by a Brownian motion:
where
We note that the ground-state wave function of a Fermionic system is
antisymmetric and changes sign. Therefore, its interpretation as a probability
distribution is somewhat problematic. In fact, mathematically, since
the Bosonic ground state is lower in energy than the Fermionic one, for
large
For the systems you will study, this is not an issue:
- Hydrogen atom: You only have one electron!
- Two-electron system (
or He): The ground-wave function is antisymmetric in the spin variables but symmetric in the space ones.
Therefore, in both cases, you are dealing with a "Bosonic" ground state.
4.3. Importance sampling
In a molecular system, the potential is far from being constant
and, in fact, diverges at the inter-particle coalescence points. Hence,
it results in very large fluctuations of the weight associated with
the potential, making the calculations impossible in practice.
Fortunately, if we multiply the Schrödinger equation by a chosen
trial wave function
Defining
The new "kinetic energy" can be simulated by the drift-diffusion
scheme presented in the previous section (VMC).
The new "potential" is the local energy, which has smaller fluctuations
when
This equation generates the N-electron density
To this aim, we use the mixed estimator of the energy:
For large
and, using that
Therefore, we can compute the energy within DMC by generating the
density
4.4. Pure Diffusion Monte Carlo
Instead of having a variable number of particles to simulate the branching process as in the Diffusion Monte Carlo (DMC) algorithm, we use variant called pure Diffusion Monte Carlo (PDMC) where the potential term is considered as a cumulative product of weights:
where
The PDMC algorithm is less stable than the DMC algorithm: it
requires to have a value of
- Start with
- Evaluate the local energy at
- Compute the contribution to the weight
- Update
- Accumulate the weighted energy
, and the weight for the normalization - Update
- If
( is an input parameter), the long projection time has been reached and we can start an new trajectory from the current position: reset and - Compute a new position
- Evaluate
and at the new position - Compute the ratio
- Draw a uniform random number
- if
, accept the move : set - else, reject the move : set
Some comments are needed:
You estimate the energy as
The result will be affected by a time-step error (the finite size of
) due to the discretization of the integral, and one has in principle to extrapolate to the limit . This amounts to fitting the energy computed for multiple values of .Here, you will be using a small enough time-step and you should not worry about the extrapolation.
- The accept/reject step (steps 9-12 in the algorithm) is in principle not needed for the correctness of the DMC algorithm. However, its use reduces significantly the time-step error.
4.5. Hydrogen atom
4.5.1. Exercise
Modify the Metropolis VMC program into a PDMC program.
In the limit
- Python
from hydrogen import * from qmc_stats import * def MonteCarlo(a, nmax, dt, Eref): # TODO # Run simulation a = 1.2 nmax = 100000 dt = 0.05 tau = 100. E_ref = -0.5 X0 = [ MonteCarlo(a, nmax, dt, E_ref) for i in range(30)] # 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}")
- Fortran
subroutine pdmc(a, nmax, dt, energy, accep, tau, E_ref) implicit none double precision, intent(in) :: a, dt, tau integer*8 , intent(in) :: nmax double precision, intent(out) :: energy, accep double precision, intent(in) :: E_ref integer*8 :: istep integer*8 :: n_accep double precision :: sq_dt, chi(3) 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 pdmc program qmc implicit none double precision, parameter :: a = 1.2d0 double precision, parameter :: dt = 0.05d0 double precision, parameter :: E_ref = -0.5d0 double precision, parameter :: tau = 100.d0 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 pdmc(a, nmax, dt, X(irun), accep(irun), tau, E_ref) 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
gfortran hydrogen.F90 qmc_stats.F90 pdmc.F90 -o pdmc ./pdmc
- C
#include <stdio.h> #include <stdlib.h> #include <stddef.h> // for size_t #include <math.h> #include <time.h> #include "hydrogen.h" #include "qmc_stats.h" void pdmc(double a, size_t nmax, double dt, double *energy, double *accep, double tau, double e_ref) { // TODO } int main(void) { #define a 1.2 #define nmax 100000 #define dt 0.05 #define tau 100.0 #define e_ref -0.5 #define nruns 30 double X[nruns]; double Y[nruns]; double ave, err; srand48(time(NULL)); for (size_t irun = 0; irun < nruns; irun++) { pdmc(a, nmax, dt, &X[irun], &Y[irun], tau, e_ref); } ave_error(X, nruns, &ave, &err); printf("E = %f +/- %f\n", ave, err); ave_error(Y, nruns, &ave, &err); printf("A = %f +/- %f\n", ave, err); return 0; }
gcc hydrogen.c qmc_stats.c pdmc.c -lm -o pdmc ./pdmc
5. Project
Change your PDMC code for one of the following:
- the Helium atom, or
- the H2 molecule at
=1.401 bohr.
And compute the ground state energy.
6. Acknowledgments
TREX : Targeting Real Chemical Accuracy at the Exascale project has received funding from the European Union’s Horizon 2020 - Research and Innovation program - under grant agreement no. 952165. The content of this document does not represent the opinion of the European Union, and the European Union is not responsible for any use that might be made of such content.