Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
afeiguin
GitHub Repository: afeiguin/comp-phys
Path: blob/master/10_01_montecarlo_integration.ipynb
374 views
Kernel: Python 3

Monte Carlo integration

Imagine that we want to measure the area of a pond with arbitrary shape. Suppose that this pond is in the middle of a field with known area AA. If we throw NN stones randomly, such that they land within the boundaries of the field, and we count the number of stones that fall in the pond NinN_{in}, the area of the pond will be approximately proportional to the fraction of stones that make a splash, multiplied by AA: Apond=NinNA.A_{pond}=\frac{N_{in}}{N}A. This simple procedure is an example of the “Monte Carlo” method.

Simple Monte Carlo integration

More generaly, imagine a rectangle of height HH in the integration interval [a,b][a,b], such that the function f(x)f(x) is within its boundaries. Compute nn pairs of random numbers (xi,yi)(x_i,y_i) such that they are uniformly distributed inside this rectangle. The fraction of points that fall within the area contained below f(x)f(x), i. e., that satisfy yif(xi)y_i \leq f(x_i) is an estimate of the ratio o fthe integral of f(x)f(x) and the area of the rectangle. Hence, the estimate of the integral will be given by: abf(x)dxI(N)=NinNH(ba).\int _a^b{f(x)dx} \simeq I(N) = \frac{N_{in}}{N}H(b-a).

Another Monte Carlo procedure is based on the definition: g=1(ba)abf(x)dx.\langle g \rangle=\frac{1}{(b-a)} \int _a^b{f(x)dx}. In order to determine this average, we sample the value of f(x)f(x): f1Ni=1Nf(xi),\langle f \rangle \simeq \frac{1}{N}\sum_{i=1}^{N}f(x_i), where the NN values xix_i are distributed unformly in the interval [a,b][a,b]. The integral will be given by I(N)=(ba)f.I(N)=(b-a) \langle f \rangle . Monte Carlo error analysis

The Monte Carlo method clearly yields approximate results. The accuracy deppends on the number of values NN that we use for the average. A possible measure of the error is the “variance” σ2\sigma^2 defined by: σ2=f2f2,\sigma ^2=\langle f^2 \rangle - \langle f \rangle ^2, where f=1Ni=1Nf(xi)\langle f \rangle = \frac{1}{N} \sum_{i=1}^N f(x_i) and f2=1Ni=1Nf(xi)2.\langle f^2 \rangle = \frac{1}{N} \sum_{i=1}^{N} f(x_i)^2. The “standard deviation” is σ\sigma. However, we should expect that the error decreases with the number of points NN, and the quantity σ\sigma defines by ([mc_sigma]) does not. Hence, this cannot be a good measure of the error.

Imagine that we perform several measurements of the integral, each of them yielding a result InI_n. These values have been obtained with different sequences of NN random numbers. According to the central limit theorem, these values whould be normally dstributed around a mean I\langle I \rangle. Suppouse that we have a set of MM of such measurements In{I_n}. A convenient measure of the differences of these measurements is the “standard deviation of the means” σM\sigma_M: σM2=I2I2,\sigma_M ^2=\langle I^2 \rangle - \langle I \rangle ^2, where I=1Mn=1MIn\langle I \rangle = \frac{1}{M} \sum_{n=1}^M I_n and I2=1Mn=1MIn2.\langle I^2 \rangle = \frac{1}{M} \sum_{n=1}^{M} I_n^2. It can be proven that σMσ/N.\sigma_M \approx \sigma/\sqrt{N}. This relation becomes exact in the limit of a very large number of measurements. Note that this expression implies that the error decreases with the square root of the number of trials, meaning that if we want to reduce the error by a factor 10, we need 100 times more points for the average.

Exercise 10.1: One dimensional integration

  1. Write a program that implements the “hit and miss” Monte Carlo integration algorithm. Find the estimate I(N)I(N) for the integral of f(x)=41x2f(x)=4\sqrt{1-x^2} as a function of NN, in the interval (0,1)(0,1). Choose H=1H=1, and sample only the xx-dependent part 1x2\sqrt{1-x^2}, and multiply the result by 4. Calculate the difference between I(N)I(N) and the exact result π\pi. This difference is a measure of the error associated with the Monte Carlo estimate. Make a log-log plot of the error as a function of NN. What is the approximate functional deppendece of the error on NN for large NN?

  2. Estimate the integral of f(x)f(x) using the simple Monte Carlo integration by averaging over NN points, using ([mc_integral2]), and compute the error as a function of NN, for NN upt to 10,000. Determine the approximate functional deppendence of the error on NN for large NN. How many trials are necessary to determine INI_N to two decimal places?

  3. Perform 10 measurements In(N)I_n(N), with N=10,000N=10,000 using different random sequences. Show in a table the values of InI_n and σ\sigma according to ([mc_integral2]) and ([mc_sigma]). Use ([mc_sigmam]) to estimate the standard deviation of the means, and compare to the values obtained from ([mc_sigma2]) using the 100,000 values.

  4. To verify that your result for the error is independent of the number of sets you used to divide your data, repeat the previous item grouping your results in 20 groups of 5,000 points each.

Exercise 10.2: Importance of randomness

To examine the effects of a poor random number generator, modify your program to use the linear congruential random number generator using the perameters a=5a=5, c=0c=0 and the seed x1=1x_1=1. Repeat the integral of the previous exercise and compare your results.

Challenge 10.1:

Exercise 10.2

%matplotlib inline import numpy as np from matplotlib import pyplot x = np.arange(0,1,0.02) pyplot.plot(x, 4*np.sqrt(1-x**2))
[<matplotlib.lines.Line2D at 0x119211390>]
Image in a Jupyter notebook
# Hit and miss Monte Carlo integration ngroups = 16 I = np.zeros(ngroups) N = np.zeros(ngroups) E = np.zeros(ngroups) n0 = 100 for i in range(ngroups): N[i] = n0 x = np.random.random(n0) y = np.random.random(n0) I[i] = 0. Nin = 0 for j in range(n0): if(y[j] < np.sqrt(1-x[j]**2)): Nin += 1 I[i] = 4.*float(Nin)/float(n0) E[i] = abs(I[i]-np.pi) print (n0,Nin,I[i],E[i]) n0 *= 2 pyplot.plot(N,E,ls='-',c='red',lw=3); pyplot.plot(N,0.8/np.sqrt(N),ls='-',c='blue',lw=3); pyplot.xscale('log') pyplot.yscale('log')
100 82 3.28 0.1384073464102067 200 163 3.26 0.11840734641020667 400 306 3.06 0.08159265358979306 800 621 3.105 0.036592653589793134 1600 1265 3.1625 0.020907346410206973 3200 2527 3.15875 0.01715734641020683 6400 5046 3.15375 0.012157346410206937 12800 10077 3.1490625 0.00746984641020676 25600 20061 3.13453125 0.007061403589792903 51200 40245 3.144140625 0.002547971410206795 102400 80438 3.142109375 0.0005167214102068662 204800 160702 3.1387109375 0.0028817160897931515 409600 321473 3.139384765625 0.0022078879647930982 819200 644142 3.145224609375 0.0036319557852069195 1638400 1287372 3.142998046875 0.0014053932852067241 3276800 2573360 3.14130859375 0.00028405983979329363
Image in a Jupyter notebook
# Simple Monte Carlo Integration ngroups = 16 I = np.zeros(ngroups) N = np.zeros(ngroups) E = np.zeros(ngroups) n0 = 100 for i in range(ngroups): N[i] = n0 r = np.random.random(n0) I[i] = 0. for j in range(n0): x = r[j] I[i] += np.sqrt(1-x**2) I[i] *= 4./float(n0) E[i] = abs(I[i]-np.pi) print (n0,I[i],E[i]) n0 *= 2 pyplot.plot(N,E,ls='-',c='red',lw=3); pyplot.plot(N,0.8/np.sqrt(N),ls='-',c='blue',lw=3); pyplot.xscale('log') pyplot.yscale('log')
100 3.083548695053129 0.05804395853666433 200 3.0169076292704244 0.12468502431936868 400 3.1806738953569527 0.03908124176715955 800 3.116184732780529 0.025407920809263906 1600 3.174486501377891 0.03289384778809801 3200 3.1495504408897914 0.007957787299998298 6400 3.14574823312641 0.004155579536616827 12800 3.1349325642964856 0.006660089293307525 25600 3.1435910835455867 0.001998429955793579 51200 3.1387551229509194 0.0028375306388737087 102400 3.1403323363592297 0.0012603172305634125 204800 3.139379645209978 0.0022130083798153066 409600 3.1408746244781374 0.0007180291116557491 819200 3.142975298734189 0.0013826451443956778 1638400 3.1414994374643403 9.321612545276636e-05 3276800 3.14129904834528 0.0002936052445132731
Image in a Jupyter notebook
n0 = 100000 I = np.zeros(n0) r = np.random.random(n0) for j in range(n0): x = r[j] I[j] = 4*np.sqrt(1-x**2) def group_measurements(ngroups): global I,n0 nmeasurements = n0//ngroups for n in range(ngroups): Ig = 0. Ig2 = 0. for i in range(n*nmeasurements,(n+1)*nmeasurements): Ig += I[i] Ig2 += I[i]**2 Ig /= nmeasurements Ig2 /= nmeasurements sigma = Ig2-Ig**2 print(Ig,Ig2,sigma) group_measurements(10) print("=============================") group_measurements(20) print("=============================") group_measurements(1)
3.1457636696861324 10.69446669333447 0.7986376278173068 3.1372366968536767 10.639779583215944 0.7975254911305765 3.150995532231518 10.719719006029514 0.7909461618865272 3.1430803624629067 10.687552183237157 0.8085980183372001 3.1549093390648597 10.725146283018782 0.7716933453001129 3.1524565362061563 10.719672155349185 0.7816899426802681 3.128145606617917 10.599685249016789 0.8143903128138117 3.1384054020541434 10.645817378212246 0.7962289105696172 3.1445580189862543 10.683852412854495 0.7956072780837395 3.1332986511743366 10.615776688479952 0.7982162510290358 ============================= 3.136743104879589 10.660760589836224 0.8216032838265797 3.154784234492714 10.728172796832657 0.7755092306288791 3.1387311928691286 10.653576800128173 0.8019432990385091 3.135742200838215 10.625982366303749 0.7931032161860561 3.143268496642032 10.687526668004205 0.8073898260219465 3.158722567820987 10.751911344054774 0.7743830835931629 3.1415836052163484 10.683297105508757 0.8137495569446074 3.1445771197094867 10.691807260965554 0.8034419991651429 3.1845050824958947 10.884185738328162 0.7431131178859758 3.125313595633803 10.566106827709524 0.798521756656033 3.1575195303371184 10.73188976937331 0.761960184912974 3.1473935420752115 10.70745454132511 0.8013684326283634 3.124137415684765 10.567418034780752 0.8071834426992712 3.1321537975510703 10.631952463252777 0.821565051739185 3.1296417982029836 10.57924507556743 0.7845872905082256 3.1471690059053077 10.712389680857035 0.8077169291260322 3.143595399324597 10.681346784020045 0.7991547493652735 3.145520638647886 10.68635804168898 0.7920579535291754 3.1259537277454843 10.587599370989583 0.8160126629836935 3.1406435746031978 10.643954005970329 0.7803119432739773 ============================= 3.1428849815338245 10.673146763274831 0.7954207561239635

Variance reduction

If the function being integrated does not fluctuate too much in the interval of integration, and does not differ much from the average value, then the standard Monte Carlo mean-value method should work well with a reasonable number of points. Otherwise, we will find that the variance is very large, meaning that some points will make small contributions, while others will make large contributions to the integral. If this is the case, the algorithm will be very inefficient. The method can be improved by splitting the function f(x)f(x) in two f(x)=f1(x)+f2(x)f(x)=f_1(x)+f_2(x), such that the integral of f1(x)f_1(x) is known, and f2(x)f_2(x) as a small variance. The “variance reduction” technique, consists then in evaluating the integral of f2(x)f_2(x) to obtain: abf(x)dx=abf1(x)dx+abf2(x)dx=abf1(x)dx+J.\int _a^b{f(x)dx}=\int _a^b {f_1(x)dx} + \int _a^b{f_2(x)dx} = \int _a^b{f_1(x)dx}+J.

Importance Sampling

Imagine that we want to sample the function f(x)=ex2f(x)=e^{-x^2} in the interval [0,1][0,1]. It is evident that most of our points will fall in the region where the value of f(x)f(x) is very small, and therefore we will need a large number of values to achieve a decent accuracy. A way to improve the measurement by reducing the variance is obtained by “importance sampling”. As the name says, the idea is to sample the regions with larger contributions to the integral. For this goal, we introduce a probability distribution P(x)P(x) normalized in the interval of integration abP(x)dx=1.\int _a^b{P(x)dx} = 1. Then, we can rewrite the integral of f(x)f(x) as I=abf(x)P(x)P(x)dxI=\int _a^b{\frac{f(x)}{P(x)}P(x)dx} We can evaluate this integral, by sampling according to the probability distribution P(x)P(x) and evaluating the sum I(N)=1Ni=1Nf(xi)P(xi).I(N)=\frac{1}{N} \sum_{i=1}^N \frac{f(x_i)}{P(x_i)}. Note that for the uniform case P(x)=1/(ba)P(x)=1/(b-a), the expression reduces to the simple Monte Carlo integral.

We are free to choose P(x)P(x) now. We wish to do it in a way to reduce and minimize the variance of the integrand f(x)/P(x)f(x)/P(x). The way to to this is picking a P(x)P(x) that mimics f(x)f(x) where f(x)f(x) is large. if we are able to determine an apropiate P(x)P(x), the integrand will be slowly varying, and hence the variance will be reduced. Another consideration is that the generation of points according to the distribution P(x)P(x) should be a simple task. As an example, let us consider again the integral I=01ex2dx.I=\int _0^1 {e^{-x^2}dx}. A reasonable choice for a weigh function is P(x)=AexP(x)=Ae^{-x}, where AA is a normalization constant.

Notice that for P(x)=f(x)P(x)=f(x) the variance is zero! This is known as the zero variance property. There is a catch, though: The probability function P(x)P(x) needs to be normalized, implying that in reality, P(x)=f(x)/f(x)dxP(x)=f(x)/\int f(x)dx, which assumes that we know in advance precisely the integral that we are trying to calculate!

Exercise 10.3: Importance sampling

  1. Choose the weight function P(x)=exP(x)=e^{-x} and evaluate the integral: 0x3/2exdx.\int _0^{\infty} {x^{3/2}e^{-x}dx}.

  2. Choose P(x)=eaxP(x)=e^{-ax} and estimate the integral 0πdxx2+cos2x.\int _0^{\pi} \frac{dx}{x^2+\cos ^2{x}}. Determine the value of aa that minimizes the variance of the integral.

pyplot.xlim(0,10) x = np.arange(0,10,0.1) pyplot.plot(x,np.exp(-x)); pyplot.plot(x,np.exp(-x**2)); pyplot.plot(x,x**1.5*np.exp(-x));
Image in a Jupyter notebook
# Trapezoidal integration def trapezoids(func, xmin, xmax, nmax): Isim = func(xmin)+func(xmax) h = (xmax-xmin)/nmax for i in range(1,nmax): x = xmin+i*h Isim += 2*func(x) Isim *= h/2 return Isim def f(x): return x**1.5*np.exp(-x) print("Trapezoids: ", trapezoids(f, 0., 20., 100000)) # Simple Monte Carlo integration n0 = 1000000 r = np.random.random(n0) Itot = np.sum(r**1.5*np.exp(-r)) print("Simple Monte Carlo: ", Itot/n0) x = -np.log(r) Itot = np.sum(x**1.5) print("Importance Sampling: ", Itot/n0)
Trapezoids: 1.3293401896452883 Simple Monte Carlo: 0.20055243365045128 Importance Sampling: 1.3288525836135752
pyplot.xlim(0,np.pi) x = np.arange(0,np.pi,0.05) pyplot.plot(x,1./(x**2+np.cos(x)**2)); pyplot.plot(x,np.exp(-x)); pyplot.plot(x,np.exp(-2*x)); pyplot.plot(x,np.exp(-0.2*x));
Image in a Jupyter notebook
# Trapezoidal integration def g(x): return 1./(x**2+np.cos(x)**2) print("Trapezoids: ", trapezoids(g, 0., np.pi, 1000000)) # Simple Monte Carlo integration n0 = 1000000 a = np.arange(0.1,2.1,0.1) I = np.arange(0.1,2.1,0.1) r = np.random.random(n0) I0 = np.sum(1./((r*np.pi)**2+np.cos(r*np.pi)**2)) print("Simple Monte Carlo: ", I0/n0*np.pi) # Importance Sampling print("Importance Sampling:") x = -np.log(r) i = 0 for ai in a: norm = (1.-np.exp(-ai*np.pi))/ai x1 = norm*x/ai Itot = 0. Nin = 0 I2 = 0. for xi in x1: if(xi <= np.pi): Nin += 1 Itot += g(xi)*np.exp(xi*ai) I2 += (g(xi)*np.exp(xi*ai))**2 Itot *= norm I2 *= norm I[i] = Itot/Nin i += 1 print(ai,Itot/Nin,np.sqrt(abs(Itot**2/Nin**2-I2/Nin))/np.sqrt(Nin)) pyplot.plot(a,I,ls='-',marker='o',c='red',lw=3);
Trapezoids: 1.5811879708476726 Simple Monte Carlo: 1.5815947622343893 Importance Sampling: 0.1 1.5244330723519488 0.00321381588763609 0.2 1.50095680220008 0.0020648357111915415 0.30000000000000004 1.4904411006780327 0.0015512787071996397 0.4 1.4968606674298741 0.0012547420671273522 0.5 1.5146208065214364 0.0010525050870957287 0.6 1.5384283496617999 0.0008903917445926586 0.7000000000000001 1.5621715230895488 0.0007325196183547532 0.8 1.5795007309626732 0.0005398433954381976 0.9 1.5848844693695763 0.00018141907194581285 1.0 1.5741096356771225 0.0004972401813311438 1.1 1.5451370063730057 0.0007367892503514896 1.2000000000000002 1.4986266667374715 0.0009147552770853103 1.3000000000000003 1.4371939198631032 0.0010517960421495326 1.4000000000000001 1.364831049309352 0.0011532665188242917 1.5000000000000002 1.2862725322212676 0.0012231156728855844 1.6 1.2052048099895123 0.0012636640779222586 1.7000000000000002 1.1250482091736054 0.0012800101146486459 1.8000000000000003 1.0480214022503385 0.0012761065837385973 1.9000000000000001 0.9756899351567899 0.0012583224714410996 2.0 0.908556238655731 0.0012275630544004439
Image in a Jupyter notebook

Exercise 10.4: The Metropolis algorithm

Use the Metropolis algorithm to sample points according to a ditribution and estimate the integral 04x2exdx,\int _0^4 {x^2e^{-x}dx}, with P(x)=exP(x)=e^{-x} for 0x40 \leq x \leq 4. Plot the number of times the walker is at points x0x_0, x1x_1, x2x_2, ... Is the integrand sampled uniformly? If not, what is the approximate region of xx where the integrand is sampled more often?

delta = 2 xmin = 0. xmax = 4. def f(x): return x**2*np.exp(-x) def P(x): global xmin, xmax if(x < xmin or x > xmax): return 0. return np.exp(-x) def metropolis(xold): global delta xtrial = np.random.random() xtrial = xold+(2*xtrial-1)*delta weight = P(xtrial)/P(xold) xnew = xold if(weight >= 1): #Accept xnew = xtrial elif(weight != 0): r = np.random.random() if(r <= weight): #Accept xnew = xtrial return xnew xwalker = (xmax+xmin)/2. for i in range(100000): xwalker = metropolis(xwalker) I0 = 0. N = 300000 x = np.zeros(N) x[0] = xwalker for i in range(1,N): for j in range(20): xwalker = metropolis(xwalker) x[i] = xwalker I0 += x[i]**2 binwidth=0.1 pyplot.hist(x,bins=np.arange(xmin-1, xmax+1, 0.1),normed=True); print("Trapezoids: ", trapezoids(f,xmin,xmax,100000)) print("Metropolis: ", I0*(1.-np.exp(-4.))/N)
('Trapezoids: ', 1.5237933888733828) ('Metropolis: ', 1.5251458730929552)
Image in a Jupyter notebook
pyplot.hist(x**2,bins=np.arange(xmin**2-1, xmax**2+1, 0.1),normed=True);
Image in a Jupyter notebook

Challenge 10.1

  • Calculate the integral 01x2dx=1/3\int_0^1 x^2 dx=1/3 using simple MC integration and importance sampling with P(x)=xP(x)=x.

  • Calculate the integral 01xdx=2/3\int_0^1 \sqrt{x}dx=2/3 using simple MC integration and P(x)=1eaxP(x)=1-e^{-ax}. Find the values of aa that minimizes the variance.