xref: /petsc/src/mat/tests/ex7.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
8*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
9*d71ae5a4SJacob Faibussowitsch {
10c4762a1bSJed Brown   Mat           C, LU;
11c4762a1bSJed Brown   MatInfo       info;
121a6d72e3SSatish Balay   PetscInt      i, j, m, n, Ii, J;
13c4762a1bSJed Brown   PetscScalar   v, one = 1.0;
14c4762a1bSJed Brown   IS            perm, iperm;
15c4762a1bSJed Brown   Vec           x, u, b;
161a6d72e3SSatish Balay   PetscReal     norm, fill;
17c4762a1bSJed Brown   MatFactorInfo luinfo;
18c4762a1bSJed Brown 
19327415f7SBarry Smith   PetscFunctionBeginUser;
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
211a6d72e3SSatish Balay 
22d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Mat test ex7 options", "Mat");
239371c9d4SSatish Balay   m    = 3;
249371c9d4SSatish Balay   n    = 3;
259371c9d4SSatish Balay   fill = 2.0;
269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-m", "Number of rows in grid", NULL, m, &m, NULL));
279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-n", "Number of columns in grid", NULL, n, &n, NULL));
289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-fill", "Expected fill ratio for factorization", NULL, fill, &fill, NULL));
291a6d72e3SSatish Balay 
30d0609cedSBarry Smith   PetscOptionsEnd();
311a6d72e3SSatish Balay 
32c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
339566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &C));
349566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
369566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
37c4762a1bSJed Brown   for (i = 0; i < m; i++) {
38c4762a1bSJed Brown     for (j = 0; j < n; j++) {
399371c9d4SSatish Balay       v  = -1.0;
409371c9d4SSatish Balay       Ii = j + n * i;
419371c9d4SSatish Balay       if (i > 0) {
429371c9d4SSatish Balay         J = Ii - n;
439371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
449371c9d4SSatish Balay       }
459371c9d4SSatish Balay       if (i < m - 1) {
469371c9d4SSatish Balay         J = Ii + n;
479371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
489371c9d4SSatish Balay       }
499371c9d4SSatish Balay       if (j > 0) {
509371c9d4SSatish Balay         J = Ii - 1;
519371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
529371c9d4SSatish Balay       }
539371c9d4SSatish Balay       if (j < n - 1) {
549371c9d4SSatish Balay         J = Ii + 1;
559371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
569371c9d4SSatish Balay       }
579371c9d4SSatish Balay       v = 4.0;
589371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
59c4762a1bSJed Brown     }
60c4762a1bSJed Brown   }
619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
639566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C, MATORDERINGRCM, &perm, &iperm));
649566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
659566063dSJacob Faibussowitsch   PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF));
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&luinfo));
68c4762a1bSJed Brown 
691a6d72e3SSatish Balay   luinfo.fill          = fill;
70c4762a1bSJed Brown   luinfo.dtcol         = 0.0;
71c4762a1bSJed Brown   luinfo.zeropivot     = 1.e-14;
72c4762a1bSJed Brown   luinfo.pivotinblocks = 1.0;
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(C, MATSOLVERPETSC, MAT_FACTOR_LU, &LU));
759566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(LU, C, perm, iperm, &luinfo));
769566063dSJacob Faibussowitsch   PetscCall(MatLUFactorNumeric(LU, C, &luinfo));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, m * n, &u));
799566063dSJacob Faibussowitsch   PetscCall(VecSet(u, one));
809566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &x));
819566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &b));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(MatMult(C, u, b));
849566063dSJacob Faibussowitsch   PetscCall(MatSolve(LU, b, x));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(VecView(b, PETSC_VIEWER_STDOUT_SELF));
879566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, -1.0, u));
909566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
919566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Norm of error %g\n", (double)norm));
92c4762a1bSJed Brown 
939566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(C, MAT_LOCAL, &info));
949566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "original matrix nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used));
959566063dSJacob Faibussowitsch   PetscCall(MatGetInfo(LU, MAT_LOCAL, &info));
969566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "factored matrix nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used));
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1019566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1029566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
1039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&LU));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
107b122ec5aSJacob Faibussowitsch   return 0;
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown /*TEST
111c4762a1bSJed Brown 
112c4762a1bSJed Brown    test:
1131a6d72e3SSatish Balay       suffix: 1
1148cc725e6SPierre Jolivet       filter: grep -v " MPI process"
1151a6d72e3SSatish Balay 
1161a6d72e3SSatish Balay    test:
1171a6d72e3SSatish Balay       suffix: 2
1181a6d72e3SSatish Balay       args: -m 1 -n 1 -fill 0.49
1198cc725e6SPierre Jolivet       filter: grep -v " MPI process"
120c4762a1bSJed Brown 
121c4762a1bSJed Brown TEST*/
122