10b4b7b1cSBarry Smith static char help[] = "Illustrate how to do one symbolic factorization and multiple numeric factorizations using same matrix nonzero structure. \n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7dc6ed827SStefano Zampini PetscInt i, rstart, rend, N = 10, num_numfac = 5, col[3], k;
8c4762a1bSJed Brown Mat A[5], F;
9c4762a1bSJed Brown Vec u, x, b;
10c4762a1bSJed Brown PetscMPIInt rank;
11c4762a1bSJed Brown PetscScalar value[3];
12c4762a1bSJed Brown PetscReal norm, tol = 100 * PETSC_MACHINE_EPSILON;
13c4762a1bSJed Brown IS perm, iperm;
14c4762a1bSJed Brown MatFactorInfo info;
15dc6ed827SStefano Zampini MatFactorType facttype = MAT_FACTOR_LU;
16dc6ed827SStefano Zampini char solvertype[64];
17dc6ed827SStefano Zampini char factortype[64];
18c4762a1bSJed Brown
19327415f7SBarry Smith PetscFunctionBeginUser;
20c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
22*c13ae4d5SJunchao Zhang PetscCall(PetscOptionsGetReal(NULL, NULL, "-tol", &tol, NULL));
23c4762a1bSJed Brown
24c4762a1bSJed Brown /* Create and assemble matrices, all have same data structure */
25c4762a1bSJed Brown for (k = 0; k < num_numfac; k++) {
269566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A[k]));
279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A[k], PETSC_DECIDE, PETSC_DECIDE, N, N));
289566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A[k]));
299566063dSJacob Faibussowitsch PetscCall(MatSetUp(A[k]));
309566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A[k], &rstart, &rend));
31c4762a1bSJed Brown
32dc6ed827SStefano Zampini value[0] = -1.0 * (k + 1);
33dc6ed827SStefano Zampini value[1] = 2.0 * (k + 1);
34dc6ed827SStefano Zampini value[2] = -1.0 * (k + 1);
35c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
369371c9d4SSatish Balay col[0] = i - 1;
379371c9d4SSatish Balay col[1] = i;
389371c9d4SSatish Balay col[2] = i + 1;
39c4762a1bSJed Brown if (i == 0) {
409566063dSJacob Faibussowitsch PetscCall(MatSetValues(A[k], 1, &i, 2, col + 1, value + 1, INSERT_VALUES));
41c4762a1bSJed Brown } else if (i == N - 1) {
429566063dSJacob Faibussowitsch PetscCall(MatSetValues(A[k], 1, &i, 2, col, value, INSERT_VALUES));
43c4762a1bSJed Brown } else {
449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A[k], 1, &i, 3, col, value, INSERT_VALUES));
45c4762a1bSJed Brown }
46c4762a1bSJed Brown }
479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A[k], MAT_FINAL_ASSEMBLY));
489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A[k], MAT_FINAL_ASSEMBLY));
499566063dSJacob Faibussowitsch PetscCall(MatSetOption(A[k], MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
50c4762a1bSJed Brown }
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* Create vectors */
539566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A[0], &x, &b));
549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u));
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* Set rhs vector b */
579566063dSJacob Faibussowitsch PetscCall(VecSet(b, 1.0));
58c4762a1bSJed Brown
59c4762a1bSJed Brown /* Get a symbolic factor F from A[0] */
609566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(solvertype, "petsc", sizeof(solvertype)));
619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", solvertype, sizeof(solvertype), NULL));
629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-mat_factor_type", MatFactorTypes, (PetscEnum *)&facttype, NULL));
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A[0], solvertype, facttype, &F));
65c4762a1bSJed Brown /* test mumps options */
669566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(F, 7, 5));
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) {
76d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_LU:
77d71ae5a4SJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A[0], perm, iperm, &info));
78d71ae5a4SJacob Faibussowitsch break;
79d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_ILU:
80d71ae5a4SJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(F, A[0], perm, iperm, &info));
81d71ae5a4SJacob Faibussowitsch break;
82d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_ICC:
83d71ae5a4SJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(F, A[0], perm, &info));
84d71ae5a4SJacob Faibussowitsch break;
85d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY:
86d71ae5a4SJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A[0], perm, &info));
87d71ae5a4SJacob Faibussowitsch break;
88d71ae5a4SJacob Faibussowitsch default:
89d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not for factor type %s", factortype);
90dc6ed827SStefano Zampini }
91c4762a1bSJed Brown
92c4762a1bSJed Brown /* Compute numeric factors using same F, then solve */
93c4762a1bSJed Brown for (k = 0; k < num_numfac; k++) {
94dc6ed827SStefano Zampini switch (facttype) {
95dc6ed827SStefano Zampini case MAT_FACTOR_LU:
96d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_ILU:
97d71ae5a4SJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A[k], &info));
98d71ae5a4SJacob Faibussowitsch break;
99dc6ed827SStefano Zampini case MAT_FACTOR_ICC:
100d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY:
101d71ae5a4SJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A[k], &info));
102d71ae5a4SJacob Faibussowitsch break;
103d71ae5a4SJacob Faibussowitsch default:
104d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not for factor type %s", factortype);
105dc6ed827SStefano Zampini }
106c4762a1bSJed Brown
107c4762a1bSJed Brown /* Solve A[k] * x = b */
1089566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x));
109c4762a1bSJed Brown
110c4762a1bSJed Brown /* Check the residual */
1119566063dSJacob Faibussowitsch PetscCall(MatMult(A[k], x, u));
1129566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, b));
1139566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_INFINITY, &norm));
11448a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the %s numfact and solve: residual %g\n", k, factortype, (double)norm));
115c4762a1bSJed Brown }
116c4762a1bSJed Brown
117c4762a1bSJed Brown /* Free data structures */
11848a46eb9SPierre Jolivet for (k = 0; k < num_numfac; k++) PetscCall(MatDestroy(&A[k]));
1199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F));
1209566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm));
1219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm));
1229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
1259566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
126b122ec5aSJacob Faibussowitsch return 0;
127c4762a1bSJed Brown }
128c4762a1bSJed Brown
129c4762a1bSJed Brown /*TEST
130c4762a1bSJed Brown
131c4762a1bSJed Brown test:
132c4762a1bSJed Brown
133c4762a1bSJed Brown test:
134c4762a1bSJed Brown suffix: 2
135c4762a1bSJed Brown args: -mat_solver_type superlu
136c4762a1bSJed Brown requires: superlu
137c4762a1bSJed Brown
138*c13ae4d5SJunchao Zhang testset:
139c4762a1bSJed Brown nsize: 2
140c4762a1bSJed Brown requires: mumps
141c4762a1bSJed Brown args: -mat_solver_type mumps
142*c13ae4d5SJunchao Zhang output_file: output/ex28_3.out
143*c13ae4d5SJunchao Zhang
144*c13ae4d5SJunchao Zhang test:
145*c13ae4d5SJunchao Zhang suffix: 3
146*c13ae4d5SJunchao Zhang requires: !__float128
147*c13ae4d5SJunchao Zhang test:
148*c13ae4d5SJunchao Zhang suffix: 3_fp128
149*c13ae4d5SJunchao Zhang requires: __float128
150*c13ae4d5SJunchao Zhang args: -tol 1e-14
151c4762a1bSJed Brown
152dc6ed827SStefano Zampini test:
153dc6ed827SStefano Zampini suffix: 4
154dc6ed827SStefano Zampini args: -mat_solver_type cusparse -mat_type aijcusparse -mat_factor_type {{lu cholesky ilu icc}separate output}
155dc6ed827SStefano Zampini requires: cuda
156dc6ed827SStefano Zampini
157c4762a1bSJed Brown TEST*/
158