-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathpi_monte_carlo_openmp.f90
62 lines (52 loc) · 2.2 KB
/
pi_monte_carlo_openmp.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
! Computes an approximation of Pi with a Monte Carlo algorithm
! OpenMP version
! Vincent Magnin, 2021-04-22
! Last modification: 2024-09-03
! MIT license
! $ gfortran -Wall -Wextra -std=f2018 -pedantic -O3 -fopenmp m_xoroshiro128plus.f90 pi_monte_carlo_openmp.f90
! $ ifx -O3 -qopenmp m_xoroshiro128plus.f90 pi_monte_carlo_openmp.f90
program pi_monte_carlo_openmp
use, intrinsic :: iso_fortran_env, only: wp=>real64, int64
use m_xoroshiro128plus
use omp_lib, only: omp_get_thread_num
implicit none
type(rng_t) :: rng ! xoroshiro128+ pseudo-random number generator
real(wp) :: x, y ! Coordinates of a point
integer(int64) :: n ! Total number of points
integer(int64) :: k = 0 ! Points into the quarter disk
integer(int64) :: i ! Loop counter
integer :: t1, t2 ! Clock ticks
real :: count_rate ! Clock ticks per second
integer :: thread ! OpenMP thread number
n = 1000000000
call system_clock(t1, count_rate)
!$OMP PARALLEL DEFAULT(NONE) SHARED(n) PRIVATE(thread, i, x, y, rng) REDUCTION(+: k)
thread = omp_get_thread_num()
! Each image have its own RNG seed, thanks to rng%jump() which
! generates non-overlapping subsequences for parallel computations:
call rng%seed([ -1337_i8, 9812374_i8 ])
! Threads are numbered from 0
if (thread+1 /= 1) then
do i = 2, thread+1
call rng%jump()
end do
end if
x = rng%U01()
!$OMP DO SCHEDULE(STATIC)
do i = 1, n
! Computing a random point (x,y) into the square 0<=x<1, 0<=y<1:
x = rng%U01()
y = rng%U01()
! Is it in the quarter disk (R=1, center=origin) ?
if ((x**2 + y**2) < 1.0_wp) k = k + 1
end do
!$OMP END DO
print '(a, i0, a, i0)', "k", thread, " = ", k
!$OMP END PARALLEL
write(*,*)
write(*, '(a, i0, a, i0)', advance='no') "4 * ", k, " / ", n
write(*, '(a, f17.15)') " = ", (4.0_wp * k) / n
call system_clock(t2)
write(*,'(a, f6.3, a)') "Execution time: ", (t2 - t1) / count_rate, " s"
write(*,'(a)') "---------------------------------------------------"
end program pi_monte_carlo_openmp