xref: /petsc/src/mat/tests/ex28.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Illustrate how to do one symbolic factorization and multiple numeric factorizations using same matrix structure. \n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
59371c9d4SSatish Balay int main(int argc, char **args) {
6dc6ed827SStefano Zampini   PetscInt      i, rstart, rend, N = 10, num_numfac = 5, col[3], k;
7c4762a1bSJed Brown   Mat           A[5], F;
8c4762a1bSJed Brown   Vec           u, x, b;
9c4762a1bSJed Brown   PetscMPIInt   rank;
10c4762a1bSJed Brown   PetscScalar   value[3];
11c4762a1bSJed Brown   PetscReal     norm, tol = 100 * PETSC_MACHINE_EPSILON;
12c4762a1bSJed Brown   IS            perm, iperm;
13c4762a1bSJed Brown   MatFactorInfo info;
14dc6ed827SStefano Zampini   MatFactorType facttype = MAT_FACTOR_LU;
15dc6ed827SStefano Zampini   char          solvertype[64];
16dc6ed827SStefano Zampini   char          factortype[64];
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* Create and assemble matrices, all have same data structure */
23c4762a1bSJed Brown   for (k = 0; k < num_numfac; k++) {
249566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &A[k]));
259566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(A[k], PETSC_DECIDE, PETSC_DECIDE, N, N));
269566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(A[k]));
279566063dSJacob Faibussowitsch     PetscCall(MatSetUp(A[k]));
289566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(A[k], &rstart, &rend));
29c4762a1bSJed Brown 
30dc6ed827SStefano Zampini     value[0] = -1.0 * (k + 1);
31dc6ed827SStefano Zampini     value[1] = 2.0 * (k + 1);
32dc6ed827SStefano Zampini     value[2] = -1.0 * (k + 1);
33c4762a1bSJed Brown     for (i = rstart; i < rend; i++) {
349371c9d4SSatish Balay       col[0] = i - 1;
359371c9d4SSatish Balay       col[1] = i;
369371c9d4SSatish Balay       col[2] = i + 1;
37c4762a1bSJed Brown       if (i == 0) {
389566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A[k], 1, &i, 2, col + 1, value + 1, INSERT_VALUES));
39c4762a1bSJed Brown       } else if (i == N - 1) {
409566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A[k], 1, &i, 2, col, value, INSERT_VALUES));
41c4762a1bSJed Brown       } else {
429566063dSJacob Faibussowitsch         PetscCall(MatSetValues(A[k], 1, &i, 3, col, value, INSERT_VALUES));
43c4762a1bSJed Brown       }
44c4762a1bSJed Brown     }
459566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A[k], MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A[k], MAT_FINAL_ASSEMBLY));
479566063dSJacob Faibussowitsch     PetscCall(MatSetOption(A[k], MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
48c4762a1bSJed Brown   }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* Create vectors */
519566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A[0], &x, &b));
529566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &u));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Set rhs vector b */
559566063dSJacob Faibussowitsch   PetscCall(VecSet(b, 1.0));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* Get a symbolic factor F from A[0] */
589566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(solvertype, "petsc", sizeof(solvertype)));
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", solvertype, sizeof(solvertype), NULL));
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-mat_factor_type", MatFactorTypes, (PetscEnum *)&facttype, NULL));
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(MatGetFactor(A[0], solvertype, facttype, &F));
63c4762a1bSJed Brown   /* test mumps options */
64dc6ed827SStefano Zampini #if defined(PETSC_HAVE_MUMPS)
659566063dSJacob Faibussowitsch   PetscCall(MatMumpsSetIcntl(F, 7, 5));
66c4762a1bSJed Brown #endif
679566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(factortype, MatFactorTypes[facttype], sizeof(factortype)));
689566063dSJacob Faibussowitsch   PetscCall(PetscStrtoupper(solvertype));
699566063dSJacob Faibussowitsch   PetscCall(PetscStrtoupper(factortype));
709566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %s %s:\n", solvertype, factortype));
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(MatFactorInfoInitialize(&info));
73c4762a1bSJed Brown   info.fill = 5.0;
749566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A[0], MATORDERINGNATURAL, &perm, &iperm));
75dc6ed827SStefano Zampini   switch (facttype) {
769371c9d4SSatish Balay   case MAT_FACTOR_LU: PetscCall(MatLUFactorSymbolic(F, A[0], perm, iperm, &info)); break;
779371c9d4SSatish Balay   case MAT_FACTOR_ILU: PetscCall(MatILUFactorSymbolic(F, A[0], perm, iperm, &info)); break;
789371c9d4SSatish Balay   case MAT_FACTOR_ICC: PetscCall(MatICCFactorSymbolic(F, A[0], perm, &info)); break;
799371c9d4SSatish Balay   case MAT_FACTOR_CHOLESKY: PetscCall(MatCholeskyFactorSymbolic(F, A[0], perm, &info)); break;
809371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not for factor type %s", factortype);
81dc6ed827SStefano Zampini   }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* Compute numeric factors using same F, then solve */
84c4762a1bSJed Brown   for (k = 0; k < num_numfac; k++) {
85dc6ed827SStefano Zampini     switch (facttype) {
86dc6ed827SStefano Zampini     case MAT_FACTOR_LU:
879371c9d4SSatish Balay     case MAT_FACTOR_ILU: PetscCall(MatLUFactorNumeric(F, A[k], &info)); break;
88dc6ed827SStefano Zampini     case MAT_FACTOR_ICC:
899371c9d4SSatish Balay     case MAT_FACTOR_CHOLESKY: PetscCall(MatCholeskyFactorNumeric(F, A[k], &info)); break;
909371c9d4SSatish Balay     default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not for factor type %s", factortype);
91dc6ed827SStefano Zampini     }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     /* Solve A[k] * x = b */
949566063dSJacob Faibussowitsch     PetscCall(MatSolve(F, b, x));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown     /* Check the residual */
979566063dSJacob Faibussowitsch     PetscCall(MatMult(A[k], x, u));
989566063dSJacob Faibussowitsch     PetscCall(VecAXPY(u, -1.0, b));
999566063dSJacob Faibussowitsch     PetscCall(VecNorm(u, NORM_INFINITY, &norm));
100*48a46eb9SPierre Jolivet     if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the %s numfact and solve: residual %g\n", k, factortype, (double)norm));
101c4762a1bSJed Brown   }
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Free data structures */
104*48a46eb9SPierre Jolivet   for (k = 0; k < num_numfac; k++) PetscCall(MatDestroy(&A[k]));
1059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
1069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1079566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
1089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1119566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
112b122ec5aSJacob Faibussowitsch   return 0;
113c4762a1bSJed Brown }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown /*TEST
116c4762a1bSJed Brown 
117c4762a1bSJed Brown    test:
118c4762a1bSJed Brown 
119c4762a1bSJed Brown    test:
120c4762a1bSJed Brown       suffix: 2
121c4762a1bSJed Brown       args: -mat_solver_type superlu
122c4762a1bSJed Brown       requires: superlu
123c4762a1bSJed Brown 
124c4762a1bSJed Brown    test:
125c4762a1bSJed Brown       suffix: 3
126c4762a1bSJed Brown       nsize: 2
127c4762a1bSJed Brown       requires: mumps
128c4762a1bSJed Brown       args: -mat_solver_type mumps
129c4762a1bSJed Brown 
130dc6ed827SStefano Zampini    test:
131dc6ed827SStefano Zampini       suffix: 4
132dc6ed827SStefano Zampini       args: -mat_solver_type cusparse -mat_type aijcusparse -mat_factor_type {{lu cholesky ilu icc}separate output}
133dc6ed827SStefano Zampini       requires: cuda
134dc6ed827SStefano Zampini 
135c4762a1bSJed Brown TEST*/
136