|
| 1 | +!*********************************************************************** |
| 2 | +!* Apache License 2.0 |
| 3 | +!* |
| 4 | +!* This file is part of the GFDL Flexible Modeling System (FMS). |
| 5 | +!* |
| 6 | +!* Licensed under the Apache License, Version 2.0 (the "License"); |
| 7 | +!* you may not use this file except in compliance with the License. |
| 8 | +!* You may obtain a copy of the License at |
| 9 | +!* |
| 10 | +!* http://www.apache.org/licenses/LICENSE-2.0 |
| 11 | +!* |
| 12 | +!* FMS is distributed in the hope that it will be useful, but WITHOUT |
| 13 | +!* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied; |
| 14 | +!* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A |
| 15 | +!* PARTICULAR PURPOSE. See the License for the specific language |
| 16 | +!* governing permissions and limitations under the License. |
| 17 | +!*********************************************************************** |
| 18 | +#ifdef SYSTEM_CLOCK |
| 19 | +#undef SYSTEM_CLOCK |
| 20 | +#endif |
| 21 | + |
| 22 | +!> @author Andrew Brooks |
| 23 | +!> @brief Test mpp_pelist_gather and mpp_pelist_scatter routines for generalized indices. |
| 24 | +program test_mpp_pelist_gatscat_gen_ind |
| 25 | + use mpp_mod |
| 26 | + use platform_mod |
| 27 | + |
| 28 | + implicit none |
| 29 | + |
| 30 | + integer :: pe, npes, root |
| 31 | + integer :: dim_order(3) |
| 32 | + integer :: perms(3,6) |
| 33 | + integer :: p |
| 34 | + |
| 35 | + call mpp_init(mpp_init_test_requests_allocated) |
| 36 | + call mpp_set_stack_size(3145746) |
| 37 | + |
| 38 | + pe = mpp_pe() |
| 39 | + npes = mpp_npes() |
| 40 | + root = mpp_root_pe() |
| 41 | + |
| 42 | + perms = reshape([ & |
| 43 | + 1, 2, 3, & |
| 44 | + 1, 3, 2, & |
| 45 | + 2, 1, 3, & |
| 46 | + 2, 3, 1, & |
| 47 | + 3, 1, 2, & |
| 48 | + 3, 2, 1 ], [3, 6]) |
| 49 | + |
| 50 | + dim_order = (/2, 3, 1/) |
| 51 | + |
| 52 | + if (pe == root) print *, '--- PELIST SCATTER/GATHER TESTS ---' |
| 53 | + |
| 54 | + do p = 1, 6 |
| 55 | + |
| 56 | + dim_order = perms(:,p) |
| 57 | + |
| 58 | + if (pe == root) then |
| 59 | + print *, '--------------------------------' |
| 60 | + print *, 'dim_order =', dim_order |
| 61 | + endif |
| 62 | + |
| 63 | + call test_scatter(npes, pe, root, dim_order) |
| 64 | + call mpp_sync() |
| 65 | + |
| 66 | + call test_gather(npes, pe, root, dim_order) |
| 67 | + call mpp_sync() |
| 68 | + |
| 69 | + enddo |
| 70 | + |
| 71 | + if (pe == root) print *, 'ALL PERMUTATIONS PASSED' |
| 72 | + |
| 73 | + call mpp_exit() |
| 74 | + |
| 75 | +contains |
| 76 | + |
| 77 | +!> @brief Map logical indices (i,j,k) to a unique scalar test value. |
| 78 | + pure real function val(i, j, k) |
| 79 | + integer, intent(in) :: i, j, k |
| 80 | + |
| 81 | + val = 100.0*i + 10.0*j + k |
| 82 | + end function val |
| 83 | + |
| 84 | +!> @brief Convert logical indices (i,j,k) into storage indices using dim_order. |
| 85 | + subroutine permute(i, j, k, dim_order, u, v, w) |
| 86 | + integer, intent(in) :: i, j, k, dim_order(3) |
| 87 | + integer, intent(out) :: u, v, w |
| 88 | + |
| 89 | + integer :: idx(3) |
| 90 | + |
| 91 | + idx = (/i, j, k/) |
| 92 | + u = idx(dim_order(1)) |
| 93 | + v = idx(dim_order(2)) |
| 94 | + w = idx(dim_order(3)) |
| 95 | + end subroutine permute |
| 96 | + |
| 97 | +!> @brief Compute 1D domain decomposition in i-direction for each MPI rank. |
| 98 | + subroutine get_decomp(pe, npes, NI, NJ, is, ie, js, je) |
| 99 | + integer, intent(in) :: pe, npes, NI, NJ |
| 100 | + integer, intent(out) :: is, ie, js, je |
| 101 | + |
| 102 | + integer :: chunk |
| 103 | + |
| 104 | + chunk = NI / npes |
| 105 | + is = pe*chunk + 1 |
| 106 | + ie = (pe+1)*chunk |
| 107 | + js = 1 |
| 108 | + je = NJ |
| 109 | + end subroutine get_decomp |
| 110 | + |
| 111 | +!> @brief Construct pelist containing all MPI ranks participating in the operation. |
| 112 | + subroutine build_pelist(npes, pelist) |
| 113 | + integer, intent(in) :: npes |
| 114 | + integer, intent(out) :: pelist(npes) |
| 115 | + |
| 116 | + integer :: i |
| 117 | + |
| 118 | + do i=0,npes-1 |
| 119 | + pelist(i+1)=i |
| 120 | + enddo |
| 121 | + end subroutine build_pelist |
| 122 | + |
| 123 | +!> @brief Allocate a 3D field in permuted layout, handling root and non-root cases. |
| 124 | + subroutine alloc_field(field, dim_order, dims_logical, pe, root, is_global) |
| 125 | + real, allocatable, intent(out) :: field(:,:,:) |
| 126 | + integer, intent(in) :: dim_order(3) |
| 127 | + integer, intent(in) :: dims_logical(3) |
| 128 | + integer, intent(in) :: pe, root |
| 129 | + logical, intent(in) :: is_global |
| 130 | + |
| 131 | + integer :: dims(3) |
| 132 | + |
| 133 | + if (is_global) then |
| 134 | + if (pe == root) then |
| 135 | + dims = dims_logical |
| 136 | + allocate(field(dims(dim_order(1)), & |
| 137 | + dims(dim_order(2)), & |
| 138 | + dims(dim_order(3)))) |
| 139 | + else |
| 140 | + allocate(field(1,1,1)) |
| 141 | + endif |
| 142 | + else |
| 143 | + dims = dims_logical |
| 144 | + allocate(field(dims(dim_order(1)), & |
| 145 | + dims(dim_order(2)), & |
| 146 | + dims(dim_order(3)))) |
| 147 | + endif |
| 148 | + end subroutine alloc_field |
| 149 | + |
| 150 | +!> @brief Populate a field with values from val(i,j,k) under the given permutation. |
| 151 | + subroutine fill_from_val(segment, is, ie, NJ, NK, dim_order) |
| 152 | + real, intent(inout) :: segment(:,:,:) |
| 153 | + integer, intent(in) :: is, ie, NJ, NK, dim_order(3) |
| 154 | + |
| 155 | + integer :: i, j, k |
| 156 | + integer :: u, v, w |
| 157 | + |
| 158 | + do i=is,ie |
| 159 | + do j=1,NJ |
| 160 | + do k=1,NK |
| 161 | + call permute(i-is+1, j, k, dim_order, u, v, w) |
| 162 | + segment(u,v,w) = val(i,j,k) |
| 163 | + enddo |
| 164 | + enddo |
| 165 | + enddo |
| 166 | + end subroutine fill_from_val |
| 167 | + |
| 168 | +!> @brief Verify field values against val(i,j,k) for local (scatter) or global (gather) domains. |
| 169 | + subroutine check_answer(field, is, ie, NI, NJ, NK, dim_order, check_global) |
| 170 | + real, intent(in) :: field(:,:,:) |
| 171 | + integer, intent(in) :: is, ie, NI, NJ, NK |
| 172 | + integer, intent(in) :: dim_order(3) |
| 173 | + logical, intent(in) :: check_global ! .true. = gather, .false. = scatter |
| 174 | + |
| 175 | + integer :: i, j, k |
| 176 | + integer :: u, v, w |
| 177 | + integer :: istart, iend, iloc |
| 178 | + |
| 179 | + ! --- choose iteration domain --- |
| 180 | + if (check_global) then |
| 181 | + istart = 1 |
| 182 | + iend = NI |
| 183 | + else |
| 184 | + istart = is |
| 185 | + iend = ie |
| 186 | + endif |
| 187 | + |
| 188 | + do i = istart, iend |
| 189 | + do j = 1, NJ |
| 190 | + do k = 1, NK |
| 191 | + |
| 192 | + if (check_global) then |
| 193 | + ! gather -> global indexing |
| 194 | + call permute(i, j, k, dim_order, u, v, w) |
| 195 | + else |
| 196 | + ! scatter -> local indexing |
| 197 | + iloc = i - is + 1 |
| 198 | + call permute(iloc, j, k, dim_order, u, v, w) |
| 199 | + endif |
| 200 | + |
| 201 | + if (field(u,v,w) /= val(i,j,k)) then |
| 202 | + print *, 'FAIL' |
| 203 | + print *, 'i,j,k=', i,j,k |
| 204 | + print *, 'u,v,w=', u,v,w |
| 205 | + print *, 'got=', field(u,v,w) |
| 206 | + print *, 'expected=', val(i,j,k) |
| 207 | + call mpp_error(FATAL,'check_answer failed') |
| 208 | + endif |
| 209 | + |
| 210 | + enddo |
| 211 | + enddo |
| 212 | + enddo |
| 213 | + end subroutine check_answer |
| 214 | + |
| 215 | +!> @brief Test mpp_scatter with pelist and permuted storage layouts. |
| 216 | + subroutine test_scatter(npes, pe, root, dim_order) |
| 217 | + integer,intent(in) :: npes, pe, root |
| 218 | + integer,intent(in) :: dim_order(3) |
| 219 | + |
| 220 | + integer :: pelist(npes) |
| 221 | + integer :: is, ie, js, je |
| 222 | + integer :: NI, NJ, NK |
| 223 | + real, allocatable :: global_perm(:,:,:) |
| 224 | + real, allocatable :: segment(:,:,:) |
| 225 | + |
| 226 | + NI = 8 |
| 227 | + NJ = 5 |
| 228 | + NK = 3 |
| 229 | + |
| 230 | + call get_decomp(pe, npes, NI, NJ, is, ie, js, je) |
| 231 | + call build_pelist(npes, pelist) |
| 232 | + |
| 233 | + call alloc_field(global_perm, dim_order, (/NI, NJ, NK/), pe, root, .true.) |
| 234 | + call alloc_field(segment, dim_order, (/ie-is+1, NJ, NK/), pe, root, .false.) |
| 235 | + |
| 236 | + segment = -2.0 |
| 237 | + |
| 238 | + if (pe == root) then |
| 239 | + call fill_from_val(global_perm, 1, NI, NJ, NK, dim_order) |
| 240 | + endif |
| 241 | + |
| 242 | + ! --- scatter --- |
| 243 | + if (pe == root) then |
| 244 | + call mpp_scatter(is, ie, js, je, NK, pelist, segment, global_perm, dim_order, .true.) |
| 245 | + else |
| 246 | + call mpp_scatter(is, ie, js, je, NK, pelist, segment, global_perm, dim_order, .false.) |
| 247 | + endif |
| 248 | + |
| 249 | + call check_answer(segment, is, ie, NI, NJ, NK, dim_order, .false.) |
| 250 | + call mpp_sync() |
| 251 | + |
| 252 | + if (pe == root) print *, 'SCATTER PASS' |
| 253 | + |
| 254 | + deallocate(segment) |
| 255 | + if (pe == root) deallocate(global_perm) |
| 256 | + end subroutine test_scatter |
| 257 | + |
| 258 | +!> @brief Test mpp_gather with pelist and permuted storage layouts. |
| 259 | + subroutine test_gather(npes, pe, root, dim_order) |
| 260 | + integer,intent(in) :: npes, pe, root |
| 261 | + integer,intent(in) :: dim_order(3) |
| 262 | + |
| 263 | + integer :: pelist(npes) |
| 264 | + integer :: is, ie, js, je |
| 265 | + integer :: NI,NJ,NK |
| 266 | + real, allocatable :: segment(:,:,:) |
| 267 | + real, allocatable :: gather_data(:,:,:) |
| 268 | + |
| 269 | + NI = 8 |
| 270 | + NJ = 5 |
| 271 | + NK = 3 |
| 272 | + |
| 273 | + call get_decomp(pe, npes, NI, NJ, is, ie, js, je) |
| 274 | + call build_pelist(npes, pelist) |
| 275 | + |
| 276 | + call alloc_field(gather_data, dim_order, (/NI, NJ, NK/), pe, root, .true.) |
| 277 | + call alloc_field(segment, dim_order, (/ie-is+1, NJ, NK/), pe, root, .false.) |
| 278 | + |
| 279 | + call fill_from_val(segment, is, ie, NJ, NK, dim_order) |
| 280 | + |
| 281 | + ! --- GATHER --- |
| 282 | + if (pe == root) then |
| 283 | + call mpp_gather(is, ie, js, je, NK, pelist, segment, gather_data, dim_order, .true.) |
| 284 | + else |
| 285 | + call mpp_gather(is, ie, js, je, NK, pelist, segment, gather_data, dim_order, .false.) |
| 286 | + endif |
| 287 | + |
| 288 | + call mpp_sync() |
| 289 | + |
| 290 | + if (pe == root) then |
| 291 | + call check_answer(gather_data, is, ie, NI, NJ, NK, dim_order, .true.) |
| 292 | + print *, 'GATHER PASS' |
| 293 | + endif |
| 294 | + |
| 295 | + deallocate(segment) |
| 296 | + if (pe == root) deallocate(gather_data) |
| 297 | + |
| 298 | + end subroutine test_gather |
| 299 | + |
| 300 | +end program test_mpp_pelist_gatscat_gen_ind |
0 commit comments