| 1 | !!!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!!!
|
|---|
| 2 | !!! Copyright (c) 2017-20, Lawrence Livermore National Security, LLC
|
|---|
| 3 | !!! and DataRaceBench project contributors. See the DataRaceBench/COPYRIGHT file for details.
|
|---|
| 4 | !!!
|
|---|
| 5 | !!! SPDX-License-Identifier: (BSD-3-Clause)
|
|---|
| 6 | !!!~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!!!
|
|---|
| 7 |
|
|---|
| 8 | !Two parallel for loops within one single parallel region,
|
|---|
| 9 | !combined with private() and reduction().
|
|---|
| 10 |
|
|---|
| 11 | !3.7969326424804763E-007 vs 3.7969326424804758E-007. There is no race condition. The minute
|
|---|
| 12 | !difference at 22nd point after decimal is due to the precision in fortran95
|
|---|
| 13 |
|
|---|
| 14 | module DRB058
|
|---|
| 15 | use omp_lib
|
|---|
| 16 | implicit none
|
|---|
| 17 |
|
|---|
| 18 | integer :: MSIZE
|
|---|
| 19 | integer :: n,m,mits
|
|---|
| 20 | integer, parameter :: dp = kind(1.0d0)
|
|---|
| 21 | real(dp), dimension(:,:), pointer :: u,f,uold
|
|---|
| 22 | real(dp) :: dx,dy,tol,relax,alpha
|
|---|
| 23 |
|
|---|
| 24 | contains
|
|---|
| 25 | subroutine initialize()
|
|---|
| 26 | integer :: i,j,xx,yy
|
|---|
| 27 |
|
|---|
| 28 | MSIZE = 200
|
|---|
| 29 | mits = 1000
|
|---|
| 30 | tol=0.0000000001
|
|---|
| 31 | relax = 1.0
|
|---|
| 32 | alpha = 0.0543
|
|---|
| 33 | n = MSIZE
|
|---|
| 34 | m = MSIZE
|
|---|
| 35 | allocate(u(MSIZE,MSIZE))
|
|---|
| 36 | allocate(f(MSIZE,MSIZE))
|
|---|
| 37 | allocate(uold(MSIZE,MSIZE))
|
|---|
| 38 |
|
|---|
| 39 | dx = 2.0D0 / DBLE(n-1)
|
|---|
| 40 | dy = 2.0D0 / DBLE(m-1)
|
|---|
| 41 |
|
|---|
| 42 | do i = 1, n
|
|---|
| 43 | do j = 1, m
|
|---|
| 44 | xx = int(-1.0 + dx * (i-1))
|
|---|
| 45 | yy = int(-1.0 + dy * (i-1))
|
|---|
| 46 | u(i,j) = 0.0
|
|---|
| 47 | f(i,j) = -1.0 * alpha * (1.0-xx*xx) * (1.0-yy*yy) - 2.0* (1.0-xx*xx) -2.0 * (1.0-yy*yy)
|
|---|
| 48 | end do
|
|---|
| 49 | end do
|
|---|
| 50 |
|
|---|
| 51 | end subroutine
|
|---|
| 52 |
|
|---|
| 53 | subroutine jacobi()
|
|---|
| 54 | integer, parameter :: dp = kind(1.0d0)
|
|---|
| 55 | real(dp) :: omega
|
|---|
| 56 | integer :: i, j, k
|
|---|
| 57 | real(dp) :: error, resid, ax, ay, b
|
|---|
| 58 |
|
|---|
| 59 | MSIZE = 200
|
|---|
| 60 | mits = 1000
|
|---|
| 61 | tol=0.0000000001
|
|---|
| 62 | relax = 1.0
|
|---|
| 63 | alpha = 0.0543
|
|---|
| 64 | n = MSIZE
|
|---|
| 65 | m = MSIZE
|
|---|
| 66 |
|
|---|
| 67 | omega = relax
|
|---|
| 68 | dx = 2.0D0 / DBLE(n-1)
|
|---|
| 69 | dy = 2.0D0 / DBLE(m-1)
|
|---|
| 70 |
|
|---|
| 71 | ax = 1.0D0 / (dx * dx); !/* X-direction coef */
|
|---|
| 72 | ay = 1.0D0 / (dy * dy); !/* Y-direction coef */
|
|---|
| 73 | b = -2.0D0 / (dx * dx) - 2.0D0 / (dy * dy) - alpha;
|
|---|
| 74 |
|
|---|
| 75 | error = 10.0 * tol
|
|---|
| 76 | k = 1
|
|---|
| 77 |
|
|---|
| 78 | do k = 1, mits
|
|---|
| 79 | error = 0.0
|
|---|
| 80 |
|
|---|
| 81 | !Copy new solution into old
|
|---|
| 82 | !$omp parallel
|
|---|
| 83 | !$omp do private(i,j)
|
|---|
| 84 | do i = 1, n
|
|---|
| 85 | do j = 1, m
|
|---|
| 86 | uold(i,j) = u(i,j)
|
|---|
| 87 | end do
|
|---|
| 88 | end do
|
|---|
| 89 | !$omp end do
|
|---|
| 90 | !$omp do private(i,j,resid) reduction(+:error)
|
|---|
| 91 | do i = 2, (n-1)
|
|---|
| 92 | do j = 2, (m-1)
|
|---|
| 93 | resid = (ax * (uold(i - 1,j) + uold(i + 1,j)) + ay * (uold(i,j - 1) + uold(i,j + 1)) + b * uold(i,j) - f(i,j)) / b
|
|---|
| 94 | u(i,j) = uold(i,j) - omega * resid
|
|---|
| 95 | error = error + resid * resid
|
|---|
| 96 | end do
|
|---|
| 97 | end do
|
|---|
| 98 | !$omp end do nowait
|
|---|
| 99 | !$omp end parallel
|
|---|
| 100 |
|
|---|
| 101 | !Error check
|
|---|
| 102 | error = sqrt(error)/(n*m)
|
|---|
| 103 | end do
|
|---|
| 104 |
|
|---|
| 105 | print*,"Total number of iterations: ",k
|
|---|
| 106 | print*,"Residual: ",error
|
|---|
| 107 |
|
|---|
| 108 | end subroutine
|
|---|
| 109 | end module
|
|---|
| 110 |
|
|---|
| 111 | program DRB058_jacobikernel_orig_no
|
|---|
| 112 | use omp_lib
|
|---|
| 113 | use DRB058
|
|---|
| 114 | implicit none
|
|---|
| 115 |
|
|---|
| 116 | call initialize()
|
|---|
| 117 | call jacobi()
|
|---|
| 118 | end program
|
|---|