Actual source code: ex7.c

  1: /*$Id: ex7.c,v 1.20 2001/08/07 21:30:08 bsmith Exp $*/

  3: static char help[] = "Tests matrix factorization.  Note that most users shouldn
  4: employ the SLES interface to the linear solvers instead of using the factorizationn
  5: routines directly.nn";

 7:  #include petscmat.h

  9: int main(int argc,char **args)
 10: {
 11:   Mat         C,LU;
 12:   MatInfo     info;
 13:   int         i,j,m = 3,n = 3,I,J,ierr;
 14:   PetscScalar v,mone = -1.0,one = 1.0;
 15:   IS          perm,iperm;
 16:   Vec         x,u,b;
 17:   PetscReal   norm;

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

 21:   /* Create the matrix for the five point stencil, YET AGAIN */
 22:   MatCreate(PETSC_COMM_SELF,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n,&C);
 23:   MatSetFromOptions(C);
 24:   for (i=0; i<m; i++) {
 25:     for (j=0; j<n; j++) {
 26:       v = -1.0;  I = j + n*i;
 27:       if (i>0)   {J = I - n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 28:       if (i<m-1) {J = I + n; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 29:       if (j>0)   {J = I - 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 30:       if (j<n-1) {J = I + 1; MatSetValues(C,1,&I,1,&J,&v,INSERT_VALUES);}
 31:       v = 4.0; MatSetValues(C,1,&I,1,&I,&v,INSERT_VALUES);
 32:     }
 33:   }
 34:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 35:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
 36:   MatGetOrdering(C,MATORDERING_RCM,&perm,&iperm);
 37:   MatView(C,PETSC_VIEWER_STDOUT_WORLD);
 38:   ISView(perm,PETSC_VIEWER_STDOUT_SELF);

 40:   MatLUFactorSymbolic(C,perm,iperm,PETSC_NULL,&LU);
 41:   MatLUFactorNumeric(C,&LU);
 42:   MatView(LU,PETSC_VIEWER_STDOUT_WORLD);

 44:   VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
 45:   VecSet(&one,u);
 46:   VecDuplicate(u,&x);
 47:   VecDuplicate(u,&b);

 49:   MatMult(C,u,b);
 50:   MatSolve(LU,b,x);

 52:   VecView(b,PETSC_VIEWER_STDOUT_SELF);
 53:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

 55:   VecAXPY(&mone,u,x);
 56:   VecNorm(x,NORM_2,&norm);
 57:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %An",norm);

 59:   MatGetInfo(C,MAT_LOCAL,&info);
 60:   PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %dn",(int)info.nz_used);
 61:   MatGetInfo(LU,MAT_LOCAL,&info);
 62:   PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %dn",(int)info.nz_used);

 64:   VecDestroy(u);
 65:   VecDestroy(b);
 66:   VecDestroy(x);
 67:   ISDestroy(perm);
 68:   ISDestroy(iperm);
 69:   MatDestroy(C);
 70:   MatDestroy(LU);

 72:   PetscFinalize();
 73:   return 0;
 74: }