1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests ILU and ICC factorization with and without matrix ordering on seqaij format, and illustrates drawing of matrix sparsity structure with MatView().\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,A; 16*c4762a1bSJed Brown PetscInt i,j,m = 5,n = 5,Ii,J,lf = 0; 17*c4762a1bSJed Brown PetscErrorCode ierr; 18*c4762a1bSJed Brown PetscBool LU=PETSC_FALSE,CHOLESKY,TRIANGULAR=PETSC_FALSE,MATDSPL=PETSC_FALSE,flg,matordering; 19*c4762a1bSJed Brown PetscScalar v; 20*c4762a1bSJed Brown IS row,col; 21*c4762a1bSJed Brown PetscViewer viewer1,viewer2; 22*c4762a1bSJed Brown MatFactorInfo info; 23*c4762a1bSJed Brown Vec x,y,b,ytmp; 24*c4762a1bSJed Brown PetscReal norm2,norm2_inplace, tol = 100.*PETSC_MACHINE_EPSILON; 25*c4762a1bSJed Brown PetscRandom rdm; 26*c4762a1bSJed Brown PetscMPIInt size; 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 29*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 30*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 31*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&C);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatSetSizes(C,m*n,m*n,m*n,m*n);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatSetUp(C);CHKERRQ(ierr); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ 44*c4762a1bSJed Brown for (i=0; i<m; i++) { 45*c4762a1bSJed Brown for (j=0; j<n; j++) { 46*c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 47*c4762a1bSJed Brown J = Ii - n; if (J>=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 48*c4762a1bSJed Brown J = Ii + n; if (J<m*n) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 49*c4762a1bSJed Brown J = Ii - 1; if (J>=0) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 50*c4762a1bSJed Brown J = Ii + 1; if (J<m*n) {ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 51*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown } 54*c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown ierr = MatIsSymmetric(C,0.0,&flg);CHKERRQ(ierr); 58*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric"); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown /* Create vectors for error checking */ 61*c4762a1bSJed Brown ierr = MatCreateVecs(C,&x,&b);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecDuplicate(x,&ytmp);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_SELF,&rdm);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = VecSetRandom(x,rdm);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatMult(C,x,b);CHKERRQ(ierr); 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-mat_ordering",&matordering);CHKERRQ(ierr); 70*c4762a1bSJed Brown if (matordering) { 71*c4762a1bSJed Brown ierr = MatGetOrdering(C,MATORDERINGRCM,&row,&col);CHKERRQ(ierr); 72*c4762a1bSJed Brown } else { 73*c4762a1bSJed Brown ierr = MatGetOrdering(C,MATORDERINGNATURAL,&row,&col);CHKERRQ(ierr); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-display_matrices",&MATDSPL);CHKERRQ(ierr); 77*c4762a1bSJed Brown if (MATDSPL) { 78*c4762a1bSJed Brown printf("original matrix:\n"); 79*c4762a1bSJed Brown ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = MatView(C,viewer1);CHKERRQ(ierr); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown /* Compute LU or ILU factor A */ 87*c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown info.fill = 1.0; 90*c4762a1bSJed Brown info.diagonal_fill = 0; 91*c4762a1bSJed Brown info.zeropivot = 0.0; 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-lu",&LU);CHKERRQ(ierr); 94*c4762a1bSJed Brown if (LU) { 95*c4762a1bSJed Brown printf("Test LU...\n"); 96*c4762a1bSJed Brown ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = MatLUFactorSymbolic(A,C,row,col,&info);CHKERRQ(ierr); 98*c4762a1bSJed Brown } else { 99*c4762a1bSJed Brown printf("Test ILU...\n"); 100*c4762a1bSJed Brown info.levels = lf; 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ILU,&A);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = MatILUFactorSymbolic(A,C,row,col,&info);CHKERRQ(ierr); 104*c4762a1bSJed Brown } 105*c4762a1bSJed Brown ierr = MatLUFactorNumeric(A,C,&info);CHKERRQ(ierr); 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown /* Solve A*y = b, then check the error */ 108*c4762a1bSJed Brown ierr = MatSolve(A,b,y);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown /* Test in-place ILU(0) and compare it with the out-place ILU(0) */ 114*c4762a1bSJed Brown if (!LU && lf==0) { 115*c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_COPY_VALUES,&A);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = MatILUFactor(A,row,col,&info);CHKERRQ(ierr); 117*c4762a1bSJed Brown /* 118*c4762a1bSJed Brown printf("In-place factored matrix:\n"); 119*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 120*c4762a1bSJed Brown */ 121*c4762a1bSJed Brown ierr = MatSolve(A,b,y);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2_inplace);CHKERRQ(ierr); 124*c4762a1bSJed Brown if (PetscAbs(norm2 - norm2_inplace) > tol) SETERRQ2(PETSC_COMM_SELF,1,"ILU(0) %g and in-place ILU(0) %g give different residuals",(double)norm2,(double)norm2_inplace); 125*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown /* Test Cholesky and ICC on seqaij matrix with matrix reordering on aij matrix C */ 129*c4762a1bSJed Brown CHOLESKY = LU; 130*c4762a1bSJed Brown if (CHOLESKY) { 131*c4762a1bSJed Brown printf("Test Cholesky...\n"); 132*c4762a1bSJed Brown lf = -1; 133*c4762a1bSJed Brown ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&A);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(A,C,row,&info);CHKERRQ(ierr); 135*c4762a1bSJed Brown } else { 136*c4762a1bSJed Brown printf("Test ICC...\n"); 137*c4762a1bSJed Brown info.levels = lf; 138*c4762a1bSJed Brown info.fill = 1.0; 139*c4762a1bSJed Brown info.diagonal_fill = 0; 140*c4762a1bSJed Brown info.zeropivot = 0.0; 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ICC,&A);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = MatICCFactorSymbolic(A,C,row,&info);CHKERRQ(ierr); 144*c4762a1bSJed Brown } 145*c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(A,C,&info);CHKERRQ(ierr); 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 148*c4762a1bSJed Brown if (lf == -1) { 149*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR);CHKERRQ(ierr); 150*c4762a1bSJed Brown if (TRIANGULAR) { 151*c4762a1bSJed Brown printf("Test MatForwardSolve...\n"); 152*c4762a1bSJed Brown ierr = MatForwardSolve(A,b,ytmp);CHKERRQ(ierr); 153*c4762a1bSJed Brown printf("Test MatBackwardSolve...\n"); 154*c4762a1bSJed Brown ierr = MatBackwardSolve(A,ytmp,y);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 157*c4762a1bSJed Brown if (norm2 > tol) { 158*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2);CHKERRQ(ierr); 159*c4762a1bSJed Brown } 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown } 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown ierr = MatSolve(A,b,y);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2);CHKERRQ(ierr); 167*c4762a1bSJed Brown if (lf == -1 && norm2 > tol) { 168*c4762a1bSJed Brown PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %d, residual %g\n",lf,norm2);CHKERRQ(ierr); 169*c4762a1bSJed Brown } 170*c4762a1bSJed Brown 171*c4762a1bSJed Brown /* Test in-place ICC(0) and compare it with the out-place ICC(0) */ 172*c4762a1bSJed Brown if (!CHOLESKY && lf==0 && !matordering) { 173*c4762a1bSJed Brown ierr = MatConvert(C,MATSBAIJ,MAT_INITIAL_MATRIX,&A);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = MatICCFactor(A,row,&info);CHKERRQ(ierr); 175*c4762a1bSJed Brown /* 176*c4762a1bSJed Brown printf("In-place factored matrix:\n"); 177*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 178*c4762a1bSJed Brown */ 179*c4762a1bSJed Brown ierr = MatSolve(A,b,y);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 181*c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm2_inplace);CHKERRQ(ierr); 182*c4762a1bSJed Brown if (PetscAbs(norm2 - norm2_inplace) > tol) SETERRQ2(PETSC_COMM_SELF,1,"ICC(0) %g and in-place ICC(0) %g give different residuals",(double)norm2,(double)norm2_inplace); 183*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 184*c4762a1bSJed Brown } 185*c4762a1bSJed Brown 186*c4762a1bSJed Brown /* Free data structures */ 187*c4762a1bSJed Brown ierr = ISDestroy(&row);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = ISDestroy(&col);CHKERRQ(ierr); 189*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 190*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer1);CHKERRQ(ierr); 191*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer2);CHKERRQ(ierr); 192*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 193*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = VecDestroy(&ytmp);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 197*c4762a1bSJed Brown ierr = PetscFinalize(); 198*c4762a1bSJed Brown return ierr; 199*c4762a1bSJed Brown } 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown 202*c4762a1bSJed Brown /*TEST 203*c4762a1bSJed Brown 204*c4762a1bSJed Brown test: 205*c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox 206*c4762a1bSJed Brown filter: grep -v "MPI processes" 207*c4762a1bSJed Brown 208*c4762a1bSJed Brown test: 209*c4762a1bSJed Brown suffix: 2 210*c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox -lu 211*c4762a1bSJed Brown 212*c4762a1bSJed Brown test: 213*c4762a1bSJed Brown suffix: 3 214*c4762a1bSJed Brown args: -mat_ordering -lu -triangular_solve 215*c4762a1bSJed Brown 216*c4762a1bSJed Brown test: 217*c4762a1bSJed Brown suffix: 4 218*c4762a1bSJed Brown 219*c4762a1bSJed Brown test: 220*c4762a1bSJed Brown suffix: 5 221*c4762a1bSJed Brown args: -lu 222*c4762a1bSJed Brown 223*c4762a1bSJed Brown test: 224*c4762a1bSJed Brown suffix: 6 225*c4762a1bSJed Brown args: -lu -triangular_solve 226*c4762a1bSJed Brown output_file: output/ex30_3.out 227*c4762a1bSJed Brown 228*c4762a1bSJed Brown TEST*/ 229