nf_mm.f90 7.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237
  1. ! StarPU --- Runtime system for heterogeneous multicore architectures.
  2. !
  3. ! Copyright (C) 2016-2020 Université de Bordeaux, CNRS (LaBRI UMR 5800), Inria
  4. !
  5. ! StarPU is free software; you can redistribute it and/or modify
  6. ! it under the terms of the GNU Lesser General Public License as published by
  7. ! the Free Software Foundation; either version 2.1 of the License, or (at
  8. ! your option) any later version.
  9. !
  10. ! StarPU is distributed in the hope that it will be useful, but
  11. ! WITHOUT ANY WARRANTY; without even the implied warranty of
  12. ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  13. !
  14. ! See the GNU Lesser General Public License in COPYING.LGPL for more details.
  15. !
  16. program nf_mm
  17. use iso_c_binding ! C interfacing module
  18. use fstarpu_mod ! StarPU interfacing module
  19. use fstarpu_mpi_mod ! StarPU-MPI interfacing module
  20. use nf_mm_cl
  21. implicit none
  22. logical, parameter :: verbose = .false.
  23. integer(c_int) :: comm_rank, comm_size, comm_world
  24. integer(c_int) :: N = 16, BS = 4, NB
  25. real(kind=c_double),allocatable,target :: A(:,:), B(:,:), C(:,:)
  26. type(c_ptr),allocatable :: dh_A(:), dh_B(:), dh_C(:,:)
  27. type(c_ptr) :: cl_mm
  28. integer(c_int) :: ncpu
  29. integer(c_int) :: ret
  30. integer(c_int) :: row, col
  31. integer(c_int) :: b_row, b_col
  32. integer(c_int) :: mr, tag, rank
  33. ret = fstarpu_init(C_NULL_PTR)
  34. if (ret == -19) then
  35. stop 77
  36. else if (ret /= 0) then
  37. stop 1
  38. end if
  39. ret = fstarpu_mpi_init(1)
  40. print *,"fstarpu_mpi_init status:", ret
  41. if (ret /= 0) then
  42. stop 1
  43. end if
  44. ! stop there if no CPU worker available
  45. ncpu = fstarpu_cpu_worker_get_count()
  46. if (ncpu == 0) then
  47. call fstarpu_shutdown()
  48. stop 77
  49. end if
  50. comm_world = fstarpu_mpi_world_comm()
  51. comm_size = fstarpu_mpi_world_size()
  52. comm_rank = fstarpu_mpi_world_rank()
  53. if (comm_size < 2) then
  54. call fstarpu_shutdown()
  55. ret = fstarpu_mpi_shutdown()
  56. stop 77
  57. end if
  58. ! TODO: process app's argc/argv
  59. NB = N/BS
  60. ! allocate and initialize codelet
  61. cl_mm = fstarpu_codelet_allocate()
  62. call fstarpu_codelet_set_name(cl_mm, c_char_"nf_mm_cl"//c_null_char)
  63. call fstarpu_codelet_add_cpu_func(cl_mm, C_FUNLOC(cl_cpu_mult))
  64. call fstarpu_codelet_add_buffer(cl_mm, FSTARPU_R)
  65. call fstarpu_codelet_add_buffer(cl_mm, FSTARPU_R)
  66. call fstarpu_codelet_add_buffer(cl_mm, FSTARPU_RW)
  67. ! allocate matrices
  68. if (comm_rank == 0) then
  69. allocate(A(N,N))
  70. allocate(B(N,N))
  71. allocate(C(N,N))
  72. end if
  73. ! init matrices
  74. if (comm_rank == 0) then
  75. do col=1,N
  76. do row=1,N
  77. if (row == col) then
  78. A(row,col) = 2
  79. else
  80. A(row,col) = 0
  81. end if
  82. B(row,col) = row*N+col
  83. C(row,col) = 0
  84. end do
  85. end do
  86. if (verbose) then
  87. print *,"A"
  88. call mat_disp(A)
  89. print *,"B"
  90. call mat_disp(B)
  91. print *,"C"
  92. call mat_disp(C)
  93. end if
  94. end if
  95. ! allocate data handles
  96. allocate(dh_A(NB))
  97. allocate(dh_B(NB))
  98. allocate(dh_C(NB,NB))
  99. ! register matrices
  100. if (comm_rank == 0) then
  101. mr = 0 ! TODO: use STARPU_MAIN_RAM constant
  102. else
  103. mr = -1
  104. end if
  105. tag = 0
  106. do b_row=1,NB
  107. if (comm_rank == 0) then
  108. call fstarpu_matrix_data_register(dh_A(b_row), mr, &
  109. c_loc( A(1+(b_row-1)*BS,1) ), N, BS, N, c_sizeof(A(1,1)))
  110. else
  111. call fstarpu_matrix_data_register(dh_A(b_row), mr, &
  112. c_null_ptr, N, BS, N, c_sizeof(A(1,1)))
  113. end if
  114. call fstarpu_mpi_data_register(dh_A(b_row), tag, 0)
  115. tag = tag+1
  116. end do
  117. do b_col=1,NB
  118. if (comm_rank == 0) then
  119. call fstarpu_matrix_data_register(dh_B(b_col), mr, &
  120. c_loc( B(1,1+(b_col-1)*BS) ), N, N, BS, c_sizeof(B(1,1)))
  121. else
  122. call fstarpu_matrix_data_register(dh_B(b_col), mr, &
  123. c_null_ptr, N, N, BS, c_sizeof(B(1,1)))
  124. end if
  125. call fstarpu_mpi_data_register(dh_B(b_col), tag, 0)
  126. tag = tag+1
  127. end do
  128. do b_col=1,NB
  129. do b_row=1,NB
  130. if (comm_rank == 0) then
  131. call fstarpu_matrix_data_register(dh_C(b_row,b_col), mr, &
  132. c_loc( C(1+(b_row-1)*BS,1+(b_col-1)*BS) ), N, BS, BS, c_sizeof(C(1,1)))
  133. else
  134. call fstarpu_matrix_data_register(dh_C(b_row,b_col), mr, &
  135. c_null_ptr, N, BS, BS, c_sizeof(C(1,1)))
  136. end if
  137. call fstarpu_mpi_data_register(dh_C(b_row,b_col), tag, 0)
  138. tag = tag+1
  139. end do
  140. end do
  141. ! distribute matrix C
  142. do b_col=1,NB
  143. do b_row=1,NB
  144. rank = modulo(b_row+b_col, comm_size)
  145. call fstarpu_mpi_data_migrate(comm_world, dh_c(b_row,b_col), rank)
  146. end do
  147. end do
  148. do b_col=1,NB
  149. do b_row=1,NB
  150. ret = fstarpu_mpi_task_insert(comm_world, (/ cl_mm, &
  151. FSTARPU_R, dh_A(b_row), &
  152. FSTARPU_R, dh_B(b_col), &
  153. FSTARPU_RW, dh_C(b_row,b_col), &
  154. C_NULL_PTR /))
  155. end do
  156. end do
  157. call fstarpu_task_wait_for_all()
  158. ! undistribute matrix C
  159. do b_col=1,NB
  160. do b_row=1,NB
  161. call fstarpu_mpi_data_migrate(comm_world, dh_c(b_row,b_col), 0)
  162. end do
  163. end do
  164. ! unregister matrices
  165. do b_row=1,NB
  166. call fstarpu_data_unregister(dh_A(b_row))
  167. end do
  168. do b_col=1,NB
  169. call fstarpu_data_unregister(dh_B(b_col))
  170. end do
  171. do b_col=1,NB
  172. do b_row=1,NB
  173. call fstarpu_data_unregister(dh_C(b_row,b_col))
  174. end do
  175. end do
  176. ! check result
  177. if (comm_rank == 0) then
  178. if (verbose) then
  179. print *,"final C"
  180. call mat_disp(C)
  181. end if
  182. do col=1,N
  183. do row=1,N
  184. if (abs(C(row,col) - 2*(row*N+col)) > 1.0) then
  185. print *, "check failed"
  186. stop 1
  187. end if
  188. end do
  189. end do
  190. end if
  191. ! free handles
  192. deallocate(dh_A)
  193. deallocate(dh_B)
  194. deallocate(dh_C)
  195. ! free matrices
  196. if (comm_rank == 0) then
  197. deallocate(A)
  198. deallocate(B)
  199. deallocate(C)
  200. end if
  201. call fstarpu_codelet_free(cl_mm)
  202. call fstarpu_shutdown()
  203. ret = fstarpu_mpi_shutdown()
  204. print *,"fstarpu_mpi_shutdown status:", ret
  205. if (ret /= 0) then
  206. stop 1
  207. end if
  208. end program nf_mm