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