Actual source code: ex28.c

  1: /*$Id: ex28.c,v 1.20 2001/03/23 23:22:29 balay Exp $*/

  3: static char help[] = "Tests MatReorderForNonzeroDiagonal()nn";

  5: #include "petscmat.h"

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

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

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

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

 36:   MatGetOrdering(A,MATORDERING_NATURAL,&rowperm,&colperm);
 37:   MatReorderForNonzeroDiagonal(A,1.e-12,rowperm,colperm);
 38:   PetscPrintf(PETSC_COMM_SELF,"column and row permsn");
 39:   ISView(rowperm,0);
 40:   ISView(colperm,0);
 41:   MatLUFactorSymbolic(A,rowperm,colperm,PETSC_NULL,&LU);
 42:   MatLUFactorNumeric(A,&LU);
 43:   MatView(LU,PETSC_VIEWER_STDOUT_SELF);
 44:   VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,4,&x);
 45:   VecSetFromOptions(x);
 46:   VecDuplicate(x,&y);
 47:   values[0]=0;values[1]=1.0;values[2]=-1.0;values[3]=1.0;
 48:   for (i=0; i<4; i++) col[i]=i;
 49:   VecSetValues(x,4,col,values,INSERT_VALUES);
 50:   VecAssemblyBegin(x);
 51:   VecAssemblyEnd(x);
 52:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

 54:   MatSolve(LU,x,y);
 55:   VecView(y,PETSC_VIEWER_STDOUT_SELF);

 57:   ISDestroy(rowperm);
 58:   ISDestroy(colperm);
 59:   MatDestroy(LU);
 60:   MatDestroy(A);
 61:   VecDestroy(x);
 62:   VecDestroy(y);
 63:   PetscFinalize();
 64:   return 0;
 65: }