xref: /petsc/src/mat/tests/ex21.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
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";
4 
5 #include <petscmat.h>
6 
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   const PetscInt    *idx;
13   PetscScalar       v;
14   const PetscScalar *values;
15 
16   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
17   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
18   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
19   n    = 2*size;
20 
21   /* create the matrix for the five point stencil, YET AGAIN*/
22   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
23   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
24   PetscCall(MatSetFromOptions(C));
25   PetscCall(MatMPIAIJSetPreallocation(C,5,NULL,5,NULL));
26   PetscCall(MatSeqAIJSetPreallocation(C,5,NULL));
27   for (i=0; i<m; i++) {
28     for (j=2*rank; j<2*rank+2; j++) {
29       v = -1.0;  Ii = j + n*i;
30       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
31       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
32       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
33       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
34       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
35     }
36   }
37   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
38   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
39   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
40   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
41   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
42   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
43 
44   PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
45   for (i=rstart; i<rend; i++) {
46     PetscCall(MatGetRow(C,i,&nz,&idx,&values));
47     PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %" PetscInt_FMT ": ",rank,i));
48     for (j=0; j<nz; j++) {
49       PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%" PetscInt_FMT " %g  ",idx[j],(double)PetscRealPart(values[j])));
50     }
51     PetscCall(PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"\n"));
52     PetscCall(MatRestoreRow(C,i,&nz,&idx,&values));
53   }
54   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
55 
56   PetscCall(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A));
57   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO));
58   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
59   PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
60   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
61 
62   PetscCall(MatDestroy(&A));
63   PetscCall(MatDestroy(&C));
64   PetscCall(PetscFinalize());
65   return 0;
66 }
67 
68 /*TEST
69 
70    test:
71       args: -mat_type seqaij
72 
73 TEST*/
74