Actual source code: ex21.c
petsc-3.10.5 2019-03-28
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>
7: int main(int argc,char **args)
8: {
9: Mat C,A;
10: PetscInt i,j,m = 3,n = 2,Ii,J,rstart,rend,nz;
11: PetscMPIInt rank,size;
12: PetscErrorCode ierr;
13: const PetscInt *idx;
14: PetscScalar v;
15: const PetscScalar *values;
17: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
19: MPI_Comm_size(PETSC_COMM_WORLD,&size);
20: n = 2*size;
22: /* create the matrix for the five point stencil, YET AGAIN*/
23: MatCreate(PETSC_COMM_WORLD,&C);
24: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);
25: MatSetFromOptions(C);
26: MatMPIAIJSetPreallocation(C,5,NULL,5,NULL);
27: MatSeqAIJSetPreallocation(C,5,NULL);
28: for (i=0; i<m; i++) {
29: for (j=2*rank; j<2*rank+2; j++) {
30: v = -1.0; Ii = j + n*i;
31: if (i>0) {J = Ii - n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
32: if (i<m-1) {J = Ii + n; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
33: if (j>0) {J = Ii - 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
34: if (j<n-1) {J = Ii + 1; MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);}
35: v = 4.0; MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);
36: }
37: }
38: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
39: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
40: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
41: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
42: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
43: MatView(C,PETSC_VIEWER_STDOUT_WORLD);
45: MatGetOwnershipRange(C,&rstart,&rend);
46: for (i=rstart; i<rend; i++) {
47: MatGetRow(C,i,&nz,&idx,&values);
48: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %D: ",rank,i);
49: for (j=0; j<nz; j++) {
50: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%D %g ",idx[j],(double)PetscRealPart(values[j]));
51: }
52: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n");
53: MatRestoreRow(C,i,&nz,&idx,&values);
54: }
55: PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);
57: MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A);
58: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO);
59: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
60: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
61: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
63: MatDestroy(&A);
64: MatDestroy(&C);
65: PetscFinalize();
66: return ierr;
67: }
70: /*TEST
72: test:
73: args: -mat_type seqaij
75: TEST*/