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; 28ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 29*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 30c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL);CHKERRQ(ierr); 33c4762a1bSJed Brown 34c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = MatSetSizes(C,m*n,m*n,m*n,m*n);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 37c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 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; 43c4762a1bSJed Brown J = Ii - n; if (J>=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 44c4762a1bSJed Brown J = Ii + n; if (J<m*n) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 45c4762a1bSJed Brown J = Ii - 1; if (J>=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 46c4762a1bSJed Brown J = Ii + 1; if (J<m*n) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 47c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown } 50c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 52c4762a1bSJed Brown 53c4762a1bSJed Brown ierr = MatIsSymmetric(C,0.0,&flg);CHKERRQ(ierr); 54*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric"); 55c4762a1bSJed Brown ierr = MatConvert(C,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sC);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Create vectors for error checking */ 58c4762a1bSJed Brown ierr = MatCreateVecs(C,&x,&b);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = VecDuplicate(x,&ytmp);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = MatMult(C,x,b);CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown ierr = MatGetOrdering(C,MATORDERINGNATURAL,&row,&col);CHKERRQ(ierr); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Compute CHOLESKY or ICC factor sA */ 69c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown info.fill = 1.0; 72c4762a1bSJed Brown info.diagonal_fill = 0; 73c4762a1bSJed Brown info.zeropivot = 0.0; 74c4762a1bSJed Brown 75c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-cholesky",&CHOLESKY);CHKERRQ(ierr); 76c4762a1bSJed Brown if (CHOLESKY) { 77c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Test CHOLESKY...\n");CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sA);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(sA,sC,row,&info);CHKERRQ(ierr); 80c4762a1bSJed Brown } else { 81c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Test ICC...\n");CHKERRQ(ierr); 82c4762a1bSJed Brown info.levels = lf; 83c4762a1bSJed Brown 84c4762a1bSJed Brown ierr = MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_ICC,&sA);CHKERRQ(ierr); 85c4762a1bSJed Brown ierr = MatICCFactorSymbolic(sA,sC,row,&info);CHKERRQ(ierr); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(sA,sC,&info);CHKERRQ(ierr); 88c4762a1bSJed Brown 89c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 90c4762a1bSJed Brown if (CHOLESKY) { 91c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR);CHKERRQ(ierr); 92c4762a1bSJed Brown if (TRIANGULAR) { 93c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Test MatForwardSolve...\n");CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = MatForwardSolve(sA,b,ytmp);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Test MatBackwardSolve...\n");CHKERRQ(ierr); 96c4762a1bSJed Brown ierr = MatBackwardSolve(sA,ytmp,y);CHKERRQ(ierr); 97c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 99c4762a1bSJed Brown if (norm2 > tol) { 100c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2);CHKERRQ(ierr); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 103c4762a1bSJed Brown } 104c4762a1bSJed Brown 105c4762a1bSJed Brown ierr = MatSolve(sA,b,y);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = MatDestroy(&sC);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = MatDestroy(&sA);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 110c4762a1bSJed Brown if (lf == -1 && norm2 > tol) { 1118fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2);CHKERRQ(ierr); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Free data structures */ 115c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 116c4762a1bSJed Brown ierr = ISDestroy(&row);CHKERRQ(ierr); 117c4762a1bSJed Brown ierr = ISDestroy(&col);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = VecDestroy(&ytmp);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 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