nf_mm.f90 7.5 KB

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