Actual source code: ex21.c
petsc-3.5.4 2015-05-23
2: static char help[] = "Tests converting a parallel AIJ formatted matrix to the parallel Row format.\n\
3: This also tests MatGetRow() and MatRestoreRow() for the parallel case.\n\n";
5: #include <petscmat.h>
9: int main(int argc,char **args)
10: {
11: Mat C,A;
12: PetscInt i,j,m = 3,n = 2,Ii,J,rstart,rend,nz;
13: PetscMPIInt rank,size;
14: PetscErrorCode ierr;
15: const PetscInt *idx;
16: PetscScalar v;
17: const PetscScalar *values;
19: PetscInitialize(&argc,&args,(char*)0,help);
20: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: n = 2*size;
24: /* create the matrix for the five point stencil, YET AGAIN*/
25: MatCreate(PETSC_COMM_WORLD,&C);
26: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
27: MatSetFromOptions(C);
28: MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);
29: MatSeqAIJSetPreallocation(C,5,NULL);
30: for (i=0; i<m; i++) {
31: for (j=2*rank; j<2*rank+2; j++) {
32: v = -1.0; Ii = j + n*i;
33: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
34: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
35: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
36: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
37: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
38: }
39: }
40: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
41: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
42: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
43: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
44: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
45: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
47: MatGetOwnershipRange(C,&rstart,&rend);
48: for (i=rstart; i<rend; i++) {
49: MatGetRow(C,i,&nz,&idx,&values);
50: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %D: ",rank,i);
51: for (j=0; j<nz; j++) {
52: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%D %g ",idx[j],(double)PetscRealPart(values[j]));
53: }
54: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n");
55: MatRestoreRow(C,i,&nz,&idx,&values);
56: }
57: PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);
59: MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
60: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
61: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
62: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
63: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
65: MatDestroy(&A);
66: MatDestroy(&C);
67: PetscFinalize();
68: return 0;
69: }