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