mpi_cholesky.c 3.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009-2012 Université de Bordeaux 1
  4. * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
  5. * Copyright (C) 2010, 2011, 2012, 2013 Centre National de la Recherche Scientifique
  6. *
  7. * StarPU is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU Lesser General Public License as published by
  9. * the Free Software Foundation; either version 2.1 of the License, or (at
  10. * your option) any later version.
  11. *
  12. * StarPU is distributed in the hope that it will be useful, but
  13. * WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  15. *
  16. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  17. */
  18. #include <starpu_mpi.h>
  19. #include "mpi_cholesky_params.h"
  20. #include "mpi_cholesky_models.h"
  21. #include "mpi_cholesky_codelets.h"
  22. void display_matrix(float ***bmat, int rank)
  23. {
  24. unsigned i,j,x,y;
  25. if (display)
  26. {
  27. printf("[%d] Input :\n", rank);
  28. for(y=0 ; y<nblocks ; y++)
  29. {
  30. for(x=0 ; x<nblocks ; x++)
  31. {
  32. printf("Block %u,%u :\n", x, y);
  33. for (j = 0; j < BLOCKSIZE; j++)
  34. {
  35. for (i = 0; i < BLOCKSIZE; i++)
  36. {
  37. if (i <= j)
  38. {
  39. printf("%2.2f\t", bmat[y][x][j +i*BLOCKSIZE]);
  40. }
  41. else
  42. {
  43. printf(".\t");
  44. }
  45. }
  46. printf("\n");
  47. }
  48. }
  49. }
  50. }
  51. }
  52. int main(int argc, char **argv)
  53. {
  54. /* create a simple definite positive symetric matrix example
  55. *
  56. * Hilbert matrix : h(i,j) = 1/(i+j+1)
  57. * */
  58. float ***bmat;
  59. int rank, nodes, ret;
  60. unsigned i,j,x,y;
  61. ret = starpu_init(NULL);
  62. STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
  63. ret = starpu_mpi_init(&argc, &argv, 1);
  64. STARPU_CHECK_RETURN_VALUE(ret, "starpu_mpi_init");
  65. MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  66. MPI_Comm_size(MPI_COMM_WORLD, &nodes);
  67. starpu_helper_cublas_init();
  68. parse_args(argc, argv, nodes);
  69. bmat = malloc(nblocks * sizeof(float *));
  70. for(x=0 ; x<nblocks ; x++)
  71. {
  72. bmat[x] = malloc(nblocks * sizeof(float *));
  73. for(y=0 ; y<nblocks ; y++)
  74. {
  75. starpu_malloc((void **)&bmat[x][y], BLOCKSIZE*BLOCKSIZE*sizeof(float));
  76. for (i = 0; i < BLOCKSIZE; i++)
  77. {
  78. for (j = 0; j < BLOCKSIZE; j++)
  79. {
  80. bmat[x][y][j +i*BLOCKSIZE] = (1.0f/(1.0f+(i+(x*BLOCKSIZE)+j+(y*BLOCKSIZE)))) + ((i+(x*BLOCKSIZE) == j+(y*BLOCKSIZE))?1.0f*size:0.0f);
  81. //mat[j +i*size] = ((i == j)?1.0f*size:0.0f);
  82. }
  83. }
  84. }
  85. }
  86. display_matrix(bmat, rank);
  87. double timing, flops;
  88. dw_cholesky(bmat, size, size/nblocks, nblocks, rank, nodes, &timing, &flops);
  89. starpu_mpi_shutdown();
  90. display_matrix(bmat, rank);
  91. int correctness;
  92. dw_cholesky_check_computation(bmat, size, rank, nodes, &correctness, &flops);
  93. for(x=0 ; x<nblocks ; x++)
  94. {
  95. for(y=0 ; y<nblocks ; y++)
  96. {
  97. starpu_free((void *)bmat[x][y]);
  98. }
  99. free(bmat[x]);
  100. }
  101. free(bmat);
  102. starpu_helper_cublas_shutdown();
  103. starpu_shutdown();
  104. assert(correctness);
  105. if (rank == 0)
  106. {
  107. fprintf(stdout, "Computation time (in ms): %2.2f\n", timing/1000);
  108. fprintf(stdout, "Synthetic GFlops : %2.2f\n", (flops/timing/1000.0f));
  109. }
  110. return 0;
  111. }