1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqsbaij format. Modified from ex30.c\n\ 3c4762a1bSJed Brown Input parameters are:\n\ 4c4762a1bSJed Brown -lf <level> : level of fill for ILU (default is 0)\n\ 5c4762a1bSJed Brown -lu : use full LU or Cholesky factorization\n\ 6c4762a1bSJed Brown -m <value>,-n <value> : grid dimensions\n\ 7c4762a1bSJed Brown Note that most users should employ the KSP interface to the\n\ 8c4762a1bSJed Brown linear solvers instead of using the factorization routines\n\ 9c4762a1bSJed Brown directly.\n\n"; 10c4762a1bSJed Brown 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown int main(int argc,char **args) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown Mat C,sC,sA; 16c4762a1bSJed Brown PetscInt i,j,m = 5,n = 5,Ii,J,lf = 0; 17c4762a1bSJed Brown PetscErrorCode ierr; 18c4762a1bSJed Brown PetscBool CHOLESKY=PETSC_FALSE,TRIANGULAR=PETSC_FALSE,flg; 19c4762a1bSJed Brown PetscScalar v; 20c4762a1bSJed Brown IS row,col; 21c4762a1bSJed Brown MatFactorInfo info; 22c4762a1bSJed Brown Vec x,y,b,ytmp; 23c4762a1bSJed Brown PetscReal norm2,tol = 100*PETSC_MACHINE_EPSILON; 24c4762a1bSJed Brown PetscRandom rdm; 25c4762a1bSJed Brown PetscMPIInt size; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 28*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 292c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL)); 33c4762a1bSJed Brown 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&C)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,m*n,m*n,m*n,m*n)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ 40c4762a1bSJed Brown for (i=0; i<m; i++) { 41c4762a1bSJed Brown for (j=0; j<n; j++) { 42c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 43*5f80ce2aSJacob Faibussowitsch J = Ii - n; if (J>=0) CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 44*5f80ce2aSJacob Faibussowitsch J = Ii + n; if (J<m*n) CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 45*5f80ce2aSJacob Faibussowitsch J = Ii - 1; if (J>=0) CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 46*5f80ce2aSJacob Faibussowitsch J = Ii + 1; if (J<m*n) CHKERRQ(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 47*5f80ce2aSJacob Faibussowitsch v = 4.0; CHKERRQ(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown } 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 52c4762a1bSJed Brown 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatIsSymmetric(C,0.0,&flg)); 542c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric"); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sC)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Create vectors for error checking */ 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(C,&x,&b)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&ytmp)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,rdm)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,x,b)); 65c4762a1bSJed Brown 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(C,MATORDERINGNATURAL,&row,&col)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Compute CHOLESKY or ICC factor sA */ 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInfoInitialize(&info)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown info.fill = 1.0; 72c4762a1bSJed Brown info.diagonal_fill = 0; 73c4762a1bSJed Brown info.zeropivot = 0.0; 74c4762a1bSJed Brown 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-cholesky",&CHOLESKY)); 76c4762a1bSJed Brown if (CHOLESKY) { 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Test CHOLESKY...\n")); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sA)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(sA,sC,row,&info)); 80c4762a1bSJed Brown } else { 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Test ICC...\n")); 82c4762a1bSJed Brown info.levels = lf; 83c4762a1bSJed Brown 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_ICC,&sA)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(sA,sC,row,&info)); 86c4762a1bSJed Brown } 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(sA,sC,&info)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 90c4762a1bSJed Brown if (CHOLESKY) { 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR)); 92c4762a1bSJed Brown if (TRIANGULAR) { 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Test MatForwardSolve...\n")); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatForwardSolve(sA,b,ytmp)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Test MatBackwardSolve...\n")); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatBackwardSolve(sA,ytmp,y)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 99c4762a1bSJed Brown if (norm2 > tol) { 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(sA,b,y)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sC)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&sA)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 110c4762a1bSJed Brown if (lf == -1 && norm2 > tol) { 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2)); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Free data structures */ 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&row)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&col)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ytmp)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 123c4762a1bSJed Brown ierr = PetscFinalize(); 124c4762a1bSJed Brown return ierr; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*TEST 128c4762a1bSJed Brown 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown output_file: output/ex128.out 131c4762a1bSJed Brown 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown suffix: 2 134c4762a1bSJed Brown args: -cholesky -triangular_solve 135c4762a1bSJed Brown 136c4762a1bSJed Brown TEST*/ 137