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 PetscFunctionBeginUser; 16 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 17 /* Create the matrix for the five point stencil, YET AGAIN */ 18 PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 19 PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 20 PetscCall(MatSetFromOptions(C)); 21 PetscCall(MatSetUp(C)); 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; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 26 if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 27 if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 28 if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 29 v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 30 } 31 } 32 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 33 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 34 PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 35 36 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 37 PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); 38 for (i=rstart; i<rend; i++) { 39 PetscCall(MatGetRow(C,i,&nz,&idx,&values)); 40 if (rank == 0) { 41 for (j=0; j<nz; j++) { 42 PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %g ",idx[j],(double)PetscRealPart(values[j]))); 43 } 44 PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n")); 45 } 46 PetscCall(MatRestoreRow(C,i,&nz,&idx,&values)); 47 } 48 49 PetscCall(MatDestroy(&C)); 50 PetscCall(PetscFinalize()); 51 return 0; 52 } 53 54 /*TEST 55 56 test: 57 58 TEST*/ 59