static char help[] = "Tests ILU and ICC factorization with matrix ordering, and illustrates drawing of matrix sparsity structure with MatView().\n\ Input parameters are:\n\ -lf : level of fill for ILU (default is 0)\n\ -lu : use full LU or Cholesky factorization\n\ -m ,-n : grid dimensions\n\ Note that most users should employ the KSP interface to the\n\ linear solvers instead of using the factorization routines\n\ directly.\n\n"; #include "petscmat.h" #undef __FUNCT__ #define __FUNCT__ "main" int main(int argc,char **args) { Mat C,A,sC,sA;; PetscInt i,j,m = 5,n = 5,Ii,J,lf = 0; PetscErrorCode ierr; PetscTruth LU=PETSC_FALSE,flg; PetscScalar v; IS row,col; PetscViewer viewer1,viewer2; MatFactorInfo info; Vec x,y,b,ytmp; PetscReal norm2,norm2_inplace; PetscRandom rdm; PetscInitialize(&argc,&args,(char *)0,help); ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(PETSC_NULL,"-lf",&lf,PETSC_NULL);CHKERRQ(ierr); ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1);CHKERRQ(ierr); ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr); ierr = MatSetSizes(C,m*n,m*n,m*n,m*n);CHKERRQ(ierr); ierr = MatSetFromOptions(C);CHKERRQ(ierr); ierr = MatSeqBDiagSetPreallocation(C,0,1,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation(C,5,PETSC_NULL);CHKERRQ(ierr); /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ for (i=0; i=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} J = Ii + n; if (J=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} J = Ii + 1; if (J 1.e-16) SETERRQ2(1,"ILU(0) %G and in-place ILU(0) %G give different residuals",norm2,norm2_inplace); ierr = MatDestroy(A);CHKERRQ(ierr); } /* Test Cholesky and ICC on seqaij matrix with matrix reordering */ if (LU){ lf = -1; ierr = MatCholeskyFactorSymbolic(C,row,&info,&A);CHKERRQ(ierr); } else { info.levels = lf; info.fill = 1.0; info.diagonal_fill = 0; info.shiftnz = 0; info.zeropivot = 0.0; ierr = MatICCFactorSymbolic(C,row,&info,&A);CHKERRQ(ierr); } ierr = MatCholeskyFactorNumeric(C,&info,&A);CHKERRQ(ierr); /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering */ if (lf == -1){ ierr = MatForwardSolve(A,b,ytmp);CHKERRQ(ierr); ierr = MatBackwardSolve(A,ytmp,y);CHKERRQ(ierr); ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (norm2 > 1.e-14){ ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%G\n",norm2);CHKERRQ(ierr); } } ierr = MatSolve(A,b,y);CHKERRQ(ierr); ierr = MatDestroy(A);CHKERRQ(ierr); ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (lf == -1 && norm2 > 1.e-14){ PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %d, residual %g\n",lf,norm2);CHKERRQ(ierr); } /* Test Cholesky and ICC on seqaij matrix without matrix reordering */ ierr = ISDestroy(row);CHKERRQ(ierr); ierr = ISDestroy(col);CHKERRQ(ierr); ierr = MatGetOrdering(C,MATORDERING_NATURAL,&row,&col);CHKERRQ(ierr); if (LU){ lf = -1; ierr = MatCholeskyFactorSymbolic(C,row,&info,&A);CHKERRQ(ierr); } else { info.levels = lf; info.fill = 1.0; info.diagonal_fill = 0; info.shiftnz = 0; info.zeropivot = 0.0; ierr = MatICCFactorSymbolic(C,row,&info,&A);CHKERRQ(ierr); } ierr = MatCholeskyFactorNumeric(C,&info,&A);CHKERRQ(ierr); /* test MatForwardSolve() and MatBackwardSolve() */ if (lf == -1){ ierr = MatForwardSolve(A,b,ytmp);CHKERRQ(ierr); ierr = MatBackwardSolve(A,ytmp,y);CHKERRQ(ierr); ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (norm2 > 1.e-14){ ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%G\n",norm2);CHKERRQ(ierr); } } /* Test MatSolve() */ ierr = MatSolve(A,b,y);CHKERRQ(ierr); ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); if (lf == -1 && norm2 > 1.e-14){ printf(" SEQAIJ: Cholesky/ICC levels %d, residual %g\n",lf,norm2);CHKERRQ(ierr); } /* Test Cholesky and ICC on seqsbaij matrix without matrix reordering */ if (LU){ ierr = MatCholeskyFactorSymbolic(sC,row,&info,&sA);CHKERRQ(ierr); } else { ierr = MatICCFactorSymbolic(sC,row,&info,&sA);CHKERRQ(ierr); } ierr = MatCholeskyFactorNumeric(sC,&info,&sA);CHKERRQ(ierr); ierr = MatEqual(A,sA,&flg);CHKERRQ(ierr); if (!flg) SETERRQ(1,"CholeskyFactors for aij and sbaij matrices are different"); ierr = MatDestroy(sC);CHKERRQ(ierr); ierr = MatDestroy(sA);CHKERRQ(ierr); ierr = MatDestroy(A);CHKERRQ(ierr); /* Free data structures */ ierr = MatDestroy(C);CHKERRQ(ierr); ierr = ISDestroy(row);CHKERRQ(ierr); ierr = ISDestroy(col);CHKERRQ(ierr); ierr = PetscViewerDestroy(viewer1);CHKERRQ(ierr); ierr = PetscViewerDestroy(viewer2);CHKERRQ(ierr); ierr = PetscRandomDestroy(rdm);CHKERRQ(ierr); ierr = VecDestroy(x);CHKERRQ(ierr); ierr = VecDestroy(y);CHKERRQ(ierr); ierr = VecDestroy(ytmp);CHKERRQ(ierr); ierr = VecDestroy(b);CHKERRQ(ierr); ierr = PetscFinalize();CHKERRQ(ierr); return 0; }