123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382 |
- /* StarPU --- Runtime system for heterogeneous multicore architectures.
- *
- * Copyright (C) 2009-2011 Université de Bordeaux 1
- * Copyright (C) 2010 Mehdi Juhoor <mjuhoor@gmail.com>
- * Copyright (C) 2010, 2011, 2012 Centre National de la Recherche Scientifique
- *
- * StarPU is free software; you can redistribute it and/or modify
- * it under the terms of the GNU Lesser General Public License as published by
- * the Free Software Foundation; either version 2.1 of the License, or (at
- * your option) any later version.
- *
- * StarPU is distributed in the hope that it will be useful, but
- * WITHOUT ANY WARRANTY; without even the implied warranty of
- * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
- *
- * See the GNU Lesser General Public License in COPYING.LGPL for more details.
- */
- #include <starpu_mpi.h>
- #include "mpi_cholesky.h"
- #include "mpi_cholesky_models.h"
- /*
- * Create the codelets
- */
- static struct starpu_codelet cl11 =
- {
- .where = STARPU_CPU|STARPU_CUDA,
- .cpu_funcs = {chol_cpu_codelet_update_u11, NULL},
- #ifdef STARPU_USE_CUDA
- .cuda_funcs = {chol_cublas_codelet_update_u11, NULL},
- #endif
- .nbuffers = 1,
- .modes = {STARPU_RW},
- .model = &chol_model_11
- };
- static struct starpu_codelet cl21 =
- {
- .where = STARPU_CPU|STARPU_CUDA,
- .cpu_funcs = {chol_cpu_codelet_update_u21, NULL},
- #ifdef STARPU_USE_CUDA
- .cuda_funcs = {chol_cublas_codelet_update_u21, NULL},
- #endif
- .nbuffers = 2,
- .modes = {STARPU_R, STARPU_RW},
- .model = &chol_model_21
- };
- static struct starpu_codelet cl22 =
- {
- .where = STARPU_CPU|STARPU_CUDA,
- .cpu_funcs = {chol_cpu_codelet_update_u22, NULL},
- #ifdef STARPU_USE_CUDA
- .cuda_funcs = {chol_cublas_codelet_update_u22, NULL},
- #endif
- .nbuffers = 3,
- .modes = {STARPU_R, STARPU_R, STARPU_RW},
- .model = &chol_model_22
- };
- /* Returns the MPI node number where data indexes index is */
- int my_distrib(int x, int y, int nb_nodes)
- {
- return (x+y) % nb_nodes;
- }
- /*
- * code to bootstrap the factorization
- * and construct the DAG
- */
- static void dw_cholesky(float ***matA, unsigned size, unsigned ld, unsigned nblocks, int rank, int nodes)
- {
- struct timeval start;
- struct timeval end;
- starpu_data_handle_t **data_handles;
- int x, y;
- /* create all the DAG nodes */
- unsigned i,j,k;
- data_handles = malloc(nblocks*sizeof(starpu_data_handle_t *));
- for(x=0 ; x<nblocks ; x++) data_handles[x] = malloc(nblocks*sizeof(starpu_data_handle_t));
- for(x = 0; x < nblocks ; x++)
- {
- for (y = 0; y < nblocks; y++)
- {
- int mpi_rank = my_distrib(x, y, nodes);
- if (mpi_rank == rank)
- {
- //fprintf(stderr, "[%d] Owning data[%d][%d]\n", rank, x, y);
- starpu_matrix_data_register(&data_handles[x][y], 0, (uintptr_t)matA[x][y],
- ld, size/nblocks, size/nblocks, sizeof(float));
- }
- /* TODO: make better test to only registering what is needed */
- else
- {
- /* I don't own that index, but will need it for my computations */
- //fprintf(stderr, "[%d] Neighbour of data[%d][%d]\n", rank, x, y);
- starpu_matrix_data_register(&data_handles[x][y], -1, (uintptr_t)NULL,
- ld, size/nblocks, size/nblocks, sizeof(float));
- }
- if (data_handles[x][y])
- {
- starpu_data_set_rank(data_handles[x][y], mpi_rank);
- starpu_data_set_tag(data_handles[x][y], (y*nblocks)+x);
- }
- }
- }
- starpu_mpi_barrier(MPI_COMM_WORLD);
- gettimeofday(&start, NULL);
- for (k = 0; k < nblocks; k++)
- {
- int prio = STARPU_DEFAULT_PRIO;
- if (!noprio) prio = STARPU_MAX_PRIO;
- starpu_mpi_insert_task(MPI_COMM_WORLD, &cl11,
- STARPU_PRIORITY, prio,
- STARPU_RW, data_handles[k][k],
- 0);
- for (j = k+1; j<nblocks; j++)
- {
- prio = STARPU_DEFAULT_PRIO;
- if (!noprio&& (j == k+1)) prio = STARPU_MAX_PRIO;
- starpu_mpi_insert_task(MPI_COMM_WORLD, &cl21,
- STARPU_PRIORITY, prio,
- STARPU_R, data_handles[k][k],
- STARPU_RW, data_handles[k][j],
- 0);
- for (i = k+1; i<nblocks; i++)
- {
- if (i <= j)
- {
- prio = STARPU_DEFAULT_PRIO;
- if (!noprio && (i == k + 1) && (j == k +1) ) prio = STARPU_MAX_PRIO;
- starpu_mpi_insert_task(MPI_COMM_WORLD, &cl22,
- STARPU_PRIORITY, prio,
- STARPU_R, data_handles[k][i],
- STARPU_R, data_handles[k][j],
- STARPU_RW, data_handles[i][j],
- 0);
- }
- }
- }
- }
- starpu_task_wait_for_all();
- for(x = 0; x < nblocks ; x++)
- {
- for (y = 0; y < nblocks; y++)
- {
- if (data_handles[x][y])
- starpu_data_unregister(data_handles[x][y]);
- }
- free(data_handles[x]);
- }
- free(data_handles);
- starpu_mpi_barrier(MPI_COMM_WORLD);
- gettimeofday(&end, NULL);
- if (rank == 0)
- {
- double timing = (double)((end.tv_sec - start.tv_sec)*1000000 + (end.tv_usec - start.tv_usec));
- fprintf(stderr, "Computation took (in ms)\n");
- fprintf(stdout, "%2.2f\n", timing/1000);
- double flop = (1.0f*size*size*size)/3.0f;
- fprintf(stderr, "Synthetic GFlops : %2.2f\n", (flop/timing/1000.0f));
- }
- }
- int main(int argc, char **argv)
- {
- /* create a simple definite positive symetric matrix example
- *
- * Hilbert matrix : h(i,j) = 1/(i+j+1)
- * */
- float ***bmat;
- int rank, nodes;
- parse_args(argc, argv);
- struct starpu_conf conf;
- starpu_conf_init(&conf);
- conf.sched_policy_name = "heft";
- conf.calibrate = 1;
- int ret = starpu_init(&conf);
- STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");
- starpu_mpi_initialize_extended(&rank, &nodes);
- starpu_helper_cublas_init();
- unsigned i,j,x,y;
- bmat = malloc(nblocks * sizeof(float *));
- for(x=0 ; x<nblocks ; x++)
- {
- bmat[x] = malloc(nblocks * sizeof(float *));
- for(y=0 ; y<nblocks ; y++)
- {
- starpu_malloc((void **)&bmat[x][y], BLOCKSIZE*BLOCKSIZE*sizeof(float));
- for (i = 0; i < BLOCKSIZE; i++)
- {
- for (j = 0; j < BLOCKSIZE; j++)
- {
- 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);
- //mat[j +i*size] = ((i == j)?1.0f*size:0.0f);
- }
- }
- }
- }
- if (display)
- {
- printf("[%d] Input :\n", rank);
- for(y=0 ; y<nblocks ; y++)
- {
- for(x=0 ; x<nblocks ; x++)
- {
- printf("Block %d,%d :\n", x, y);
- for (j = 0; j < BLOCKSIZE; j++)
- {
- for (i = 0; i < BLOCKSIZE; i++)
- {
- if (i <= j)
- {
- printf("%2.2f\t", bmat[y][x][j +i*BLOCKSIZE]);
- }
- else
- {
- printf(".\t");
- }
- }
- printf("\n");
- }
- }
- }
- }
- dw_cholesky(bmat, size, size/nblocks, nblocks, rank, nodes);
- starpu_mpi_shutdown();
- if (display)
- {
- printf("[%d] Results :\n", rank);
- for(y=0 ; y<nblocks ; y++)
- {
- for(x=0 ; x<nblocks ; x++)
- {
- printf("Block %d,%d :\n", x, y);
- for (j = 0; j < BLOCKSIZE; j++)
- {
- for (i = 0; i < BLOCKSIZE; i++)
- {
- if (i <= j)
- {
- printf("%2.2f\t", bmat[y][x][j +i*BLOCKSIZE]);
- }
- else
- {
- printf(".\t");
- }
- }
- printf("\n");
- }
- }
- }
- }
- float *rmat = malloc(size*size*sizeof(float));
- for(x=0 ; x<nblocks ; x++)
- {
- for(y=0 ; y<nblocks ; y++)
- {
- for (i = 0; i < BLOCKSIZE; i++)
- {
- for (j = 0; j < BLOCKSIZE; j++)
- {
- rmat[j+(y*BLOCKSIZE)+(i+(x*BLOCKSIZE))*size] = bmat[x][y][j +i*BLOCKSIZE];
- }
- }
- }
- }
- fprintf(stderr, "[%d] compute explicit LLt ...\n", rank);
- for (j = 0; j < size; j++)
- {
- for (i = 0; i < size; i++)
- {
- if (i > j)
- {
- rmat[j+i*size] = 0.0f; // debug
- }
- }
- }
- float *test_mat = malloc(size*size*sizeof(float));
- STARPU_ASSERT(test_mat);
- SSYRK("L", "N", size, size, 1.0f,
- rmat, size, 0.0f, test_mat, size);
- fprintf(stderr, "[%d] comparing results ...\n", rank);
- if (display)
- {
- for (j = 0; j < size; j++)
- {
- for (i = 0; i < size; i++)
- {
- if (i <= j)
- {
- printf("%2.2f\t", test_mat[j +i*size]);
- }
- else
- {
- printf(".\t");
- }
- }
- printf("\n");
- }
- }
- int correctness = 1;
- for(x = 0; x < nblocks ; x++)
- {
- for (y = 0; y < nblocks; y++)
- {
- int mpi_rank = my_distrib(x, y, nodes);
- if (mpi_rank == rank)
- {
- for (i = (size/nblocks)*x ; i < (size/nblocks)*x+(size/nblocks); i++)
- {
- for (j = (size/nblocks)*y ; j < (size/nblocks)*y+(size/nblocks); j++)
- {
- if (i <= j)
- {
- float orig = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);
- float err = abs(test_mat[j +i*size] - orig);
- if (err > 0.00001)
- {
- fprintf(stderr, "[%d] Error[%d, %d] --> %2.2f != %2.2f (err %2.2f)\n", rank, i, j, test_mat[j +i*size], orig, err);
- correctness = 0;
- break;
- }
- }
- }
- }
- }
- }
- }
- for(x=0 ; x<nblocks ; x++)
- {
- for(y=0 ; y<nblocks ; y++)
- {
- starpu_free((void *)bmat[x][y]);
- }
- free(bmat[x]);
- }
- free(bmat);
- free(rmat);
- free(test_mat);
- starpu_helper_cublas_shutdown();
- starpu_shutdown();
- assert(correctness);
- return 0;
- }
|