| 123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370 | /* StarPU --- Runtime system for heterogeneous multicore architectures. * * Copyright (C) 2009-2013  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 "cholesky.h"#include "../sched_ctx_utils/sched_ctx_utils.h"/* *	Create the codelets */static struct starpu_codelet cl11 ={	.type = STARPU_SEQ,	.cpu_funcs = {chol_cpu_codelet_update_u11, NULL},#ifdef STARPU_USE_CUDA	.cuda_funcs = {chol_cublas_codelet_update_u11, NULL},#elif defined(STARPU_SIMGRID)	.cuda_funcs = {(void*)1, NULL},#endif	.nbuffers = 1,	.modes = {STARPU_RW},	.model = &chol_model_11};static struct starpu_codelet cl21 ={	.type = STARPU_SEQ,	.cpu_funcs = {chol_cpu_codelet_update_u21, NULL},#ifdef STARPU_USE_CUDA	.cuda_funcs = {chol_cublas_codelet_update_u21, NULL},#elif defined(STARPU_SIMGRID)	.cuda_funcs = {(void*)1, NULL},#endif	.nbuffers = 2,	.modes = {STARPU_R, STARPU_RW},	.model = &chol_model_21};static struct starpu_codelet cl22 ={	.type = STARPU_SEQ,	.max_parallelism = INT_MAX,	.cpu_funcs = {chol_cpu_codelet_update_u22, NULL},#ifdef STARPU_USE_CUDA	.cuda_funcs = {chol_cublas_codelet_update_u22, NULL},#elif defined(STARPU_SIMGRID)	.cuda_funcs = {(void*)1, NULL},#endif	.nbuffers = 3,	.modes = {STARPU_R, STARPU_R, STARPU_RW},	.model = &chol_model_22};/* *	code to bootstrap the factorization *	and construct the DAG */static void callback_turn_spmd_on(void *arg STARPU_ATTRIBUTE_UNUSED){	cl22.type = STARPU_SPMD;}static int _cholesky(starpu_data_handle_t dataA, unsigned nblocks){	int ret;	double start;	double end;	unsigned i,j,k;	unsigned long n = starpu_matrix_get_nx(dataA);	unsigned long nn = n/nblocks;	int prio_level = noprio?STARPU_DEFAULT_PRIO:STARPU_MAX_PRIO;	start = starpu_timing_now();	if (bound || bound_lp || bound_mps)		starpu_bound_start(bound_deps, 0);	/* create all the DAG nodes */	for (k = 0; k < nblocks; k++)	{                starpu_data_handle_t sdatakk = starpu_data_get_sub_data(dataA, 2, k, k);                ret = starpu_insert_task(&cl11,					 STARPU_PRIORITY, prio_level,					 STARPU_RW, sdatakk,					 STARPU_CALLBACK, (k == 3*nblocks/4)?callback_turn_spmd_on:NULL,					 STARPU_FLOPS, (double) FLOPS_SPOTRF(nn),					 0);		if (ret == -ENODEV) return 77;		STARPU_CHECK_RETURN_VALUE(ret, "starpu_insert_task");		for (j = k+1; j<nblocks; j++)		{                        starpu_data_handle_t sdatakj = starpu_data_get_sub_data(dataA, 2, k, j);                        ret = starpu_insert_task(&cl21,						 STARPU_PRIORITY, (j == k+1)?prio_level:STARPU_DEFAULT_PRIO,						 STARPU_R, sdatakk,						 STARPU_RW, sdatakj,						 STARPU_FLOPS, (double) FLOPS_STRSM(nn, nn),						 0);			if (ret == -ENODEV) return 77;			STARPU_CHECK_RETURN_VALUE(ret, "starpu_insert_task");			for (i = k+1; i<nblocks; i++)			{				if (i <= j)                                {					starpu_data_handle_t sdataki = starpu_data_get_sub_data(dataA, 2, k, i);					starpu_data_handle_t sdataij = starpu_data_get_sub_data(dataA, 2, i, j);					ret = starpu_insert_task(&cl22,								 STARPU_PRIORITY, ((i == k+1) && (j == k+1))?prio_level:STARPU_DEFAULT_PRIO,								 STARPU_R, sdataki,								 STARPU_R, sdatakj,								 STARPU_RW, sdataij,								 STARPU_FLOPS, (double) FLOPS_SGEMM(nn, nn, nn),								 0);					if (ret == -ENODEV) return 77;					STARPU_CHECK_RETURN_VALUE(ret, "starpu_insert_task");                                }			}		}	}	starpu_task_wait_for_all();	if (bound || bound_lp || bound_mps)		starpu_bound_stop();	end = starpu_timing_now();	double timing = end - start;	double flop = FLOPS_SPOTRF(n);	if(with_ctxs || with_noctxs || chole1 || chole2)		update_sched_ctx_timing_results((flop/timing/1000.0f), (timing/1000000.0f));	else	{		FPRINTF(stderr, "Computation took (in ms)\n");		FPRINTF(stdout, "%2.2f\n", timing/1000);			FPRINTF(stderr, "Synthetic GFlops : %2.2f\n", (flop/timing/1000.0f));		if (bound_lp)		{			FILE *f = fopen("cholesky.lp", "w");			starpu_bound_print_lp(f);		}		if (bound_mps)		{			FILE *f = fopen("cholesky.mps", "w");			starpu_bound_print_mps(f);		}		if (bound)		{			double res;			starpu_bound_compute(&res, NULL, 0);			FPRINTF(stderr, "Theoretical GFlops: %2.2f\n", (flop/res/1000000.0f));		}	}	return 0;}static int cholesky(float *matA, unsigned size, unsigned ld, unsigned nblocks){	starpu_data_handle_t dataA;	/* monitor and partition the A matrix into blocks :	 * one block is now determined by 2 unsigned (i,j) */	starpu_matrix_data_register(&dataA, STARPU_MAIN_RAM, (uintptr_t)matA, ld, size, size, sizeof(float));	struct starpu_data_filter f =	{		.filter_func = starpu_matrix_filter_vertical_block,		.nchildren = nblocks	};	struct starpu_data_filter f2 =	{		.filter_func = starpu_matrix_filter_block,		.nchildren = nblocks	};	starpu_data_map_filters(dataA, 2, &f, &f2);	int ret = _cholesky(dataA, nblocks);	starpu_data_unpartition(dataA, STARPU_MAIN_RAM);	starpu_data_unregister(dataA);	return ret;}static void execute_cholesky(unsigned size, unsigned nblocks){	int ret;	float *mat = NULL;	unsigned i,j;#ifndef STARPU_SIMGRID	starpu_malloc((void **)&mat, (size_t)size*size*sizeof(float));	for (i = 0; i < size; i++)	{		for (j = 0; j < size; j++)		{			mat[j +i*size] = (1.0f/(1.0f+i+j)) + ((i == j)?1.0f*size:0.0f);			/* mat[j +i*size] = ((i == j)?1.0f*size:0.0f); */		}	}#endif/* #define PRINT_OUTPUT */#ifdef PRINT_OUTPUT	FPRINTF(stdout, "Input :\n");	for (j = 0; j < size; j++)	{		for (i = 0; i < size; i++)		{			if (i <= j)			{				FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);			}			else			{				FPRINTF(stdout, ".\t");			}		}		FPRINTF(stdout, "\n");	}#endif	ret = cholesky(mat, size, size, nblocks);#ifdef PRINT_OUTPUT	FPRINTF(stdout, "Results :\n");	for (j = 0; j < size; j++)	{		for (i = 0; i < size; i++)		{			if (i <= j)			{				FPRINTF(stdout, "%2.2f\t", mat[j +i*size]);			}			else			{				FPRINTF(stdout, ".\t");				mat[j+i*size] = 0.0f; /* debug */			}		}		FPRINTF(stdout, "\n");	}#endif	if (check)	{		FPRINTF(stderr, "compute explicit LLt ...\n");		for (j = 0; j < size; j++)		{			for (i = 0; i < size; i++)			{				if (i > j)				{					mat[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,					mat, size, 0.0f, test_mat, size);		FPRINTF(stderr, "comparing results ...\n");#ifdef PRINT_OUTPUT		for (j = 0; j < size; j++)		{			for (i = 0; i < size; i++)			{				if (i <= j)				{					FPRINTF(stdout, "%2.2f\t", test_mat[j +i*size]);				}				else				{					FPRINTF(stdout, ".\t");				}			}			FPRINTF(stdout, "\n");		}#endif		for (j = 0; j < size; j++)		{			for (i = 0; i < size; i++)			{				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, "Error[%u, %u] --> %2.2f != %2.2f (err %2.2f)\n", i, j, test_mat[j +i*size], orig, err);	                                        assert(0);	                                }	                        }			}	        }		free(test_mat);	}	starpu_free(mat);}int main(int argc, char **argv){	/* create a simple definite positive symetric matrix example	 *	 *	Hilbert matrix : h(i,j) = 1/(i+j+1)	 * */	parse_args(argc, argv);	if(with_ctxs || with_noctxs || chole1 || chole2)		parse_args_ctx(argc, argv);	int ret;	ret = starpu_init(NULL);	if (ret == -ENODEV)                return 77;        STARPU_CHECK_RETURN_VALUE(ret, "starpu_init");	starpu_cublas_init();	if(with_ctxs)	{		construct_contexts(execute_cholesky);		start_2benchs(execute_cholesky);	}	else if(with_noctxs)		start_2benchs(execute_cholesky);	else if(chole1)		start_1stbench(execute_cholesky);	else if(chole2)		start_2ndbench(execute_cholesky);	else		execute_cholesky(size, nblocks);	starpu_cublas_shutdown();	starpu_shutdown();	return ret;}
 |