1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C,A; 9c4762a1bSJed Brown PetscInt i,j,m = 5,n = 5,Ii,J; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown IS perm,iperm; 13c4762a1bSJed Brown 14c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 15*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); 16c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 17c4762a1bSJed Brown for (i=0; i<m; i++) { 18c4762a1bSJed Brown for (j=0; j<n; j++) { 19c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 20*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 21*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 22*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 23*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 24*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 25c4762a1bSJed Brown } 26c4762a1bSJed Brown } 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 29c4762a1bSJed Brown 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATSAME,MAT_INITIAL_MATRIX,&A)); 31c4762a1bSJed Brown 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGND,&perm,&iperm)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(perm,PETSC_VIEWER_STDOUT_SELF)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(iperm,PETSC_VIEWER_STDOUT_SELF)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 36c4762a1bSJed Brown 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 41c4762a1bSJed Brown ierr = PetscFinalize(); 42c4762a1bSJed Brown return ierr; 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45c4762a1bSJed Brown /*TEST 46c4762a1bSJed Brown 47c4762a1bSJed Brown test: 48c4762a1bSJed Brown filter: grep -v "MPI processes" 49c4762a1bSJed Brown 50c4762a1bSJed Brown TEST*/ 51