xref: /petsc/src/mat/tests/ex14.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
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   PetscErrorCode    ierr;
11   PetscMPIInt       rank;
12   const PetscInt    *idx;
13   PetscScalar       v;
14   const PetscScalar *values;
15 
16   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
17   /* Create the matrix for the five point stencil, YET AGAIN */
18   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
19   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
20   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
21   ierr = MatSetUp(C);CHKERRQ(ierr);
22   for (i=0; i<m; i++) {
23     for (j=0; j<n; j++) {
24       v = -1.0;  Ii = j + n*i;
25       if (i>0)   {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
26       if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
27       if (j>0)   {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
28       if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);}
29       v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr);
30     }
31   }
32   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
33   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
35 
36   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
37   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
38   for (i=rstart; i<rend; i++) {
39     ierr = MatGetRow(C,i,&nz,&idx,&values);CHKERRQ(ierr);
40     if (!rank) {
41       for (j=0; j<nz; j++) {
42         ierr = PetscPrintf(PETSC_COMM_SELF,"%D %g ",idx[j],(double)PetscRealPart(values[j]));CHKERRQ(ierr);
43       }
44       ierr = PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr);
45     }
46     ierr = MatRestoreRow(C,i,&nz,&idx,&values);CHKERRQ(ierr);
47   }
48 
49   ierr = MatDestroy(&C);CHKERRQ(ierr);
50   ierr = PetscFinalize();
51   return ierr;
52 }
53 
54 
55 /*TEST
56 
57    test:
58 
59 TEST*/
60