Actual source code: ex7.c

  2: static char help[] = "Tests matrix factorization.  Note that most users should\n\
  3: employ the KSP  interface to the linear solvers instead of using the factorization\n\
  4: routines directly.\n\n";

 6:  #include petscmat.h

 10: int main(int argc,char **args)
 11: {
 12:   Mat            C,LU;
 13:   MatInfo        info;
 14:   PetscInt       i,j,m = 3,n = 3,Ii,J;
 16:   PetscScalar    v,one = 1.0;
 17:   IS             perm,iperm;
 18:   Vec            x,u,b;
 19:   PetscReal      norm;
 20:   MatFactorInfo  luinfo;

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

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

 44:   MatFactorInfoInitialize(&luinfo);
 45:   luinfo.fill = 2.0;
 46:   luinfo.dtcol = 0.0;
 47:   luinfo.shiftnz = 0.0;
 48:   luinfo.zeropivot = 1.e-14;
 49:   luinfo.pivotinblocks = 1.0;
 50:   MatLUFactorSymbolic(C,perm,iperm,&luinfo,&LU);
 51:   MatLUFactorNumeric(C,&luinfo,&LU);
 52:   MatView(LU,PETSC_VIEWER_STDOUT_WORLD);

 54:   VecCreateSeq(PETSC_COMM_SELF,m*n,&u);
 55:   VecSet(u,one);
 56:   VecDuplicate(u,&x);
 57:   VecDuplicate(u,&b);

 59:   MatMult(C,u,b);
 60:   MatSolve(LU,b,x);

 62:   VecView(b,PETSC_VIEWER_STDOUT_SELF);
 63:   VecView(x,PETSC_VIEWER_STDOUT_SELF);

 65:   VecAXPY(x,-1.0,u);
 66:   VecNorm(x,NORM_2,&norm);
 67:   PetscPrintf(PETSC_COMM_SELF,"Norm of error %A\n",norm);

 69:   MatGetInfo(C,MAT_LOCAL,&info);
 70:   PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %D\n",(PetscInt)info.nz_used);
 71:   MatGetInfo(LU,MAT_LOCAL,&info);
 72:   PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %D\n",(PetscInt)info.nz_used);

 74:   VecDestroy(u);
 75:   VecDestroy(b);
 76:   VecDestroy(x);
 77:   ISDestroy(perm);
 78:   ISDestroy(iperm);
 79:   MatDestroy(C);
 80:   MatDestroy(LU);

 82:   PetscFinalize();
 83:   return 0;
 84: }