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