SpM Handbook 1.2.4
Loading...
Searching...
No Matches
example_mdof2.c
Go to the documentation of this file.
1/**
2 *
3 * @file example_mdof2.c
4 *
5 * Example to show how to use the SPM library with a 1-based constant multi-dof
6 * sparse matrix allocated through the library but initialized by the user.
7 *
8 * @copyright 2020-2024 Bordeaux INP, CNRS (LaBRI UMR 5800), Inria,
9 * Univ. Bordeaux. All rights reserved.
10 *
11 * @version 1.2.4
12 * @author Mathieu Faverge
13 * @author Tony Delarue
14 * @date 2024-06-26
15 *
16 * @ingroup examples_c
17 * @code
18 **/
19#include <spm.h>
20
21/**
22 * @brief Create an IJV symmetric laplacian matrix for a regular cube of
23 * dimensions (dim1 x dim2 x dim3) and with 1-based values as in Fortran arrays.
24 */
25static inline void
26spm_example_fill_ijv_spm( spm_int_t *colptr,
27 spm_int_t *rowptr,
28 double *values,
29 spm_int_t dim1,
30 spm_int_t dim2,
31 spm_int_t dim3,
32 spm_int_t dof )
33{
34 spm_int_t l, i, j, k, v, d;
35 spm_int_t mdim1, mdim2, mdim3;
36 double alpha = 8.;
37 double beta = 0.5;
38 double value;
39
40 mdim1 = dim1 - 1;
41 mdim2 = dim2 - 1;
42 mdim3 = dim3 - 1;
43 l = 0;
44 v = 0;
45 for( i = 0; i < dim1; i++ )
46 {
47 for( j = 0; j < dim2; j++ )
48 {
49 for( k = 0; k < dim3; k++ )
50 {
51 /*
52 * Start with the diagonal element connected with its two
53 * neigbours in the three directions
54 */
55 rowptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
56 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
57 l++;
58
59 /*
60 * Fill the element matrix
61 */
62 {
63 value = 6.;
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.; }
70
71 /* Scale by alpha */
72 value *= alpha;
73
74 for( d=0; d<(dof*dof); d++ ) {
75 values[v] = value;
76 v++;
77 }
78 }
79
80 /* Add connexions */
81 if ( i < mdim1 ) {
82 rowptr[l] = (i+1) + dim1 * j + dim1 * dim2 * k + 1;
83 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
84 l++;
85
86 for( d=0; d<(dof*dof); d++ ) {
87 values[v] = - beta;
88 v++;
89 }
90 }
91
92 if ( j < mdim2 ) {
93 rowptr[l] = i + dim1 * (j+1) + dim1 * dim2 * k + 1;
94 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
95 l++;
96
97 for( d=0; d<(dof*dof); d++ ) {
98 values[v] = - beta;
99 v++;
100 }
101 }
102
103 if ( k < mdim3 ) {
104 rowptr[l] = i + dim1 * j + dim1 * dim2 * (k+1) + 1;
105 colptr[l] = i + dim1 * j + dim1 * dim2 * k + 1;
106 l++;
107
108 for( d=0; d<(dof*dof); d++ ) {
109 values[v] = - beta;
110 v++;
111 }
112 }
113 }
114 }
115 }
116
117 assert( l == ((2*(dim1)-1) * dim2 * dim3 + (dim2-1)*dim1*dim3 + dim2*dim1*(dim3-1)) );
118}
119
120/**
121 * @brief Create a laplacian manually
122 */
123void
124spm_example_create_laplacian( spmatrix_t *spm )
125{
126 spm_int_t dim1, dim2, dim3, dof;
127 spm_int_t *colptr;
128 spm_int_t *rowptr;
129 double *values;
130 spm_int_t n, nnz;
131
132 dim1 = 6;
133 dim2 = 6;
134 dim3 = 6;
135 dof = 3;
136 n = dim1 * dim2 * dim3;
137 nnz = dim1 * dim2 * dim3; /* Diagonal elements */
138 nnz += (dim1 - 1) * dim2 * dim3; /* Connexions along first axe */
139 nnz += dim1 * (dim2 - 1) * dim3; /* Connexions along second axe */
140 nnz += dim1 * dim2 * (dim3 - 1); /* Connexions along third axe */
141
142 /*
143 * Allocate the pointers of the spm.
144 */
145 colptr = malloc( nnz * sizeof(spm_int_t) );
146 rowptr = malloc( nnz * sizeof(spm_int_t) );
147 values = malloc( nnz * dof * dof * sizeof(double) );
148
149 /*
150 * Generate the user matrix.
151 */
152 spm_example_fill_ijv_spm( colptr, rowptr, values,
153 dim1, dim2, dim3, dof );
154
155 /*
156 * Create the associated spm data structure.
157 */
158 spmInit( spm );
159
160 /*
161 * Set the fields according to what we want to generate.
162 * All non computed fields must be set.
163 */
164 spm->baseval = 1;
165 spm->mtxtype = SpmSymmetric;
166 spm->flttype = SpmDouble;
167 spm->fmttype = SpmIJV;
168 spm->n = n;
169 spm->nnz = nnz;
170 spm->layout = SpmColMajor;
171 spm->dof = dof;
172 spm->colptr = colptr;
173 spm->rowptr = rowptr;
174 spm->values = values;
175
176 /* The full matrix is replicated on all nodes */
177 spm->replicated = 1;
178
179 /*
180 * Update the computed fields.
181 */
183}
184
185int main (int argc, char **argv)
186{
187 spmatrix_t spm;
188 spm_int_t nrhs = 10;
189 double alpha = 2.;
190 spm_int_t size, ldb, ldx;
191 double epsilon, norm;
192 void *x, *b;
193 int rc = 0;
194
195#if defined(SPM_WITH_MPI)
196 MPI_Init( &argc, &argv );
197#endif
198
199 /*
200 * Create a user made sparse matrix.
201 */
202 spm_example_create_laplacian( &spm );
203
204 /*
205 * Print basic information about the sparse matrix
206 */
207 spmPrintInfo( &spm, stdout );
208
209 /*
210 * Possibly convert the SPM into another format (CSC may be needed for a few routines)
211 */
212 spmConvert( SpmCSC, &spm );
213
214 /*
215 * Compute the frobenius norm of the spm
216 */
217 norm = spmNorm( SpmFrobeniusNorm, &spm );
218 fprintf( stdout, " || A ||_f: %e\n", norm );
219
220 /*
221 * Scale the sparse matrix.
222 */
223 if ( norm > 0. ) {
224 spmScal( 1. / norm, &spm );
225 }
226
227 /*
228 * Create a random matrix x to test products (multiple right hand side).
229 * Note that you can get the size to allocate with spm_size_of() that
230 * returns the size of each value.
231 */
232 ldb = spm.nexp > 1 ? spm.nexp : 1;
233 ldx = spm.nexp > 1 ? spm.nexp : 1;
234 size = spm_size_of( spm.flttype ) * spm.nexp * nrhs;
235 x = malloc( size );
236 rc = spmGenMat( SpmRhsRndX, nrhs, &spm, &alpha, 24356, x, ldx );
237 if ( rc != SPM_SUCCESS ) {
238 free( x );
239 spmExit( &spm );
240 return rc;
241 }
242
243 /*
244 * Compute b = A * x
245 */
246 b = malloc( size );
247 spmMatMat( SpmNoTrans, nrhs, 1., &spm, x, ldx, 0., b, ldb );
248
249 /*
250 * The following routines helps to check results obtained out of solvers by
251 * computing b = b - A * x with a given precision.
252 */
253 if ( (spm.flttype == SpmDouble ) ||
254 (spm.flttype == SpmComplex64) )
255 {
256 epsilon = 1e-15;
257 }
258 else {
259 epsilon = 1e-7;
260 }
261 epsilon = epsilon * spm.nnzexp / spm.gN;
262 rc = spmCheckAxb( epsilon, nrhs, &spm, NULL, 1, b, ldb, x, ldx );
263
264 free( x );
265 free( b );
266
267 /* Warning: Free the pointers allocated by the user */
268 spmExit( &spm );
269
270#if defined(SPM_WITH_MPI)
271 MPI_Finalize();
272#endif
273
274 (void)argc;
275 (void)argv;
276 return rc;
277}
278/**
279 * @endcode
280 */
static size_t spm_size_of(spm_coeftype_t type)
Double datatype that is not converted through precision generator functions.
Definition datatypes.h:209
@ SPM_SUCCESS
Definition const.h:82
@ SpmNoTrans
Definition const.h:155
@ SpmColMajor
Definition const.h:148
@ SpmComplex64
Definition const.h:66
@ SpmDouble
Definition const.h:64
@ SpmSymmetric
Definition const.h:166
@ SpmCSC
Definition const.h:73
@ SpmIJV
Definition const.h:75
@ SpmFrobeniusNorm
Definition const.h:200
spm_coeftype_t flttype
Definition spm.h:67
spm_int_t nexp
Definition spm.h:78
void * values
Definition spm.h:92
spm_int_t nnz
Definition spm.h:75
spm_int_t nnzexp
Definition spm.h:80
spm_int_t * rowptr
Definition spm.h:90
spm_fmttype_t fmttype
Definition spm.h:69
spm_int_t dof
Definition spm.h:82
spm_int_t n
Definition spm.h:73
spm_mtxtype_t mtxtype
Definition spm.h:65
spm_layout_t layout
Definition spm.h:87
spm_int_t baseval
Definition spm.h:71
spm_int_t * colptr
Definition spm.h:89
spm_int_t gN
Definition spm.h:72
int replicated
Definition spm.h:98
void spmExit(spmatrix_t *spm)
Cleanup the spm structure but do not free the spm pointer.
Definition spm.c:224
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.
Definition datatypes.h:70
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.
Definition spm.c:1423
int spmConvert(int ofmttype, spmatrix_t *ospm)
Convert the storage format of the spm.
Definition spm.c:418
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.
Definition spm.c:1327
void spmScal(double alpha, spmatrix_t *spm)
Scale the spm.
Definition spm.c:1481
double spmNorm(spm_normtype_t ntype, const spmatrix_t *spm)
Compute the norm of the spm.
Definition spm.c:514
void spmInit(spmatrix_t *spm)
Init the spm structure.
Definition spm.c:152
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.
Definition spm.c:1732
void spmPrintInfo(const spmatrix_t *spm, FILE *f)
Print basic informations about the spm matrix into a given stream.
Definition spm.c:1026
The sparse matrix data structure.
Definition spm.h:64