Actual source code: ex28.c

petsc-3.4.5 2014-06-29
  2: static char help[] = "Tests MatReorderForNonzeroDiagonal()\n\n";

  4: #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat            A,LU;
 11:   Vec            x,y;
 12:   PetscInt       nnz[4]={2,1,1,1},col[4],i;
 14:   PetscScalar    values[4];
 15:   IS             rowperm,colperm;

 17:   PetscInitialize(&argc,&args,(char*)0,help);

 19:   MatCreateSeqAIJ(PETSC_COMM_WORLD,4,4,2,nnz,&A);

 21:   /* build test matrix */
 22:   values[0]=1.0;values[1]=-1.0;
 23:   col[0]   =0;col[1]=2; i=0;
 24:   MatSetValues(A,1,&i,2,col,values,INSERT_VALUES);
 25:   values[0]=1.0;
 26:   col[0]   =1;i=1;
 27:   MatSetValues(A,1,&i,1,col,values,INSERT_VALUES);
 28:   values[0]=-1.0;
 29:   col[0]   =3;i=2;
 30:   MatSetValues(A,1,&i,1,col,values,INSERT_VALUES);
 31:   values[0]=1.0;
 32:   col[0]   =2;i=3;

 34:   MatSetValues(A,1,&i,1,col,values,INSERT_VALUES);
 35:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 36:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 37:   MatView(A,PETSC_VIEWER_STDOUT_SELF);

 39:   MatGetOrdering(A,MATORDERINGNATURAL,&rowperm,&colperm);
 40:   MatReorderForNonzeroDiagonal(A,1.e-12,rowperm,colperm);
 41:   PetscPrintf(PETSC_COMM_SELF,"column and row perms\n");
 42:   ISView(rowperm,0);
 43:   ISView(colperm,0);
 44:   MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&LU);
 45:   MatLUFactorSymbolic(LU,A,rowperm,colperm,NULL);
 46:   MatLUFactorNumeric(LU,A,NULL);
 47:   MatView(LU,PETSC_VIEWER_STDOUT_SELF);
 48:   VecCreate(PETSC_COMM_WORLD,&x);
 49:   VecSetSizes(x,PETSC_DECIDE,4);
 50:   VecSetFromOptions(x);
 51:   VecDuplicate(x,&y);

 53:   values[0]=0;values[1]=1.0;values[2]=-1.0;values[3]=1.0;
 54:   for (i=0; i<4; i++) col[i]=i;
 55:   VecSetValues(x,4,col,values,INSERT_VALUES);
 56:   VecAssemblyBegin(x);
 57:   VecAssemblyEnd(x);
 58:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

 60:   MatSolve(LU,x,y);
 61:   VecView(y,PETSC_VIEWER_STDOUT_SELF);

 63:   ISDestroy(&rowperm);
 64:   ISDestroy(&colperm);
 65:   MatDestroy(&LU);
 66:   MatDestroy(&A);
 67:   VecDestroy(&x);
 68:   VecDestroy(&y);
 69:   PetscFinalize();
 70:   return 0;
 71: }