1 2 static char help[] = "Tests matrix factorization. Note that most users should\n\ 3 employ the KSP interface to the linear solvers instead of using the factorization\n\ 4 routines directly.\n\n"; 5 6 #include <petscmat.h> 7 8 int main(int argc,char **args) 9 { 10 Mat C,LU; 11 MatInfo info; 12 PetscInt i,j,m,n,Ii,J; 13 PetscScalar v,one = 1.0; 14 IS perm,iperm; 15 Vec x,u,b; 16 PetscReal norm,fill; 17 MatFactorInfo luinfo; 18 19 PetscFunctionBeginUser; 20 PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 21 22 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Mat test ex7 options","Mat"); 23 m = 3; n = 3; fill = 2.0; 24 PetscCall(PetscOptionsInt("-m","Number of rows in grid",NULL,m,&m,NULL)); 25 PetscCall(PetscOptionsInt("-n","Number of columns in grid",NULL,n,&n,NULL)); 26 PetscCall(PetscOptionsReal("-fill","Expected fill ratio for factorization",NULL,fill,&fill,NULL)); 27 28 PetscOptionsEnd(); 29 30 /* Create the matrix for the five point stencil, YET AGAIN */ 31 PetscCall(MatCreate(PETSC_COMM_SELF,&C)); 32 PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 33 PetscCall(MatSetFromOptions(C)); 34 PetscCall(MatSetUp(C)); 35 for (i=0; i<m; i++) { 36 for (j=0; j<n; j++) { 37 v = -1.0; Ii = j + n*i; 38 if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 39 if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 40 if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 41 if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 42 v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 43 } 44 } 45 PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 46 PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 47 PetscCall(MatGetOrdering(C,MATORDERINGRCM,&perm,&iperm)); 48 PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 49 PetscCall(ISView(perm,PETSC_VIEWER_STDOUT_SELF)); 50 51 PetscCall(MatFactorInfoInitialize(&luinfo)); 52 53 luinfo.fill = fill; 54 luinfo.dtcol = 0.0; 55 luinfo.zeropivot = 1.e-14; 56 luinfo.pivotinblocks = 1.0; 57 58 PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&LU)); 59 PetscCall(MatLUFactorSymbolic(LU,C,perm,iperm,&luinfo)); 60 PetscCall(MatLUFactorNumeric(LU,C,&luinfo)); 61 62 PetscCall(VecCreateSeq(PETSC_COMM_SELF,m*n,&u)); 63 PetscCall(VecSet(u,one)); 64 PetscCall(VecDuplicate(u,&x)); 65 PetscCall(VecDuplicate(u,&b)); 66 67 PetscCall(MatMult(C,u,b)); 68 PetscCall(MatSolve(LU,b,x)); 69 70 PetscCall(VecView(b,PETSC_VIEWER_STDOUT_SELF)); 71 PetscCall(VecView(x,PETSC_VIEWER_STDOUT_SELF)); 72 73 PetscCall(VecAXPY(x,-1.0,u)); 74 PetscCall(VecNorm(x,NORM_2,&norm)); 75 PetscCall(PetscPrintf(PETSC_COMM_SELF,"Norm of error %g\n",(double)norm)); 76 77 PetscCall(MatGetInfo(C,MAT_LOCAL,&info)); 78 PetscCall(PetscPrintf(PETSC_COMM_SELF,"original matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used)); 79 PetscCall(MatGetInfo(LU,MAT_LOCAL,&info)); 80 PetscCall(PetscPrintf(PETSC_COMM_SELF,"factored matrix nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used)); 81 82 PetscCall(VecDestroy(&u)); 83 PetscCall(VecDestroy(&b)); 84 PetscCall(VecDestroy(&x)); 85 PetscCall(ISDestroy(&perm)); 86 PetscCall(ISDestroy(&iperm)); 87 PetscCall(MatDestroy(&C)); 88 PetscCall(MatDestroy(&LU)); 89 90 PetscCall(PetscFinalize()); 91 return 0; 92 } 93 94 /*TEST 95 96 test: 97 suffix: 1 98 filter: grep -v " MPI process" 99 100 test: 101 suffix: 2 102 args: -m 1 -n 1 -fill 0.49 103 filter: grep -v " MPI process" 104 105 TEST*/ 106