xref: /petsc/src/mat/tests/ex14.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests sequential and parallel MatGetRow() and MatRestoreRow().\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat               C;
9   PetscInt          i,j,m = 5,n = 5,Ii,J,nz,rstart,rend;
10   PetscMPIInt       rank;
11   const PetscInt    *idx;
12   PetscScalar       v;
13   const PetscScalar *values;
14 
15   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
16   /* Create the matrix for the five point stencil, YET AGAIN */
17   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
18   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
19   PetscCall(MatSetFromOptions(C));
20   PetscCall(MatSetUp(C));
21   for (i=0; i<m; i++) {
22     for (j=0; j<n; j++) {
23       v = -1.0;  Ii = j + n*i;
24       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
25       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
26       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
27       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
28       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
29     }
30   }
31   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
32   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
33   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
34 
35   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
36   PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
37   for (i=rstart; i<rend; i++) {
38     PetscCall(MatGetRow(C,i,&nz,&idx,&values));
39     if (rank == 0) {
40       for (j=0; j<nz; j++) {
41         PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %g ",idx[j],(double)PetscRealPart(values[j])));
42       }
43       PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n"));
44     }
45     PetscCall(MatRestoreRow(C,i,&nz,&idx,&values));
46   }
47 
48   PetscCall(MatDestroy(&C));
49   PetscCall(PetscFinalize());
50   return 0;
51 }
52 
53 /*TEST
54 
55    test:
56 
57 TEST*/
58