| 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 | !*****************************************************************************
|
|---|
| 9 | !
|
|---|
| 10 | ! adi.F90: This file is part of the PolyBench/Fortran 1.0 test suite.
|
|---|
| 11 | !
|
|---|
| 12 | ! Contact: Louis-Noel Pouchet <pouchet@cse.ohio-state.edu>
|
|---|
| 13 | ! Web address: http://polybench.sourceforge.net
|
|---|
| 14 | !
|
|---|
| 15 | !*****************************************************************************
|
|---|
| 16 |
|
|---|
| 17 | !!No data race pairs.
|
|---|
| 18 |
|
|---|
| 19 | ! Include polybench common header.
|
|---|
| 20 | #include "polybench/fpolybench.h"
|
|---|
| 21 |
|
|---|
| 22 | ! Include benchmark-specific header.
|
|---|
| 23 | ! Default data type is double, default size is 10x1024x1024.
|
|---|
| 24 | #include "polybench/adi.h"
|
|---|
| 25 |
|
|---|
| 26 | program adi
|
|---|
| 27 | use omp_lib
|
|---|
| 28 | implicit none
|
|---|
| 29 |
|
|---|
| 30 | POLYBENCH_2D_ARRAY_DECL(x,DATA_TYPE,500, 500)
|
|---|
| 31 | POLYBENCH_2D_ARRAY_DECL(a,DATA_TYPE,500, 500)
|
|---|
| 32 | POLYBENCH_2D_ARRAY_DECL(b,DATA_TYPE,500, 500)
|
|---|
| 33 | polybench_declare_prevent_dce_vars
|
|---|
| 34 | polybench_declare_instruments
|
|---|
| 35 |
|
|---|
| 36 | ! Allocation of Arrays
|
|---|
| 37 | POLYBENCH_ALLOC_2D_ARRAY(x, 500, 500)
|
|---|
| 38 | POLYBENCH_ALLOC_2D_ARRAY(a, 500, 500)
|
|---|
| 39 | POLYBENCH_ALLOC_2D_ARRAY(b, 500, 500)
|
|---|
| 40 |
|
|---|
| 41 | ! Initialization
|
|---|
| 42 | call init_array(500, x, a, b)
|
|---|
| 43 |
|
|---|
| 44 | ! Kernel Execution
|
|---|
| 45 | polybench_start_instruments
|
|---|
| 46 |
|
|---|
| 47 | call kernel_adi(10, 500, x, a, b)
|
|---|
| 48 |
|
|---|
| 49 | polybench_stop_instruments
|
|---|
| 50 | polybench_print_instruments
|
|---|
| 51 |
|
|---|
| 52 | ! Prevent dead-code elimination. All live-out data must be printed
|
|---|
| 53 | ! by the function call in argument.
|
|---|
| 54 | polybench_prevent_dce(print_array(500, x));
|
|---|
| 55 |
|
|---|
| 56 | ! Deallocation of Arrays
|
|---|
| 57 | POLYBENCH_DEALLOC_ARRAY(x)
|
|---|
| 58 | POLYBENCH_DEALLOC_ARRAY(a)
|
|---|
| 59 | POLYBENCH_DEALLOC_ARRAY(b)
|
|---|
| 60 |
|
|---|
| 61 | contains
|
|---|
| 62 |
|
|---|
| 63 | subroutine init_array(n, x, a, b)
|
|---|
| 64 | implicit none
|
|---|
| 65 |
|
|---|
| 66 | DATA_TYPE, dimension(n, n) :: a
|
|---|
| 67 | DATA_TYPE, dimension(n, n) :: x
|
|---|
| 68 | DATA_TYPE, dimension(n, n) :: b
|
|---|
| 69 | integer :: n
|
|---|
| 70 | integer :: i, j, c1, c2
|
|---|
| 71 |
|
|---|
| 72 | if (n >= 1) then
|
|---|
| 73 | !$omp parallel do private(c2)
|
|---|
| 74 | do c1 = 1, n-1
|
|---|
| 75 | do c2 = 1, n-1
|
|---|
| 76 | x(c1,c2) = (DBLE(c1) * DBLE(c2+1) + 1.0D0) / DBLE(n)
|
|---|
| 77 | a(c1,c2) = (DBLE(c1) * DBLE(c2+2) + 2.0D0) / DBLE(n)
|
|---|
| 78 | b(c1,c2) = (DBLE(c1) * DBLE(c2+3) + 3.0D0) / DBLE(n)
|
|---|
| 79 | end do
|
|---|
| 80 | end do
|
|---|
| 81 | !$omp end parallel do
|
|---|
| 82 | end if
|
|---|
| 83 | end subroutine
|
|---|
| 84 |
|
|---|
| 85 | ! DCE code. Must scan the entire live-out data.
|
|---|
| 86 | ! Can be used also to check the correctness of the output.
|
|---|
| 87 |
|
|---|
| 88 | subroutine print_array(n, x)
|
|---|
| 89 | implicit none
|
|---|
| 90 |
|
|---|
| 91 | DATA_TYPE, dimension(n, n) :: x
|
|---|
| 92 | integer :: n
|
|---|
| 93 | integer :: i, j
|
|---|
| 94 |
|
|---|
| 95 | do i = 1, n
|
|---|
| 96 | do j = 1, n
|
|---|
| 97 | write(0, DATA_PRINTF_MODIFIER) x(i, j)
|
|---|
| 98 | if (mod((i * 500 + j), 20) == 0) then
|
|---|
| 99 | write(0, *)
|
|---|
| 100 | end if
|
|---|
| 101 | end do
|
|---|
| 102 | end do
|
|---|
| 103 | write(0, *)
|
|---|
| 104 | end subroutine
|
|---|
| 105 |
|
|---|
| 106 | ! Main computational kernel. The whole function will be timed,
|
|---|
| 107 | ! including the call and return.
|
|---|
| 108 |
|
|---|
| 109 | subroutine kernel_adi(tsteps, n, x, a, b)
|
|---|
| 110 | implicit none
|
|---|
| 111 |
|
|---|
| 112 | DATA_TYPE, dimension(n, n) :: a
|
|---|
| 113 | DATA_TYPE, dimension(n, n) :: x
|
|---|
| 114 | DATA_TYPE, dimension(n, n) :: b
|
|---|
| 115 | integer :: n, tsteps
|
|---|
| 116 | !integer :: i1, i2, t
|
|---|
| 117 | integer :: c0,c2,c8
|
|---|
| 118 |
|
|---|
| 119 | do c0 = 1, 10
|
|---|
| 120 | !$omp parallel do private (c8)
|
|---|
| 121 | do c2 = 1, 500
|
|---|
| 122 | do c8 = 2, 500
|
|---|
| 123 | b(c2,c8) = b(c2,c8)-a(c2,c8)*a(c2,c8)/b(c2,c8-1)
|
|---|
| 124 | end do
|
|---|
| 125 |
|
|---|
| 126 | do c8 = 2, 500
|
|---|
| 127 | x(c2,c8) = x(c2,c8)-x(c2,c8-1)*a(c2,c8)/b(c2,c8-1)
|
|---|
| 128 | end do
|
|---|
| 129 |
|
|---|
| 130 | do c8 = 1, 498
|
|---|
| 131 | x(c2,500-c8) = (X(c2,500-c8) - X(c2,500-c8-1) * A(c2,500 - c8 - 1)) / B(c2,500 - 1 - c8)
|
|---|
| 132 | end do
|
|---|
| 133 | end do
|
|---|
| 134 | !$omp end parallel do
|
|---|
| 135 |
|
|---|
| 136 | !$omp parallel do
|
|---|
| 137 | do c2 = 1, 500
|
|---|
| 138 | x(c2,499) = x(c2,499)/b(c2,499)
|
|---|
| 139 | end do
|
|---|
| 140 | !$omp end parallel do
|
|---|
| 141 |
|
|---|
| 142 | !$omp parallel do private(c8)
|
|---|
| 143 | do c2 = 1, 500
|
|---|
| 144 | do c8 = 2, 500
|
|---|
| 145 | b(c8,c2) = b(c8,c2)-a(c8,c2)*a(c8,c2)/b(c8-1,c2)
|
|---|
| 146 | end do
|
|---|
| 147 |
|
|---|
| 148 | do c8 = 2, 500
|
|---|
| 149 | x(c8,c2) = x(c8,c2)-x(c8-1,c2)*a(c8,c2)/b(c8-1,c2)
|
|---|
| 150 | end do
|
|---|
| 151 |
|
|---|
| 152 | do c8 = 1, 498
|
|---|
| 153 | X(500 - c8,c2) = (X(500 - c8,c2) - X(500 - c8 - 1,c2) * A(500 - 1 - c8,c2)) / B(500 - c8,c2)
|
|---|
| 154 | end do
|
|---|
| 155 | end do
|
|---|
| 156 | !$omp end parallel do
|
|---|
| 157 |
|
|---|
| 158 | !$omp parallel do
|
|---|
| 159 | do c2 = 1, 500
|
|---|
| 160 | x(500-1,c2) = x(499,c2)/b(499,c2)
|
|---|
| 161 | end do
|
|---|
| 162 | !$omp end parallel do
|
|---|
| 163 | end do
|
|---|
| 164 |
|
|---|
| 165 | end subroutine
|
|---|
| 166 |
|
|---|
| 167 | end program
|
|---|