nf_mm.f90 7.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238
  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_size, comm_rank
  24. integer(c_int), target :: comm_world
  25. integer(c_int) :: N = 16, BS = 4, NB
  26. real(kind=c_double),allocatable,target :: A(:,:), B(:,:), C(:,:)
  27. type(c_ptr),allocatable :: dh_A(:), dh_B(:), dh_C(:,:)
  28. type(c_ptr) :: cl_mm
  29. integer(c_int) :: ncpu
  30. integer(c_int) :: ret
  31. integer(c_int) :: row, col
  32. integer(c_int) :: b_row, b_col
  33. integer(c_int) :: mr, tag, rank
  34. ret = fstarpu_init(C_NULL_PTR)
  35. if (ret == -19) then
  36. stop 77
  37. else if (ret /= 0) then
  38. stop 1
  39. end if
  40. ret = fstarpu_mpi_init(1)
  41. print *,"fstarpu_mpi_init status:", ret
  42. if (ret /= 0) then
  43. stop 1
  44. end if
  45. ! stop there if no CPU worker available
  46. ncpu = fstarpu_cpu_worker_get_count()
  47. if (ncpu == 0) then
  48. call fstarpu_shutdown()
  49. stop 77
  50. end if
  51. comm_world = fstarpu_mpi_world_comm()
  52. comm_size = fstarpu_mpi_world_size()
  53. comm_rank = fstarpu_mpi_world_rank()
  54. if (comm_size < 2) then
  55. call fstarpu_shutdown()
  56. ret = fstarpu_mpi_shutdown()
  57. stop 77
  58. end if
  59. ! TODO: process app's argc/argv
  60. NB = N/BS
  61. ! allocate and initialize codelet
  62. cl_mm = fstarpu_codelet_allocate()
  63. call fstarpu_codelet_set_name(cl_mm, c_char_"nf_mm_cl"//c_null_char)
  64. call fstarpu_codelet_add_cpu_func(cl_mm, C_FUNLOC(cl_cpu_mult))
  65. call fstarpu_codelet_add_buffer(cl_mm, FSTARPU_R)
  66. call fstarpu_codelet_add_buffer(cl_mm, FSTARPU_R)
  67. call fstarpu_codelet_add_buffer(cl_mm, FSTARPU_RW)
  68. ! allocate matrices
  69. if (comm_rank == 0) then
  70. allocate(A(N,N))
  71. allocate(B(N,N))
  72. allocate(C(N,N))
  73. end if
  74. ! init matrices
  75. if (comm_rank == 0) then
  76. do col=1,N
  77. do row=1,N
  78. if (row == col) then
  79. A(row,col) = 2
  80. else
  81. A(row,col) = 0
  82. end if
  83. B(row,col) = row*N+col
  84. C(row,col) = 0
  85. end do
  86. end do
  87. if (verbose) then
  88. print *,"A"
  89. call mat_disp(A)
  90. print *,"B"
  91. call mat_disp(B)
  92. print *,"C"
  93. call mat_disp(C)
  94. end if
  95. end if
  96. ! allocate data handles
  97. allocate(dh_A(NB))
  98. allocate(dh_B(NB))
  99. allocate(dh_C(NB,NB))
  100. ! register matrices
  101. if (comm_rank == 0) then
  102. mr = 0 ! TODO: use STARPU_MAIN_RAM constant
  103. else
  104. mr = -1
  105. end if
  106. tag = 0
  107. do b_row=1,NB
  108. if (comm_rank == 0) then
  109. call fstarpu_matrix_data_register(dh_A(b_row), mr, &
  110. c_loc( A(1+(b_row-1)*BS,1) ), N, BS, N, c_sizeof(A(1,1)))
  111. else
  112. call fstarpu_matrix_data_register(dh_A(b_row), mr, &
  113. c_null_ptr, N, BS, N, c_sizeof(A(1,1)))
  114. end if
  115. call fstarpu_mpi_data_register(dh_A(b_row), tag, 0)
  116. tag = tag+1
  117. end do
  118. do b_col=1,NB
  119. if (comm_rank == 0) then
  120. call fstarpu_matrix_data_register(dh_B(b_col), mr, &
  121. c_loc( B(1,1+(b_col-1)*BS) ), N, N, BS, c_sizeof(B(1,1)))
  122. else
  123. call fstarpu_matrix_data_register(dh_B(b_col), mr, &
  124. c_null_ptr, N, N, BS, c_sizeof(B(1,1)))
  125. end if
  126. call fstarpu_mpi_data_register(dh_B(b_col), tag, 0)
  127. tag = tag+1
  128. end do
  129. do b_col=1,NB
  130. do b_row=1,NB
  131. if (comm_rank == 0) then
  132. call fstarpu_matrix_data_register(dh_C(b_row,b_col), mr, &
  133. c_loc( C(1+(b_row-1)*BS,1+(b_col-1)*BS) ), N, BS, BS, c_sizeof(C(1,1)))
  134. else
  135. call fstarpu_matrix_data_register(dh_C(b_row,b_col), mr, &
  136. c_null_ptr, N, BS, BS, c_sizeof(C(1,1)))
  137. end if
  138. call fstarpu_mpi_data_register(dh_C(b_row,b_col), tag, 0)
  139. tag = tag+1
  140. end do
  141. end do
  142. ! distribute matrix C
  143. do b_col=1,NB
  144. do b_row=1,NB
  145. rank = modulo(b_row+b_col, comm_size)
  146. call fstarpu_mpi_data_migrate(comm_world, dh_c(b_row,b_col), rank)
  147. end do
  148. end do
  149. do b_col=1,NB
  150. do b_row=1,NB
  151. call fstarpu_mpi_task_insert((/ c_loc(comm_world), cl_mm, &
  152. FSTARPU_R, dh_A(b_row), &
  153. FSTARPU_R, dh_B(b_col), &
  154. FSTARPU_RW, dh_C(b_row,b_col), &
  155. C_NULL_PTR /))
  156. end do
  157. end do
  158. call fstarpu_task_wait_for_all()
  159. ! undistribute matrix C
  160. do b_col=1,NB
  161. do b_row=1,NB
  162. call fstarpu_mpi_data_migrate(comm_world, dh_c(b_row,b_col), 0)
  163. end do
  164. end do
  165. ! unregister matrices
  166. do b_row=1,NB
  167. call fstarpu_data_unregister(dh_A(b_row))
  168. end do
  169. do b_col=1,NB
  170. call fstarpu_data_unregister(dh_B(b_col))
  171. end do
  172. do b_col=1,NB
  173. do b_row=1,NB
  174. call fstarpu_data_unregister(dh_C(b_row,b_col))
  175. end do
  176. end do
  177. ! check result
  178. if (comm_rank == 0) then
  179. if (verbose) then
  180. print *,"final C"
  181. call mat_disp(C)
  182. end if
  183. do col=1,N
  184. do row=1,N
  185. if (abs(C(row,col) - 2*(row*N+col)) > 1.0) then
  186. print *, "check failed"
  187. stop 1
  188. end if
  189. end do
  190. end do
  191. end if
  192. ! free handles
  193. deallocate(dh_A)
  194. deallocate(dh_B)
  195. deallocate(dh_C)
  196. ! free matrices
  197. if (comm_rank == 0) then
  198. deallocate(A)
  199. deallocate(B)
  200. deallocate(C)
  201. end if
  202. call fstarpu_codelet_free(cl_mm)
  203. call fstarpu_shutdown()
  204. ret = fstarpu_mpi_shutdown()
  205. print *,"fstarpu_mpi_shutdown status:", ret
  206. if (ret /= 0) then
  207. stop 1
  208. end if
  209. end program nf_mm