1c4762a1bSJed Brown 2c4762a1bSJed 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\ 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,A; 16c4762a1bSJed Brown PetscInt i,j,m = 5,n = 5,Ii,J,lf = 0; 17c4762a1bSJed Brown PetscBool LU=PETSC_FALSE,CHOLESKY,TRIANGULAR=PETSC_FALSE,MATDSPL=PETSC_FALSE,flg,matordering; 18c4762a1bSJed Brown PetscScalar v; 19c4762a1bSJed Brown IS row,col; 20c4762a1bSJed Brown PetscViewer viewer1,viewer2; 21c4762a1bSJed Brown MatFactorInfo info; 22c4762a1bSJed Brown Vec x,y,b,ytmp; 23c4762a1bSJed Brown PetscReal norm2,norm2_inplace, tol = 100.*PETSC_MACHINE_EPSILON; 24c4762a1bSJed Brown PetscRandom rdm; 25c4762a1bSJed Brown PetscMPIInt size; 26c4762a1bSJed Brown 27*327415f7SBarry Smith PetscFunctionBeginUser; 289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 30be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-lf",&lf,NULL)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,0,0,400,400,&viewer1)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawOpen(PETSC_COMM_SELF,0,0,400,0,400,400,&viewer2)); 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&C)); 399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m*n,m*n,m*n,m*n)); 409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 419566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Create matrix C in seqaij format and sC in seqsbaij. (This is five-point stencil with some extra elements) */ 44c4762a1bSJed Brown for (i=0; i<m; i++) { 45c4762a1bSJed Brown for (j=0; j<n; j++) { 46c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 479566063dSJacob Faibussowitsch J = Ii - n; if (J>=0) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 489566063dSJacob Faibussowitsch J = Ii + n; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 499566063dSJacob Faibussowitsch J = Ii - 1; if (J>=0) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 509566063dSJacob Faibussowitsch J = Ii + 1; if (J<m*n) PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); 519566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 52c4762a1bSJed Brown } 53c4762a1bSJed Brown } 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(C,0.0,&flg)); 5828b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_SUP,"C is non-symmetric"); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create vectors for error checking */ 619566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C,&x,&b)); 629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&ytmp)); 649566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 659566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 669566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rdm)); 679566063dSJacob Faibussowitsch PetscCall(MatMult(C,x,b)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mat_ordering",&matordering)); 70c4762a1bSJed Brown if (matordering) { 719566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C,MATORDERINGRCM,&row,&col)); 72c4762a1bSJed Brown } else { 739566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(C,MATORDERINGNATURAL,&row,&col)); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-display_matrices",&MATDSPL)); 77c4762a1bSJed Brown if (MATDSPL) { 78c4762a1bSJed Brown printf("original matrix:\n"); 799566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_INFO)); 809566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF)); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_SELF)); 829566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF)); 839566063dSJacob Faibussowitsch PetscCall(MatView(C,viewer1)); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Compute LU or ILU factor A */ 879566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown info.fill = 1.0; 90c4762a1bSJed Brown info.diagonal_fill = 0; 91c4762a1bSJed Brown info.zeropivot = 0.0; 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-lu",&LU)); 94c4762a1bSJed Brown if (LU) { 95c4762a1bSJed Brown printf("Test LU...\n"); 969566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&A)); 979566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(A,C,row,col,&info)); 98c4762a1bSJed Brown } else { 99c4762a1bSJed Brown printf("Test ILU...\n"); 100c4762a1bSJed Brown info.levels = lf; 101c4762a1bSJed Brown 1029566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ILU,&A)); 1039566063dSJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(A,C,row,col,&info)); 104c4762a1bSJed Brown } 1059566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(A,C,&info)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* Solve A*y = b, then check the error */ 1089566063dSJacob Faibussowitsch PetscCall(MatSolve(A,b,y)); 1099566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 1109566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 1119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* Test in-place ILU(0) and compare it with the out-place ILU(0) */ 114c4762a1bSJed Brown if (!LU && lf==0) { 1159566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&A)); 1169566063dSJacob Faibussowitsch PetscCall(MatILUFactor(A,row,col,&info)); 117c4762a1bSJed Brown /* 118c4762a1bSJed Brown printf("In-place factored matrix:\n"); 1199566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_SELF)); 120c4762a1bSJed Brown */ 1219566063dSJacob Faibussowitsch PetscCall(MatSolve(A,b,y)); 1229566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 1239566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2_inplace)); 12408401ef6SPierre Jolivet PetscCheck(PetscAbs(norm2 - norm2_inplace) <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ILU(0) %g and in-place ILU(0) %g give different residuals",(double)norm2,(double)norm2_inplace); 1259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* Test Cholesky and ICC on seqaij matrix with matrix reordering on aij matrix C */ 129c4762a1bSJed Brown CHOLESKY = LU; 130c4762a1bSJed Brown if (CHOLESKY) { 131c4762a1bSJed Brown printf("Test Cholesky...\n"); 132c4762a1bSJed Brown lf = -1; 1339566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_CHOLESKY,&A)); 1349566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(A,C,row,&info)); 135c4762a1bSJed Brown } else { 136c4762a1bSJed Brown printf("Test ICC...\n"); 137c4762a1bSJed Brown info.levels = lf; 138c4762a1bSJed Brown info.fill = 1.0; 139c4762a1bSJed Brown info.diagonal_fill = 0; 140c4762a1bSJed Brown info.zeropivot = 0.0; 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_ICC,&A)); 1439566063dSJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(A,C,row,&info)); 144c4762a1bSJed Brown } 1459566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(A,C,&info)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* test MatForwardSolve() and MatBackwardSolve() with matrix reordering on aij matrix C */ 148c4762a1bSJed Brown if (lf == -1) { 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-triangular_solve",&TRIANGULAR)); 150c4762a1bSJed Brown if (TRIANGULAR) { 151c4762a1bSJed Brown printf("Test MatForwardSolve...\n"); 1529566063dSJacob Faibussowitsch PetscCall(MatForwardSolve(A,b,ytmp)); 153c4762a1bSJed Brown printf("Test MatBackwardSolve...\n"); 1549566063dSJacob Faibussowitsch PetscCall(MatBackwardSolve(A,ytmp,y)); 1559566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 1569566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 157c4762a1bSJed Brown if (norm2 > tol) { 1589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatForwardSolve and BackwardSolve: Norm of error=%g\n",(double)norm2)); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(MatSolve(A,b,y)); 1649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1659566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 1669566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 167c4762a1bSJed Brown if (lf == -1 && norm2 > tol) { 1689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " reordered SEQAIJ: Cholesky/ICC levels %" PetscInt_FMT ", residual %g\n",lf,(double)norm2)); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* Test in-place ICC(0) and compare it with the out-place ICC(0) */ 172c4762a1bSJed Brown if (!CHOLESKY && lf==0 && !matordering) { 1739566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATSBAIJ,MAT_INITIAL_MATRIX,&A)); 1749566063dSJacob Faibussowitsch PetscCall(MatICCFactor(A,row,&info)); 175c4762a1bSJed Brown /* 176c4762a1bSJed Brown printf("In-place factored matrix:\n"); 1779566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 178c4762a1bSJed Brown */ 1799566063dSJacob Faibussowitsch PetscCall(MatSolve(A,b,y)); 1809566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,x)); 1819566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2_inplace)); 18208401ef6SPierre Jolivet PetscCheck(PetscAbs(norm2 - norm2_inplace) <= tol,PETSC_COMM_SELF,PETSC_ERR_PLIB,"ICC(0) %g and in-place ICC(0) %g give different residuals",(double)norm2,(double)norm2_inplace); 1839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Free data structures */ 1879566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&col)); 1899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1909566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer1)); 1919566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer2)); 1929566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ytmp)); 1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1979566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 198b122ec5aSJacob Faibussowitsch return 0; 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /*TEST 202c4762a1bSJed Brown 203c4762a1bSJed Brown test: 204c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox 2058cc725e6SPierre Jolivet filter: grep -v " MPI process" 206c4762a1bSJed Brown 207c4762a1bSJed Brown test: 208c4762a1bSJed Brown suffix: 2 209c4762a1bSJed Brown args: -mat_ordering -display_matrices -nox -lu 210c4762a1bSJed Brown 211c4762a1bSJed Brown test: 212c4762a1bSJed Brown suffix: 3 213c4762a1bSJed Brown args: -mat_ordering -lu -triangular_solve 214c4762a1bSJed Brown 215c4762a1bSJed Brown test: 216c4762a1bSJed Brown suffix: 4 217c4762a1bSJed Brown 218c4762a1bSJed Brown test: 219c4762a1bSJed Brown suffix: 5 220c4762a1bSJed Brown args: -lu 221c4762a1bSJed Brown 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 6 224c4762a1bSJed Brown args: -lu -triangular_solve 225c4762a1bSJed Brown output_file: output/ex30_3.out 226c4762a1bSJed Brown 227c4762a1bSJed Brown TEST*/ 228