1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests matrix ordering routines.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown extern PetscErrorCode MatGetOrdering_myordering(Mat,MatOrderingType,IS*,IS*); 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat C,Cperm; 10c4762a1bSJed Brown PetscInt i,j,m = 5,n = 5,Ii,J,ncols; 11c4762a1bSJed Brown PetscErrorCode ierr; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown PetscMPIInt size; 14c4762a1bSJed Brown IS rperm,cperm,icperm; 15c4762a1bSJed Brown const PetscInt *rperm_ptr,*cperm_ptr,*cols; 16c4762a1bSJed Brown const PetscScalar *vals; 17c4762a1bSJed Brown PetscBool TestMyorder=PETSC_FALSE; 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN */ 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&C)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 26c4762a1bSJed Brown for (i=0; i<m; i++) { 27c4762a1bSJed Brown for (j=0; j<n; j++) { 28c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 29*5f80ce2aSJacob Faibussowitsch if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 30*5f80ce2aSJacob Faibussowitsch if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 31*5f80ce2aSJacob Faibussowitsch if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 32*5f80ce2aSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 33*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown } 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 38c4762a1bSJed Brown 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(C,MATORDERINGND,&rperm,&cperm)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&rperm)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cperm)); 43c4762a1bSJed Brown 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(C,MATORDERINGRCM,&rperm,&cperm)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&rperm)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cperm)); 48c4762a1bSJed Brown 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(C,MATORDERINGQMD,&rperm,&cperm)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&rperm)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cperm)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* create Cperm = rperm*C*icperm */ 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-testmyordering",&TestMyorder,NULL)); 56c4762a1bSJed Brown if (TestMyorder) { 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering_myordering(C,MATORDERINGQMD,&rperm,&cperm)); 58c4762a1bSJed Brown printf("myordering's rperm:\n"); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISView(rperm,PETSC_VIEWER_STDOUT_SELF)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISInvertPermutation(cperm,PETSC_DECIDE,&icperm)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(rperm,&rperm_ptr)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(icperm,&cperm_ptr)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,m*n,m*n,5,NULL,&Cperm)); 64c4762a1bSJed Brown for (i=0; i<m*n; i++) { 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(C,rperm_ptr[i],&ncols,&cols,&vals)); 66c4762a1bSJed Brown for (j=0; j<ncols; j++) { 67c4762a1bSJed Brown /* printf(" (%d %d %g)\n",i,cperm_ptr[cols[j]],vals[j]); */ 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(Cperm,1,&i,1,&cperm_ptr[cols[j]],&vals[j],INSERT_VALUES)); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown } 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Cperm,MAT_FINAL_ASSEMBLY)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Cperm,MAT_FINAL_ASSEMBLY)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(rperm,&rperm_ptr)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(icperm,&cperm_ptr)); 75c4762a1bSJed Brown 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&rperm)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&cperm)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&icperm)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cperm)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 83c4762a1bSJed Brown ierr = PetscFinalize(); 84c4762a1bSJed Brown return ierr; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87c4762a1bSJed Brown #include <petsc/private/matimpl.h> 88c4762a1bSJed Brown /* This is modified from MatGetOrdering_Natural() */ 89c4762a1bSJed Brown PetscErrorCode MatGetOrdering_myordering(Mat mat,MatOrderingType type,IS *irow,IS *icol) 90c4762a1bSJed Brown { 91c4762a1bSJed Brown PetscInt n,i,*ii; 92c4762a1bSJed Brown PetscBool done; 93c4762a1bSJed Brown MPI_Comm comm; 94c4762a1bSJed Brown 95c4762a1bSJed Brown PetscFunctionBegin; 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)mat,&comm)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&n,NULL,NULL,&done)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,NULL,NULL,&done)); 99c4762a1bSJed Brown if (done) { /* matrix may be "compressed" in symbolic factorization, due to i-nodes or block storage */ 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&ii)); 101c4762a1bSJed Brown for (i=0; i<n; i++) ii[i] = n-i-1; /* replace your index here */ 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_COPY_VALUES,irow)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,ii,PETSC_OWN_POINTER,icol)); 104c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"MatRestoreRowIJ fails!"); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISSetPermutation(*irow)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISSetPermutation(*icol)); 107c4762a1bSJed Brown PetscFunctionReturn(0); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /*TEST 111c4762a1bSJed Brown 112c4762a1bSJed Brown test: 113c4762a1bSJed Brown 114c4762a1bSJed Brown TEST*/ 115