1*0b4b7b1cSBarry 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 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)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Create and assemble matrices, all have same data structure */ 24c4762a1bSJed Brown for (k = 0; k < num_numfac; k++) { 259566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A[k])); 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A[k], PETSC_DECIDE, PETSC_DECIDE, N, N)); 279566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A[k])); 289566063dSJacob Faibussowitsch PetscCall(MatSetUp(A[k])); 299566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A[k], &rstart, &rend)); 30c4762a1bSJed Brown 31dc6ed827SStefano Zampini value[0] = -1.0 * (k + 1); 32dc6ed827SStefano Zampini value[1] = 2.0 * (k + 1); 33dc6ed827SStefano Zampini value[2] = -1.0 * (k + 1); 34c4762a1bSJed Brown for (i = rstart; i < rend; i++) { 359371c9d4SSatish Balay col[0] = i - 1; 369371c9d4SSatish Balay col[1] = i; 379371c9d4SSatish Balay col[2] = i + 1; 38c4762a1bSJed Brown if (i == 0) { 399566063dSJacob Faibussowitsch PetscCall(MatSetValues(A[k], 1, &i, 2, col + 1, value + 1, INSERT_VALUES)); 40c4762a1bSJed Brown } else if (i == N - 1) { 419566063dSJacob Faibussowitsch PetscCall(MatSetValues(A[k], 1, &i, 2, col, value, INSERT_VALUES)); 42c4762a1bSJed Brown } else { 439566063dSJacob Faibussowitsch PetscCall(MatSetValues(A[k], 1, &i, 3, col, value, INSERT_VALUES)); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown } 469566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A[k], MAT_FINAL_ASSEMBLY)); 479566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A[k], MAT_FINAL_ASSEMBLY)); 489566063dSJacob Faibussowitsch PetscCall(MatSetOption(A[k], MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Create vectors */ 529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A[0], &x, &b)); 539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &u)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Set rhs vector b */ 569566063dSJacob Faibussowitsch PetscCall(VecSet(b, 1.0)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* Get a symbolic factor F from A[0] */ 599566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(solvertype, "petsc", sizeof(solvertype))); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-mat_solver_type", solvertype, sizeof(solvertype), NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEnum(NULL, NULL, "-mat_factor_type", MatFactorTypes, (PetscEnum *)&facttype, NULL)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A[0], solvertype, facttype, &F)); 64c4762a1bSJed Brown /* test mumps options */ 65dc6ed827SStefano Zampini #if defined(PETSC_HAVE_MUMPS) 669566063dSJacob Faibussowitsch PetscCall(MatMumpsSetIcntl(F, 7, 5)); 67c4762a1bSJed Brown #endif 689566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(factortype, MatFactorTypes[facttype], sizeof(factortype))); 699566063dSJacob Faibussowitsch PetscCall(PetscStrtoupper(solvertype)); 709566063dSJacob Faibussowitsch PetscCall(PetscStrtoupper(factortype)); 719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %s %s:\n", solvertype, factortype)); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&info)); 74c4762a1bSJed Brown info.fill = 5.0; 759566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A[0], MATORDERINGNATURAL, &perm, &iperm)); 76dc6ed827SStefano Zampini switch (facttype) { 77d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_LU: 78d71ae5a4SJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F, A[0], perm, iperm, &info)); 79d71ae5a4SJacob Faibussowitsch break; 80d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_ILU: 81d71ae5a4SJacob Faibussowitsch PetscCall(MatILUFactorSymbolic(F, A[0], perm, iperm, &info)); 82d71ae5a4SJacob Faibussowitsch break; 83d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_ICC: 84d71ae5a4SJacob Faibussowitsch PetscCall(MatICCFactorSymbolic(F, A[0], perm, &info)); 85d71ae5a4SJacob Faibussowitsch break; 86d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY: 87d71ae5a4SJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F, A[0], perm, &info)); 88d71ae5a4SJacob Faibussowitsch break; 89d71ae5a4SJacob Faibussowitsch default: 90d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not for factor type %s", factortype); 91dc6ed827SStefano Zampini } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Compute numeric factors using same F, then solve */ 94c4762a1bSJed Brown for (k = 0; k < num_numfac; k++) { 95dc6ed827SStefano Zampini switch (facttype) { 96dc6ed827SStefano Zampini case MAT_FACTOR_LU: 97d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_ILU: 98d71ae5a4SJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F, A[k], &info)); 99d71ae5a4SJacob Faibussowitsch break; 100dc6ed827SStefano Zampini case MAT_FACTOR_ICC: 101d71ae5a4SJacob Faibussowitsch case MAT_FACTOR_CHOLESKY: 102d71ae5a4SJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F, A[k], &info)); 103d71ae5a4SJacob Faibussowitsch break; 104d71ae5a4SJacob Faibussowitsch default: 105d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Not for factor type %s", factortype); 106dc6ed827SStefano Zampini } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Solve A[k] * x = b */ 1099566063dSJacob Faibussowitsch PetscCall(MatSolve(F, b, x)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Check the residual */ 1129566063dSJacob Faibussowitsch PetscCall(MatMult(A[k], x, u)); 1139566063dSJacob Faibussowitsch PetscCall(VecAXPY(u, -1.0, b)); 1149566063dSJacob Faibussowitsch PetscCall(VecNorm(u, NORM_INFINITY, &norm)); 11548a46eb9SPierre Jolivet if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT "-the %s numfact and solve: residual %g\n", k, factortype, (double)norm)); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Free data structures */ 11948a46eb9SPierre Jolivet for (k = 0; k < num_numfac; k++) PetscCall(MatDestroy(&A[k])); 1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1219566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1229566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1269566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 127b122ec5aSJacob Faibussowitsch return 0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /*TEST 131c4762a1bSJed Brown 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown 134c4762a1bSJed Brown test: 135c4762a1bSJed Brown suffix: 2 136c4762a1bSJed Brown args: -mat_solver_type superlu 137c4762a1bSJed Brown requires: superlu 138c4762a1bSJed Brown 139c4762a1bSJed Brown test: 140c4762a1bSJed Brown suffix: 3 141c4762a1bSJed Brown nsize: 2 142c4762a1bSJed Brown requires: mumps 143c4762a1bSJed Brown args: -mat_solver_type mumps 144c4762a1bSJed Brown 145dc6ed827SStefano Zampini test: 146dc6ed827SStefano Zampini suffix: 4 147dc6ed827SStefano Zampini args: -mat_solver_type cusparse -mat_type aijcusparse -mat_factor_type {{lu cholesky ilu icc}separate output} 148dc6ed827SStefano Zampini requires: cuda 149dc6ed827SStefano Zampini 150c4762a1bSJed Brown TEST*/ 151