void lis_sort_id(LIS_INT is, LIS_INT ie, LIS_INT *i1, LIS_SCALAR *d1) { LIS_INT i,j; LIS_INT p,t1; LIS_INT v; LIS_SCALAR s1; if( ie <= is ) return; p = (is+ie)/2; v = i1[p]; lis_swap(i1[p],i1[ie],t1); lis_swap(d1[p],d1[ie],s1); i = is; j = ie; while(i<=j) { while(i1[i] < v) { i++; } while(i1[j] > v) { j--; } if( i<=j ) { lis_swap(i1[i],i1[j],t1); lis_swap(d1[i],d1[j],s1); i++; j--; } } lis_sort_id(is,j ,i1,d1); lis_sort_id(i ,ie,i1,d1); }
LIS_INT lis_matrix_sort_csr(LIS_MATRIX A) { LIS_INT i,n; LIS_DEBUG_FUNC_IN; if( !A->is_sorted ) { n = A->n; if( A->is_splited ) { #ifdef _OPENMP #pragma omp parallel for private(i) #endif for(i=0;i<n;i++) { lis_sort_id(A->L->ptr[i],A->L->ptr[i+1]-1,A->L->index,A->L->value); lis_sort_id(A->U->ptr[i],A->U->ptr[i+1]-1,A->U->index,A->U->value); } } else { #ifdef _OPENMP #pragma omp parallel for private(i) #endif for(i=0;i<n;i++) { lis_sort_id(A->ptr[i],A->ptr[i+1]-1,A->index,A->value); } } A->is_sorted = LIS_TRUE; } LIS_DEBUG_FUNC_OUT; return LIS_SUCCESS; }
LIS_INT main(LIS_INT argc, char* argv[]) { LIS_MATRIX A,A0; LIS_VECTOR b,x; LIS_SCALAR *value; LIS_INT nprocs,my_rank; int int_nprocs,int_my_rank; LIS_INT nthreads,maxthreads; LIS_INT gn,nnz,np; LIS_INT i,j,k,si,sj,sk,ii,jj,ctr; LIS_INT l,m,n,nn; LIS_INT is,ie; LIS_INT err,iter,matrix_type,storage,ss,se; LIS_INT *ptr,*index; double time,time2,nnzs,nnzap,nnzt; LIS_SCALAR val; double commtime,comptime,flops; LIS_DEBUG_FUNC_IN; lis_initialize(&argc, &argv); #ifdef USE_MPI MPI_Comm_size(MPI_COMM_WORLD,&int_nprocs); MPI_Comm_rank(MPI_COMM_WORLD,&int_my_rank); nprocs = int_nprocs; my_rank = int_my_rank; #else nprocs = 1; my_rank = 0; #endif if( argc < 5 ) { if( my_rank==0 ) { printf("Usage: %s l m n iter [matrix_type]\n", argv[0]); } CHKERR(1); } l = atoi(argv[1]); m = atoi(argv[2]); n = atoi(argv[3]); iter = atoi(argv[4]); if (argv[5] == NULL) { storage = 0; } else { storage = atoi(argv[5]); } if( iter<=0 ) { #ifdef _LONG__LONG if( my_rank==0 ) printf("iter=%lld <= 0\n",iter); #else if( my_rank==0 ) printf("iter=%d <= 0\n",iter); #endif CHKERR(1); } if( l<=0 || m<=0 || n<=0 ) { #ifdef _LONG__LONG if( my_rank==0 ) printf("l=%lld <=0, m=%lld <=0 or n=%lld <=0\n",l,m,n); #else if( my_rank==0 ) printf("l=%d <=0, m=%d <=0 or n=%d <=0\n",l,m,n); #endif CHKERR(1); } if( storage<0 || storage>11 ) { #ifdef _LONG__LONG if( my_rank==0 ) printf("matrix_type=%lld < 0 or matrix_type=%lld > 11\n",storage,storage); #else if( my_rank==0 ) printf("matrix_type=%d < 0 or matrix_type=%d > 11\n",storage,storage); #endif CHKERR(1); } if( my_rank==0 ) { printf("\n"); #ifdef _LONG__LONG printf("number of processes = %lld\n",nprocs); #else printf("number of processes = %d\n",nprocs); #endif } #ifdef _OPENMP nthreads = omp_get_num_procs(); maxthreads = omp_get_max_threads(); if( my_rank==0 ) { #ifdef _LONG__LONG printf("max number of threads = %lld\n", nthreads); printf("number of threads = %lld\n", maxthreads); #else printf("max number of threads = %d\n", nthreads); printf("number of threads = %d\n", maxthreads); #endif } #else nthreads = 1; maxthreads = 1; #endif /* create matrix and vectors */ nn = l*m*n; err = lis_matrix_create(LIS_COMM_WORLD,&A0); err = lis_matrix_set_size(A0,0,nn); CHKERR(err); ptr = (LIS_INT *)malloc((A0->n+1)*sizeof(LIS_INT)); if( ptr==NULL ) CHKERR(1); index = (LIS_INT *)malloc(27*A0->n*sizeof(LIS_INT)); if( index==NULL ) CHKERR(1); value = (LIS_SCALAR *)malloc(27*A0->n*sizeof(LIS_SCALAR)); if( value==NULL ) CHKERR(1); lis_matrix_get_range(A0,&is,&ie); ctr = 0; for(ii=is;ii<ie;ii++) { i = ii/(m*n); j = (ii - i*m*n)/n; k = ii - i*m*n - j*n; for(si=-1;si<=1;si++) { if( i+si>-1 && i+si<l ) { for(sj=-1;sj<=1;sj++) { if( j+sj>-1 && j+sj<m ) { for(sk=-1;sk<=1;sk++) { if( k+sk>-1 && k+sk<n ) { jj = ii + si*m*n + sj*n + sk; index[ctr] = jj; if( jj==ii ) { value[ctr++] = 26.0;} else { value[ctr++] = -1.0;} } } } } } } ptr[ii-is+1] = ctr; } ptr[0] = 0; err = lis_matrix_set_csr(ptr[ie-is],ptr,index,value,A0); CHKERR(err); err = lis_matrix_assemble(A0); CHKERR(err); n = A0->n; gn = A0->gn; nnz = A0->nnz; np = A0->np-n; #ifdef USE_MPI MPI_Allreduce(&nnz,&i,1,LIS_MPI_INT,MPI_SUM,A0->comm); nnzap = (double)i / (double)nprocs; nnzt = ((double)nnz -nnzap)*((double)nnz -nnzap); nnz = i; MPI_Allreduce(&nnzt,&nnzs,1,MPI_DOUBLE,MPI_SUM,A0->comm); nnzs = (nnzs / (double)nprocs)/nnzap; MPI_Allreduce(&np,&i,1,LIS_MPI_INT,MPI_SUM,A0->comm); np = i; #endif if( my_rank==0 ) { #ifdef _LONG__LONG printf("matrix size = %lld x %lld (%lld nonzero entries)\n",gn,gn,nnz); printf("number of iterations = %lld\n\n",iter); #else printf("matrix size = %d x %d (%d nonzero entries)\n",gn,gn,nnz); printf("number of iterations = %d\n\n",iter); #endif } err = lis_vector_duplicate(A0,&x); if( err ) CHKERR(err); err = lis_vector_duplicate(A0,&b); if( err ) CHKERR(err); lis_matrix_get_range(A0,&is,&ie); for(i=0;i<n;i++) { err = lis_vector_set_value(LIS_INS_VALUE,i+is,1.0,x); } for(i=0;i<n;i++) { lis_sort_id(A0->ptr[i],A0->ptr[i+1]-1,A0->index,A0->value); } /* MPI version of VBR is not implemented. DNS is also excluded to reduce memory usage. */ if (storage==0) { ss = 1; se = 11; } else { ss = storage; se = storage+1; } for (matrix_type=ss;matrix_type<se;matrix_type++) { if ( nprocs>1 && matrix_type==9 ) continue; lis_matrix_duplicate(A0,&A); lis_matrix_set_type(A,matrix_type); err = lis_matrix_convert(A0,A); if( err ) CHKERR(err); comptime = 0.0; commtime = 0.0; for(i=0;i<iter;i++) { #ifdef USE_MPI MPI_Barrier(A->comm); time = lis_wtime(); lis_send_recv(A->commtable,x->value); commtime += lis_wtime() - time; #endif time2 = lis_wtime(); lis_matvec(A,x,b); comptime += lis_wtime() - time2; } lis_vector_nrm2(b,&val); if( my_rank==0 ) { flops = 2.0*nnz*iter*1.0e-6 / comptime; #ifdef USE_MPI #ifdef _LONG__DOUBLE #ifdef _LONG__LONG printf("matrix_type = %2lld (%s), computation = %e sec, %8.3f MFLOPS, communication = %e sec, communication/computation = %3.3f %%, 2-norm = %Le\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,commtime,commtime/comptime*100,val); #else printf("matrix_type = %2d (%s), computation = %e sec, %8.3f MFLOPS, communication = %e sec, communication/computation = %3.3f %%, 2-norm = %Le\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,commtime,commtime/comptime*100,val); #endif #else #ifdef _LONG__LONG printf("matrix_type = %2lld (%s), computation = %e sec, %8.3f MFLOPS, communication = %e sec, communication/computation = %3.3f %%, 2-norm = %e\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,commtime,commtime/comptime*100,val); #else printf("matrix_type = %2d (%s), computation = %e sec, %8.3f MFLOPS, communication = %e sec, communication/computation = %3.3f %%, 2-norm = %e\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,commtime,commtime/comptime*100,val); #endif #endif #else #ifdef _LONG__DOUBLE #ifdef _LONG__LONG printf("matrix_type = %2lld (%s), computation = %e sec, %8.3f MFLOPS, 2-norm = %Le\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,val); #else printf("matrix_type = %2d (%s), computation = %e sec, %8.3f MFLOPS, 2-norm = %Le\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,val); #endif #else #ifdef _LONG__LONG printf("matrix_type = %2lld (%s), computation = %e sec, %8.3f MFLOPS, 2-norm = %e\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,val); #else printf("matrix_type = %2d (%s), computation = %e sec, %8.3f MFLOPS, 2-norm = %e\n",matrix_type,lis_storagename2[matrix_type-1],comptime,flops,val); #endif #endif #endif } lis_matrix_destroy(A); } lis_matrix_destroy(A0); lis_vector_destroy(b); lis_vector_destroy(x); lis_finalize(); LIS_DEBUG_FUNC_OUT; return 0; }