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 PetscBool CHOLESKY=PETSC_FALSE,TRIANGULAR=PETSC_FALSE,flg; 18c4762a1bSJed Brown PetscScalar v; 19c4762a1bSJed Brown IS row,col; 20c4762a1bSJed Brown MatFactorInfo info; 21c4762a1bSJed Brown Vec x,y,b,ytmp; 22c4762a1bSJed Brown PetscReal norm2,tol = 100*PETSC_MACHINE_EPSILON; 23c4762a1bSJed Brown PetscRandom rdm; 24c4762a1bSJed Brown PetscMPIInt size; 25c4762a1bSJed Brown 269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 28*be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&C)); 349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m*n,m*n,m*n,m*n)); 359566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 369566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ 39c4762a1bSJed Brown for (i=0; i<m; i++) { 40c4762a1bSJed Brown for (j=0; j<n; j++) { 41c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 429566063dSJacob Faibussowitsch J = Ii - n; if (J>=0) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 439566063dSJacob Faibussowitsch J = Ii + n; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 449566063dSJacob Faibussowitsch J = Ii - 1; if (J>=0) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 459566063dSJacob Faibussowitsch J = Ii + 1; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 469566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(C,0.0,&flg)); 5328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric"); 549566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sC)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Create vectors for error checking */ 579566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C,&x,&b)); 589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&ytmp)); 609566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 619566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 629566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 639566063dSJacob Faibussowitsch PetscCall(MatMult(C,x,b)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C,MATORDERINGNATURAL,&row,&col)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* Compute CHOLESKY or ICC factor sA */ 689566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown info.fill = 1.0; 71c4762a1bSJed Brown info.diagonal_fill = 0; 72c4762a1bSJed Brown info.zeropivot = 0.0; 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-cholesky",&CHOLESKY)); 75c4762a1bSJed Brown if (CHOLESKY) { 769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test CHOLESKY...\n")); 779566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&sA)); 789566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(sA,sC,row,&info)); 79c4762a1bSJed Brown } else { 809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test ICC...\n")); 81c4762a1bSJed Brown info.levels = lf; 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(MatGetFactor(sC,MATSOLVERPETSC,MAT_FACTOR_ICC,&sA)); 849566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(sA,sC,row,&info)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(sA,sC,&info)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 89c4762a1bSJed Brown if (CHOLESKY) { 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR)); 91c4762a1bSJed Brown if (TRIANGULAR) { 929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test MatForwardSolve...\n")); 939566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(sA,b,ytmp)); 949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Test MatBackwardSolve...\n")); 959566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(sA,ytmp,y)); 969566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 979566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 98c4762a1bSJed Brown if (norm2 > tol) { 999566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2)); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(MatSolve(sA,b,y)); 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sC)); 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 1079566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 1089566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 109c4762a1bSJed Brown if (lf == -1 && norm2 > tol) { 1109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Free data structures */ 1149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 1179566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 1219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1229566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 123b122ec5aSJacob Faibussowitsch return 0; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown /*TEST 127c4762a1bSJed Brown 128c4762a1bSJed Brown test: 129c4762a1bSJed Brown output_file: output/ex128.out 130c4762a1bSJed Brown 131c4762a1bSJed Brown test: 132c4762a1bSJed Brown suffix: 2 133c4762a1bSJed Brown args: -cholesky -triangular_solve 134c4762a1bSJed Brown 135c4762a1bSJed Brown TEST*/ 136