Actual source code: samgpetsctools.c
1: #include "global.h"
2: #include petscksp.h
3: #include samgfunc.h
4: #include petscfunc.h
5: #include externc.h
8: void USER_coo(int * i,int * ndim, double * x, double * y, double * z)
9: {
10: printf("in user_coo");
11: }
14: static double Machine_Precision_Eps = 2.e-16;
15: /* ------------------------------------------------------------------- */
18: /*..SamgGetGrid - This routine gets an array of grids
19: INPUT: levels: number of levels created by SAMG
20: numnodes: number of nodes on finest grid
21: numnonzeros: number of nonzeros on coarsest grid
22: OUTPUT: grid : array of grids ..*/
23: PetscErrorCode SamgGetGrid(int levels, int numnodes, int numnonzero,
24: GridCtx* grid, void* ctx)
25: {
26: int k;
27: int ia_shift[MAX_LEVELS], ja_shift[MAX_LEVELS], nnu_cg, nna_cg;
28: int iw_shift[MAX_LEVELS], jw_shift[MAX_LEVELS], rows_weights,
29: nna_weights, dummy;
31: MatInfo info;
33: /*..Get coarse grid operators..*/
34: /*....Initialize ia_shift, ja_shift, nnu_cg and nna_cg....*/
35: ia_shift[1] = 1;
36: ja_shift[1] = 1;
37: nnu_cg = numnodes;
38: nna_cg = numnonzero;
40: for (k=2;k<=levels;k++){ /*....We do not get the finest level matrix....*/
41: /*....Update ia_shift and ja_shift values with nna_cg and nnu_cg
42: from previous loop....*/
43: ia_shift[k] = ia_shift[k-1] + nna_cg ;
44: ja_shift[k] = ja_shift[k-1] + nnu_cg ;
46: /*....Get coarse grid matrix on level k....*/
47: SamgGetCoarseMat(k, ia_shift[k], ja_shift[k], &(grid[k].A),
48: PETSC_NULL);
50: /*....Get size and number of nonzeros of coarse grid matrix on
51: level k, i.e. get new nna_cg and nnu_cg values....*/
52: MatGetSize(grid[k].A, &nnu_cg, &nnu_cg);
53: MatGetInfo(grid[k].A, MAT_LOCAL, &info);
54: nna_cg = int(info.nz_used);
55: }
56:
57: /*..Get interpolation operators..*/
58: /*....Initialize iw_shift, jw_shift and nna_weights....*/
59: iw_shift[0] = 1;
60: jw_shift[0] = 1;
61: nna_weights = 0;
62: rows_weights = numnodes;
64: for (k=1;k<=levels-1;k++){/*....There's NO interpolation operator
65: associated to the coarsest level....*/
66: /*....Update iw_shift with nna_weights value from
67: previous loop....*/
68: iw_shift[k] = iw_shift[k-1] + nna_weights ;
69: /*....Update jw_shift with rows_weights value from
70: current loop....*/
71: jw_shift[k] = jw_shift[k-1] + rows_weights ;
72:
73: /*....Get interpolation from level k+1 to level k....*/
74: SamgGetInterpolation(k, iw_shift[k], jw_shift[k],
75: &(grid[k].Interp), PETSC_NULL) ;
77: /*....Get number of collumns and number of nonzeros of
78: interpolation associated to level k. NOTE: The
79: number of collums at this loop equals the number of
80: rows at the next loop...*/
81: MatGetSize(grid[k].Interp, &dummy, &rows_weights);
82: MatGetInfo(grid[k].Interp, MAT_LOCAL, &info);
83: nna_weights = int(info.nz_used);
84: }
86: return 0;
87: }
88: /* ------------------------------------------------------------------- */
91: /*..SamgGetCoarseMat - This routine gets the coarse level matrix on the
92: level specified at input.
93: WARNING: This routine does not work to get the fine level matrix,
94: i.e. the value of k at input should be at least 2
95: INPUT: level: current grid level
96: ia_shift: shift to apply on ia_cg elements
97: ja_shift: shift to apply on ja_cg elements
98: OUTPUT: coarsemat: coarse level matrix
99: ..*/
101: PetscErrorCode SamgGetCoarseMat(int level, int ia_shift, int ja_shift,
102: Mat* coarsemat, void* ctx)
103: {
104: int nnu_k, nna_k; /* size and non-zeros of operator on level k */
105: int *ia_k, *ja_k; /* coarse grid matrix in skyline format */
106: double *a_k;
107: int *nnz_per_row; /* integer vector to hold the number of nonzeros */
108: /* of each row. This vector will be used to */
109: /* allocate memory for the matrix, and to store */
110: /* elements in the matrix */
112: int I;
114: /*..Get size (nnu_k) and number of non-zeros (nna_k) of operator
115: on level k..*/
116: SAMGPETSC_get_dim_operator(&level, &nnu_k, &nna_k);
118: /*..Now that nnu_cg and nna_cg are known, we can allocate memory for
119: coarse level matrix in compresses skyline format..*/
120: PetscMalloc(nna_k * sizeof(double),&a_k);
121: PetscMalloc((nnu_k+1) * sizeof(int),&ia_k);
122: PetscMalloc(nna_k * sizeof(int),&ja_k);
124: /*..Get coarse grid matrix in skyline format..*/
125: SAMGPETSC_get_operator(&level, a_k, ia_k, ja_k);
127: /*..Apply shift on each of the ia_cg and ja_cg elements..*/
128: SAMGPETSC_apply_shift(ia_k, &nnu_k, &ia_shift,
129: ja_k, &nna_k, &ja_shift);
131: PetscMalloc(nnu_k * sizeof(int),&nnz_per_row);
133: /*..The numbero f nonzeros entries in row I can be calculated as
134: ia[I+1] - 1 - ia[I] + 1 = ia[I+1] - ia[I] ..*/
135: for (I=0;I<nnu_k;I++)
136: nnz_per_row[I] = ia_k[I+1] - ia_k[I];
138: /*..Allocate (create) SeqAIJ matrix for use within PETSc..*/
139: MatCreate(PETSC_COMM_WORLD,nnu_k,nnu_k,nnu_k,nnu_k,coarsemat);
140: MatSetType(*coarsemat,MATSEQAIJ);
141: MatSeqAIJSetPreallocation(*coarsemat,0,nnz_per_row);
143: /*..Store coarse grid matrix in Petsc Mat object..*/
144: for (I=0;I<nnu_k;I++){
145: MatSetValues(*coarsemat,
146: 1, /* number of rows */
147: &I, /* pointer to global row number */
148: nnz_per_row[I], /* number of collums = number of nonzero ... */
149: /* entries in row I */
150: &(ja_k[ ia_k[I] ]),
151: /* vector global column indices */
152: (PetscScalar *) &(a_k[ ia_k[I] ]),
153: /* vector of coefficients */
154: INSERT_VALUES);
155: }
157: MatAssemblyBegin(*coarsemat,MAT_FINAL_ASSEMBLY);
158: MatAssemblyEnd(*coarsemat,MAT_FINAL_ASSEMBLY);
160: /*..Free memory required by storage in skyline format..*/
161: PetscFree(a_k);
162: PetscFree(ia_k);
163: PetscFree(ja_k);
164: /*..Free memory for auxilary array..*/
165: PetscFree(nnz_per_row);
166:
167: return 0;
168: }
169: /* ------------------------------------------------------------------- */
172: /*..SamgGetInterpolation - Get interpolation operator that interpolates
173: from level k+1 (coarse grid) to k (fine grid), where the input level
174: equals k ..*/
175: /*..WARNING: This routine assumes that the input value for level is strictly
176: smaller than the number of levels created..*/
177: /*..Implementation notes
178: o) The interpolation is a rectangular matrix with
179: number of rows equal to fine grid dim
180: cols coarse.
181: ..*/
182: PetscErrorCode SamgGetInterpolation(int level, int iw_shift, int jw_shift,
183: Mat* interpolation, void* ctx)
184: {
185: int rows_weights, cols_weights, nna_weights;
186: int *iweights, *jweights;
187: double *weights;
188: int *nnz_per_row; /* integer vector to hold the number of nonzeros */
189: /* of each row. This vector will be used to */
190: /* allocate memory for the matrix, and to store */
191: /* elements in the matrix */
193: int I, coarser_level=level+1, dummy;
195: /*..Get number of rows and number of nonzeros of interpolation operator..*/
196: SAMGPETSC_get_dim_interpol(&level, &rows_weights, &nna_weights);
198: /*..Get number of cols of interpolation operator. NOTE: The number of
199: collums of the interpolation on level k equals the size of
200: the coarse grid matrix on the next coarsest grid.
201: SAMGPETSC_get_dim_interpol does not allow to get the number of
202: collumns of next to coarsest grid..*/
203: SAMGPETSC_get_dim_operator(&coarser_level, &cols_weights, &dummy);
205: /*..Now that nnu_weights and nna_weights are known, we can allocate
206: memory for interpolation operator in compresses skyline format..*/
207: PetscMalloc(nna_weights * sizeof(double),&weights);
208:
209: PetscMalloc((rows_weights+1) * sizeof(int),&iweights);
210:
211: PetscMalloc(nna_weights * sizeof(int),&jweights);
212:
214: /*..Get interpolation operator in compressed skyline format..*/
215: SAMGPETSC_get_interpol(&level, weights, iweights, jweights);
217: /*..Apply shift on each of the ia_cg and ja_cg elements..*/
218: SAMGPETSC_apply_shift(iweights, &rows_weights, &iw_shift,
219: jweights, &nna_weights, &jw_shift);
221: PetscMalloc(rows_weights * sizeof(int),&nnz_per_row);
222:
224: /*..The numbero f nonzeros entries in row I can be calculated as
225: ia[I+1] - 1 - ia[I] + 1 = ia[I+1] - ia[I] ..*/
226: for (I=0;I<rows_weights;I++)
227: nnz_per_row[I] = iweights[I+1] - iweights[I];
229: /*..Allocate (create) SeqAIJ matrix for use within PETSc..*/
230: MatCreate(PETSC_COMM_WORLD,rows_weights,cols_weights,rows_weights,cols_weights,interpolation);
231: MatSetType(*interpolation,MATSEQAIJ);
232: MatSeqAIJSetPreallocation(*interpolation,0,nnz_per_row);
234: /*..Store coarse grid matrix in Petsc Mat object..*/
235: for (I=0;I<rows_weights;I++){
236: MatSetValues(*interpolation,
237: 1, /* number of rows */
238: &I, /* pointer to global row number */
239: nnz_per_row[I], /* number of collums = number of nonzero ... */
240: /* entries in row I */
241: &(jweights[ iweights[I] ]),
242: /* vector global column indices */
243: (PetscScalar *) &(weights[ iweights[I] ]),
244: /* vector of coefficients */
245: INSERT_VALUES);
246: }
248: MatAssemblyBegin(*interpolation,MAT_FINAL_ASSEMBLY);
249: MatAssemblyEnd(*interpolation,MAT_FINAL_ASSEMBLY);
251: /*..Free memory required by storage in skyline format..*/
252: PetscFree(weights);
253: PetscFree(iweights);
254: PetscFree(jweights);
255: /*..Free memory for auxilary array..*/
256: PetscFree(nnz_per_row);
258: return 0;
259: }
260: /* ------------------------------------------------------------------- */
263: /*..Write coarser grid operators constructed by SAMG to ASCII file..*/
265: PetscErrorCode SamgPetscWriteOperator(const int numnodes, const double* Asky,
266: const int* ia, const int* ja, int extension)
267: {
269: static char filename[80];
270: int I,j,j1,j2;
271: FILE *output;
272:
273: /*..Switch arrays iacopy and jacopy to C conventions..*/
274: // for (j=0;j<=numnodes;j++)
275: // iacopy[j]--;
276: // for (j=0;j<ia[numnodes];j++)
277: // jacopy[j]--;
278:
279: /*....Write matrix to file....*/
280: sprintf(filename,"coarsemat.%02u", extension);
281: output=fopen(filename,"w");
282: fprintf(output, "%% \n");
283: fprintf(output, "%% %d %d \n", numnodes, ia[numnodes] );
285: for (I=0;I<numnodes;I++){
286: j1 = ia[I];
287: j2 = ia[I+1] - 1;
288: for (j=j1;j<=j2;j++){
289: fprintf(output, "%d %d %22.18e\n", I+1, ja[j]+1,
290: Asky[j] );
291: // printf("%d %d %e \n", I+1, ja[j]+1,
292: // Asky[j] );
293: }
294: }
295: fclose(output);
296:
297: return 0;
298: }
299: /* ------------------------------------------------------------------- */
302: /*..Write interpolation operators constructed by SAMG to ASCII file..*/
304: PetscErrorCode SamgPetscWriteInterpol(const int numrows, const double* weights,
305: const int* iweights, const int* jweights, int extension)
306: {
308: static char filename[80];
309: int I,j,j1,j2,numcols,numnonzero;
310: FILE *output;
311:
312: /*..Set number of nonzeros..*/
313: numnonzero = iweights[numrows];
315: /*..Determine numcols as the maximum ja value +1..*/
316: numcols = jweights[0];
317: for (j=0;j<numnonzero;j++){
318: if (jweights[j] > numcols) numcols = jweights[j];
319: }
320: numcols++;
322: /*..Write interpolation operator from grid k+1 (coarse grid) grid to k
323: (finer grid) to file..*/
324: sprintf(filename,"interpol.%02u%02u",
325: extension+1, extension);
326: output=fopen(filename,"w");
327: fprintf(output, "%% \n%% %d %d %d \n", numrows, numcols, iweights[numrows] );
328: for (I=0;I<numrows;I++){
329: j1 = iweights[I];
330: j2 = iweights[I+1] - 1;
331: for (j=j1;j<=j2;j++){
332: fprintf(output, "%d %d %22.18e\n", I+1, jweights[j]+1,
333: weights[j] );
334: // printf("%d %d %e \n", I+1, jweights[j]+1,
335: // weights[j] );
336: }
337: }
338: fclose(output);
340: return 0;
341: }
342: /* ------------------------------------------------------------------- */
345: /*..SamgCheckGalerkin - This routine offers a check on the correctness
346: of how SAMG interpolation and coarse grid operators are parsed to
347: PETSc. This routine computes I^H_h A^h I^h_H by PETSc matrix - matrix
348: multiplications, and compares this product with A^H..*/
349:
350: PetscErrorCode SamgCheckGalerkin(int levels, Mat A, GridCtx* grid,
351: void* ctx)
352: {
353: Mat FineLevelMatrix, Restriction, HalfGalerkin, Galerkin, Diff;
354: double normdiff;
356: int k;
358: for (k=1;k<=levels-1;k++){
359: if (k==1)
360: FineLevelMatrix = A;
361: else
362: FineLevelMatrix = grid[k].A;
363: /*....Compute A^h I^h_H....*/
364: MatMatMult(FineLevelMatrix, grid[k].Interp, &HalfGalerkin);
365: /*....Get I^h_H....*/
366: MatTranspose(grid[k].Interp,&Restriction);
367: /*....Compute I^H_h A^h I^h_H....*/
368: MatMatMult(Restriction, HalfGalerkin, &Galerkin);
369: /*....Compute A^H - I^H_h A^h I^h_H....*/
370: MatSubstract(grid[k+1].A, Galerkin, &Diff);
371: /*....Compute || A^H - I^H_h A^h I^h_H||_{\infty}....*/
372: MatNorm(Diff,NORM_INFINITY,&normdiff);
374: printf("SamgCheckGalerkin :: || A^H - I^H_h A^h I^h_H||_{infty} on level %8d = %e\n",
375: k+1, normdiff);
377: MatDestroy(Restriction);
378: MatDestroy(HalfGalerkin);
379: MatDestroy(Galerkin);
380: MatDestroy(Diff);
381: }
382: return 0;
383: }
384: /* ------------------------------------------------------------------- */
387: /*..MatSubstract - Computes the difference Term1 - Term2
388: INPUT: Term1, Term2 : The input matrices
389: OUTPUT: Prod: the difference
390: NOTE: Memory needed by difference has to freed outside this routine!
391: ..*/
393: PetscErrorCode MatSubstract(Mat Term1, Mat Term2, Mat* Diff)
394: {
395: Vec col_vec1, col_vec2, diff_vec;
397: int rows1, cols1, rows2, cols2, col, row;
398: static PetscScalar dminusone = -1.;
399: PetscScalar matrix_element ;
400: PetscScalar *vec_getvalues ;
401: double inf_norm_diff_vec = 0.0 ;
402: double Zero_Element_Treshold = 0.0 ;
404: /*..Get sizes of terms..*/
405: MatGetSize(Term1, &rows1, &cols1);
406: MatGetSize(Term2, &rows2, &cols2);
408: /*..Check input..*/
409: if ( (cols1 != rows1) || (cols2 != rows2) ){
410: SETERRQ(PETSC_ERR_ARG_SIZ,"Error in MatMatMult: cols1 <> rows1 or cols1 <> rows1");
411: }
413: /*..Create difference of 2 SeqAIJ matrices..*/
414: MatCreate(PETSC_COMM_WORLD,rows1,cols1,rows1,cols1,Diff);
415: MatSetType(*Diff,MATSEQAIJ);
416: MatSeqAIJSetPreallocation(*Diff,0,PETSC_NULL);
418: /*..Create vectors..*/
419: VecCreate(MPI_COMM_WORLD,&col_vec1);
420: VecSetSizes(col_vec1,PETSC_DECIDE,rows1);
421: VecSetType(col_vec1,VECSEQ);
422: VecDuplicate(col_vec1, &col_vec2);
423: VecDuplicate(col_vec1, &diff_vec);
425: for (col=0;col<cols1;col++){
426: /*..Get collumns..*/
427: MatGetColumnVector(Term1,col_vec1,col);
428: MatGetColumnVector(Term2,col_vec2,col);
429: /*..Substract collumns..*/
430: VecWAXPY(&dminusone, col_vec2, col_vec1, diff_vec);
431:
432:
433: /*..Compute norm..*/
434: VecNorm( diff_vec, NORM_INFINITY, &inf_norm_diff_vec );
435:
436: /*..Set threshold..*/
437: Zero_Element_Treshold = inf_norm_diff_vec * Machine_Precision_Eps ;
439: /*..Get Term1(:,col) - Term2(:,col) values..*/
440: VecGetArray(diff_vec, &vec_getvalues);
442: for (row=0;row<rows1;row++){
443: matrix_element = vec_getvalues[row];
444: if ( PetscAbsScalar( matrix_element ) >= Zero_Element_Treshold ) {
445: MatSetValue(*Diff, row, col,matrix_element, INSERT_VALUES );
446: }
447: }
448: VecRestoreArray(diff_vec, &vec_getvalues);
449: }
451: MatAssemblyBegin(*Diff,MAT_FINAL_ASSEMBLY);
452: MatAssemblyEnd(*Diff,MAT_FINAL_ASSEMBLY);
454: VecDestroy(col_vec1);
455: VecDestroy(col_vec2);
456: VecDestroy(diff_vec);
458: return 0;
459: }
460: /* ------------------------------------------------------------------- */