xref: /petsc/src/mat/tests/ex14.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
17c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
189566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
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;
259566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
269566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
279566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
289566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));}
299566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES));
30c4762a1bSJed Brown     }
31c4762a1bSJed Brown   }
329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
349566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
379566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C,&rstart,&rend));
38c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
399566063dSJacob Faibussowitsch     PetscCall(MatGetRow(C,i,&nz,&idx,&values));
40dd400576SPatrick Sanan     if (rank == 0) {
41c4762a1bSJed Brown       for (j=0; j<nz; j++) {
429566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %g ",idx[j],(double)PetscRealPart(values[j])));
43c4762a1bSJed Brown       }
449566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n"));
45c4762a1bSJed Brown     }
469566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(C,i,&nz,&idx,&values));
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
509566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
51b122ec5aSJacob Faibussowitsch   return 0;
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*TEST
55c4762a1bSJed Brown 
56c4762a1bSJed Brown    test:
57c4762a1bSJed Brown 
58c4762a1bSJed Brown TEST*/
59