cholesky_common.jl 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170
  1. # StarPU --- Runtime system for heterogeneous multicore architectures.
  2. #
  3. # Copyright (C) 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. # Standard kernels for the Cholesky factorization
  17. # U22 is the gemm update
  18. # U21 is the trsm update
  19. # U11 is the cholesky factorization
  20. @target STARPU_CPU+STARPU_CUDA
  21. @codelet function u11(sub11 :: Matrix{Float32}) :: Nothing
  22. nx :: Int32 = width(sub11)
  23. ld :: Int32 = ld(sub11)
  24. for z in 0:nx-1
  25. lambda11 :: Float32 = sqrt(sub11[z+1,z+1])
  26. sub11[z+1,z+1] = lambda11
  27. alpha ::Float32 = 1.0f0 / lambda11
  28. X :: Vector{Float32} = view(sub11, z+2:z+2+(nx-z-2), z+1)
  29. STARPU_SSCAL(nx-z-1, alpha, X, 1)
  30. alpha = -1.0f0
  31. A :: Matrix{Float32} = view(sub11, z+2:z+2+(nx-z-2), z+2:z+2+(nx-z-2))
  32. STARPU_SSYR("L", nx-z-1, alpha, X, 1, A, ld)
  33. end
  34. return
  35. end
  36. @target STARPU_CPU+STARPU_CUDA
  37. @codelet function u21(sub11 :: Matrix{Float32},
  38. sub21 :: Matrix{Float32}) :: Nothing
  39. ld11 :: Int32 = ld(sub11)
  40. ld21 :: Int32 = ld(sub21)
  41. nx21 :: Int32 = width(sub21)
  42. ny21 :: Int32 = height(sub21)
  43. alpha :: Float32 = 1.0f0
  44. STARPU_STRSM("R", "L", "T", "N", nx21, ny21, alpha, sub11, ld11, sub21, ld21)
  45. return
  46. end
  47. @target STARPU_CPU+STARPU_CUDA
  48. @codelet function u22(left :: Matrix{Float32},
  49. right :: Matrix{Float32},
  50. center :: Matrix{Float32}) :: Nothing
  51. dx :: Int32 = width(center)
  52. dy :: Int32 = height(center)
  53. dz :: Int32 = width(left)
  54. ld21 :: Int32 = ld(left)
  55. ld12 :: Int32 = ld(center)
  56. ld22 :: Int32 = ld(right)
  57. alpha :: Float32 = -1.0f0
  58. beta :: Float32 = 1.0f0
  59. STARPU_SGEMM("N", "T", dy, dx, dz, alpha, left, ld21, right, ld12, beta, center, ld22)
  60. return
  61. end
  62. @inline function tag11(k)
  63. return starpu_tag_t((UInt64(1)<<60) | UInt64(k))
  64. end
  65. @inline function tag21(k, j)
  66. return starpu_tag_t((UInt64(3)<<60) | (UInt64(k)<<32) | UInt64(j))
  67. end
  68. @inline function tag22(k, i, j)
  69. return starpu_tag_t((UInt64(4)<<60) | (UInt64(k)<<32) | (UInt64(i)<<16) | UInt64(j))
  70. end
  71. function check(mat::Matrix{Float32})
  72. size_p = size(mat, 1)
  73. for i in 1:size_p
  74. for j in 1:size_p
  75. if j > i
  76. mat[i, j] = 0.0f0
  77. end
  78. end
  79. end
  80. test_mat ::Matrix{Float32} = zeros(Float32, size_p, size_p)
  81. syrk!('L', 'N', 1.0f0, mat, 0.0f0, test_mat)
  82. for i in 1:size_p
  83. for j in 1:size_p
  84. if j <= i
  85. orig = (1.0f0/(1.0f0+(i-1)+(j-1))) + ((i == j) ? 1.0f0*size_p : 0.0f0)
  86. err = abs(test_mat[i,j] - orig) / orig
  87. if err > 0.0001
  88. got = test_mat[i,j]
  89. expected = orig
  90. error("[$i, $j] -> $got != $expected (err $err)")
  91. end
  92. end
  93. end
  94. end
  95. println(stderr, "Verification successful !")
  96. end
  97. function clean_tags(nblocks)
  98. for k in 1:nblocks
  99. starpu_tag_remove(tag11(k))
  100. for m in k+1:nblocks
  101. starpu_tag_remove(tag21(k, m))
  102. for n in k+1:nblocks
  103. if n <= m
  104. starpu_tag_remove(tag22(k, m, n))
  105. end
  106. end
  107. end
  108. end
  109. end
  110. function main(size_p :: Int, nblocks :: Int; verify = false, verbose = false)
  111. mat :: Matrix{Float32} = zeros(Float32, size_p, size_p)
  112. # create a simple definite positive symetric matrix
  113. # Hilbert matrix h(i,j) = 1/(i+j+1)
  114. for i in 1:size_p
  115. for j in 1:size_p
  116. mat[i, j] = 1.0f0 / (1.0f0+(i-1)+(j-1)) + ((i == j) ? 1.0f0*size_p : 0.0f0)
  117. end
  118. end
  119. if verbose
  120. display(mat)
  121. end
  122. starpu_memory_pin(mat)
  123. t_start = time_ns()
  124. cholesky(mat, size_p, nblocks)
  125. t_end = time_ns()
  126. starpu_memory_unpin(mat)
  127. flop = (1.0*size_p*size_p*size_p)/3.0
  128. time_ms = (t_end-t_start) / 1e6
  129. gflops = flop/(time_ms*1000)/1000
  130. println("$size_p\t$time_ms\t$gflops")
  131. clean_tags(nblocks)
  132. if verbose
  133. display(mat)
  134. end
  135. if verify
  136. check(mat)
  137. end
  138. end