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