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