26spm_example_fill_ijv_spm(
spm_int_t *colptr,
45 for( i = 0; i < dim1; i++ )
47 for( j = 0; j < dim2; j++ )
49 for( k = 0; k < dim3; k++ )
55 rowptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
56 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
64 if ( i == 0 ) { value -= 1.; }
65 if ( i == mdim1 ) { value -= 1.; }
66 if ( j == 0 ) { value -= 1.; }
67 if ( j == mdim2 ) { value -= 1.; }
68 if ( k == 0 ) { value -= 1.; }
69 if ( k == mdim3 ) { value -= 1.; }
74 for( d=0; d<(dof*dof); d++ ) {
82 rowptr[l] = (i+1) + dim1 * j + dim1 * dim2 * k + 1;
83 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
86 for( d=0; d<(dof*dof); d++ ) {
93 rowptr[l] = i + dim1 * (j+1) + dim1 * dim2 * k + 1;
94 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
97 for( d=0; d<(dof*dof); d++ ) {
104 rowptr[l] = i + dim1 * j + dim1 * dim2 * (k+1) + 1;
105 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
108 for( d=0; d<(dof*dof); d++ ) {
117 assert( l == ((2*(dim1)-1) * dim2 * dim3 + (dim2-1)*dim1*dim3 + dim2*dim1*(dim3-1)) );
124spm_example_create_laplacian(
spmatrix_t *spm )
136 n = dim1 * dim2 * dim3;
137 nnz = dim1 * dim2 * dim3;
138 nnz += (dim1 - 1) * dim2 * dim3;
139 nnz += dim1 * (dim2 - 1) * dim3;
140 nnz += dim1 * dim2 * (dim3 - 1);
145 colptr = malloc( nnz *
sizeof(
spm_int_t) );
146 rowptr = malloc( nnz *
sizeof(
spm_int_t) );
147 values = malloc( nnz * dof * dof *
sizeof(
double) );
152 spm_example_fill_ijv_spm( colptr, rowptr, values,
153 dim1, dim2, dim3, dof );
185int main (
int argc,
char **argv)
191 double epsilon, norm;
195#if defined(SPM_WITH_MPI)
196 MPI_Init( &argc, &argv );
202 spm_example_create_laplacian( &spm );
218 fprintf( stdout,
" || A ||_f: %e\n", norm );
236 rc =
spmGenMat( SpmRhsRndX, nrhs, &spm, &alpha, 24356, x, ldx );
261 epsilon = epsilon * spm.
nnzexp / spm.
gN;
262 rc =
spmCheckAxb( epsilon, nrhs, &spm, NULL, 1, b, ldb, x, ldx );
270#if defined(SPM_WITH_MPI)
static size_t spm_size_of(spm_coeftype_t type)
Double datatype that is not converted through precision generator functions.
void spmExit(spmatrix_t *spm)
Cleanup the spm structure but do not free the spm pointer.
void spmUpdateComputedFields(spmatrix_t *spm)
Update all the computed fields based on the static values stored.
int spm_int_t
The main integer datatype used in spm arrays.
int spmCheckAxb(double eps, spm_int_t nrhs, const spmatrix_t *spm, void *opt_X0, spm_int_t opt_ldx0, void *B, spm_int_t ldb, const void *X, spm_int_t ldx)
Check the backward error, and the forward error if x0 is provided.
int spmConvert(int ofmttype, spmatrix_t *ospm)
Convert the storage format of the spm.
int spmMatMat(spm_trans_t trans, spm_int_t n, double alpha, const spmatrix_t *A, const void *B, spm_int_t ldb, double beta, void *C, spm_int_t ldc)
Compute a matrix-matrix product.
void spmScal(double alpha, spmatrix_t *spm)
Scale the spm.
double spmNorm(spm_normtype_t ntype, const spmatrix_t *spm)
Compute the norm of the spm.
void spmInit(spmatrix_t *spm)
Init the spm structure.
int spmGenMat(spm_rhstype_t type, spm_int_t nrhs, const spmatrix_t *spm, void *alpha, unsigned long long int seed, void *A, spm_int_t lda)
Generate a set of vectors associated to a given matrix.
void spmPrintInfo(const spmatrix_t *spm, FILE *f)
Print basic informations about the spm matrix into a given stream.
The sparse matrix data structure.