Actual source code: ex87.c
1: #include <petscksp.h>
2: #include <petsc/private/petscimpl.h>
4: static char help[] = "Solves a saddle-point linear system using PCHPDDM.\n\n";
6: static PetscErrorCode MatAndISLoad(const char *prefix, const char *identifier, Mat A, IS is, Mat N, PetscMPIInt rank, PetscMPIInt size);
8: int main(int argc, char **args)
9: {
10: Vec b, x; /* computed solution and RHS */
11: Mat A[4], aux[2], S; /* linear system matrix */
12: KSP ksp, *subksp; /* linear solver context */
13: PC pc;
14: IS is[2];
15: PetscMPIInt rank, size;
16: PetscInt m, M, n, N, id = 0;
17: PetscViewer viewer;
18: const char *const system[] = {"elasticity", "stokes"};
19: /* "elasticity":
20: * 2D linear elasticity with rubber-like and steel-like material coefficients, i.e., Poisson's ratio \in {0.4999, 0.35} and Young's modulus \in {0.01 GPa, 200.0 GPa}
21: * discretized by order 2 (resp. 0) Lagrange finite elements in displacements (resp. pressure) on a triangle mesh
22: * "stokes":
23: * 2D lid-driven cavity with constant viscosity
24: * discretized by order 2 (resp. 1) Lagrange finite elements, i.e., lowest-order Taylor--Hood finite elements, in velocities (resp. pressure) on a triangle mesh
25: * if the option -empty_A11 is not set (or set to false), a pressure with a zero mean-value is computed
26: */
27: char dir[PETSC_MAX_PATH_LEN], prefix[PETSC_MAX_PATH_LEN];
28: PetscBool flg[4] = {PETSC_FALSE, PETSC_FALSE, PETSC_FALSE, PETSC_FALSE};
30: PetscFunctionBeginUser;
31: PetscCall(PetscInitialize(&argc, &args, NULL, help));
32: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
33: PetscCheck(size == 4, PETSC_COMM_WORLD, PETSC_ERR_USER, "This example requires 4 processes");
34: PetscCall(PetscOptionsGetEList(NULL, NULL, "-system", system, PETSC_STATIC_ARRAY_LENGTH(system), &id, NULL));
35: if (id == 1) PetscCall(PetscOptionsGetBool(NULL, NULL, "-empty_A11", flg, NULL));
36: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
37: for (PetscInt i = 0; i < 2; ++i) {
38: PetscCall(MatCreate(PETSC_COMM_WORLD, A + (i ? 3 : 0)));
39: PetscCall(ISCreate(PETSC_COMM_SELF, is + i));
40: PetscCall(MatCreate(PETSC_COMM_SELF, aux + i));
41: }
42: PetscCall(PetscStrncpy(dir, ".", sizeof(dir)));
43: PetscCall(PetscOptionsGetString(NULL, NULL, "-load_dir", dir, sizeof(dir), NULL));
44: /* loading matrices and auxiliary data for the diagonal blocks */
45: PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%s/%s", dir, id == 1 ? "B" : "A"));
46: PetscCall(MatAndISLoad(prefix, "00", A[0], is[0], aux[0], rank, size));
47: PetscCall(MatAndISLoad(prefix, "11", A[3], is[1], aux[1], rank, size));
48: /* loading the off-diagonal block with a coherent row/column layout */
49: PetscCall(MatCreate(PETSC_COMM_WORLD, A + 2));
50: PetscCall(MatGetLocalSize(A[0], &n, NULL));
51: PetscCall(MatGetSize(A[0], &N, NULL));
52: PetscCall(MatGetLocalSize(A[3], &m, NULL));
53: PetscCall(MatGetSize(A[3], &M, NULL));
54: PetscCall(MatSetSizes(A[2], m, n, M, N));
55: PetscCall(MatSetUp(A[2]));
56: PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%s/%s10.dat", dir, id == 1 ? "B" : "A"));
57: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, prefix, FILE_MODE_READ, &viewer));
58: PetscCall(MatLoad(A[2], viewer));
59: PetscCall(PetscViewerDestroy(&viewer));
60: /* transposing the off-diagonal block */
61: PetscCall(PetscOptionsGetBool(NULL, NULL, "-transpose", flg + 1, NULL));
62: PetscCall(PetscOptionsGetBool(NULL, NULL, "-permute", flg + 2, NULL));
63: PetscCall(PetscOptionsGetBool(NULL, NULL, "-explicit", flg + 3, NULL));
64: if (flg[1]) {
65: if (flg[2]) {
66: PetscCall(MatTranspose(A[2], MAT_INITIAL_MATRIX, A + 1));
67: PetscCall(MatDestroy(A + 2));
68: }
69: if (!flg[3]) PetscCall(MatCreateTranspose(A[2 - flg[2]], A + 1 + flg[2]));
70: else PetscCall(MatTranspose(A[2 - flg[2]], MAT_INITIAL_MATRIX, A + 1 + flg[2]));
71: } else {
72: if (flg[2]) {
73: PetscCall(MatHermitianTranspose(A[2], MAT_INITIAL_MATRIX, A + 1));
74: PetscCall(MatDestroy(A + 2));
75: }
76: if (!flg[3]) PetscCall(MatCreateHermitianTranspose(A[2 - flg[2]], A + 1 + flg[2]));
77: else PetscCall(MatHermitianTranspose(A[2 - flg[2]], MAT_INITIAL_MATRIX, A + 1 + flg[2]));
78: }
79: if (flg[0]) PetscCall(MatDestroy(A + 3));
80: /* global coefficient matrix */
81: PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, NULL, 2, NULL, A, &S));
82: PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
83: PetscCall(KSPSetOperators(ksp, S, S));
84: PetscCall(KSPGetPC(ksp, &pc));
85: /* outer preconditioner */
86: PetscCall(PCSetType(pc, PCFIELDSPLIT));
87: PetscCall(PCFieldSplitSetType(pc, PC_COMPOSITE_SCHUR));
88: PetscCall(PCFieldSplitSetSchurPre(pc, PC_FIELDSPLIT_SCHUR_PRE_SELF, NULL));
89: PetscCall(PCSetFromOptions(pc));
90: PetscCall(PCSetUp(pc));
91: PetscCall(PCFieldSplitGetSubKSP(pc, &n, &subksp));
92: PetscCall(KSPGetPC(subksp[0], &pc));
93: /* inner preconditioner associated to top-left block */
94: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
95: PetscCall(PCSetType(pc, PCHPDDM));
96: PetscCall(PCHPDDMSetAuxiliaryMat(pc, is[0], aux[0], NULL, NULL));
97: #endif
98: PetscCall(PCSetFromOptions(pc));
99: PetscCall(KSPGetPC(subksp[1], &pc));
100: /* inner preconditioner associated to Schur complement, which will be set internally to a PCKSP */
101: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
102: PetscCall(PCSetType(pc, PCHPDDM));
103: if (!flg[0]) PetscCall(PCHPDDMSetAuxiliaryMat(pc, is[1], aux[1], NULL, NULL));
104: #endif
105: PetscCall(PCSetFromOptions(pc));
106: PetscCall(PetscFree(subksp));
107: PetscCall(KSPSetFromOptions(ksp));
108: PetscCall(MatCreateVecs(S, &b, &x));
109: PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%s/rhs_%s.dat", dir, id == 1 ? "B" : "A"));
110: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, prefix, FILE_MODE_READ, &viewer));
111: PetscCall(VecLoad(b, viewer));
112: PetscCall(PetscViewerDestroy(&viewer));
113: PetscCall(KSPSolve(ksp, b, x));
114: flg[1] = PETSC_FALSE;
115: PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewer", flg + 1, NULL));
116: if (flg[1]) PetscCall(PCView(pc, PETSC_VIEWER_STDOUT_WORLD));
117: flg[1] = PETSC_FALSE;
118: PetscCall(PetscOptionsGetBool(NULL, NULL, "-successive_solves", flg + 1, NULL));
119: if (flg[1]) {
120: KSPConvergedReason reason[2];
121: PetscInt iterations[2];
122: PetscCall(KSPGetConvergedReason(ksp, reason));
123: PetscCall(KSPGetTotalIterations(ksp, iterations));
124: PetscCall(KSPMonitorCancel(ksp));
125: PetscCall(PetscOptionsClearValue(NULL, "-ksp_monitor"));
126: PetscCall(PetscObjectStateIncrease((PetscObject)S));
127: PetscCall(KSPSetUp(ksp));
128: PetscCall(KSPGetPC(ksp, &pc));
129: PetscCall(PCFieldSplitGetSubKSP(pc, &n, &subksp));
130: PetscCall(KSPGetPC(subksp[0], &pc));
131: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
132: PetscCall(PCHPDDMSetAuxiliaryMat(pc, is[0], aux[0], NULL, NULL));
133: #endif
134: PetscCall(PCSetFromOptions(pc));
135: PetscCall(KSPGetPC(subksp[1], &pc));
136: #if defined(PETSC_HAVE_HPDDM) && defined(PETSC_HAVE_DYNAMIC_LIBRARIES) && defined(PETSC_USE_SHARED_LIBRARIES)
137: if (!flg[0]) PetscCall(PCHPDDMSetAuxiliaryMat(pc, is[1], aux[1], NULL, NULL));
138: #endif
139: PetscCall(PCSetFromOptions(pc));
140: PetscCall(PetscFree(subksp));
141: PetscCall(KSPSolve(ksp, b, x));
142: PetscCall(KSPGetConvergedReason(ksp, reason + 1));
143: PetscCall(KSPGetTotalIterations(ksp, iterations + 1));
144: iterations[1] -= iterations[0];
145: PetscCheck(reason[0] == reason[1] && PetscAbs(iterations[0] - iterations[1]) <= 3, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Successive calls to KSPSolve() did not converge for the same reason (%s v. %s) or with the same number of iterations (+/- 3, %" PetscInt_FMT " v. %" PetscInt_FMT ")", KSPConvergedReasons[reason[0]], KSPConvergedReasons[reason[1]], iterations[0], iterations[1]);
146: }
147: PetscCall(VecDestroy(&x));
148: PetscCall(VecDestroy(&b));
149: PetscCall(KSPDestroy(&ksp));
150: PetscCall(MatDestroy(&S));
151: PetscCall(MatDestroy(A + 1));
152: PetscCall(MatDestroy(A + 2));
153: for (PetscInt i = 0; i < 2; ++i) {
154: PetscCall(MatDestroy(A + (i ? 3 : 0)));
155: PetscCall(MatDestroy(aux + i));
156: PetscCall(ISDestroy(is + i));
157: }
158: PetscCall(PetscFinalize());
159: return 0;
160: }
162: PetscErrorCode MatAndISLoad(const char *prefix, const char *identifier, Mat A, IS is, Mat aux, PetscMPIInt rank, PetscMPIInt size)
163: {
164: IS sizes;
165: const PetscInt *idx;
166: PetscViewer viewer;
167: char name[PETSC_MAX_PATH_LEN];
169: PetscFunctionBeginUser;
170: PetscCall(PetscSNPrintf(name, sizeof(name), "%s%s_sizes_%d_%d.dat", prefix, identifier, rank, size));
171: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
172: PetscCall(ISCreate(PETSC_COMM_SELF, &sizes));
173: PetscCall(ISLoad(sizes, viewer));
174: PetscCall(ISGetIndices(sizes, &idx));
175: PetscCall(MatSetSizes(A, idx[0], idx[1], idx[2], idx[3]));
176: PetscCall(MatSetUp(A));
177: PetscCall(ISRestoreIndices(sizes, &idx));
178: PetscCall(ISDestroy(&sizes));
179: PetscCall(PetscViewerDestroy(&viewer));
180: PetscCall(PetscSNPrintf(name, sizeof(name), "%s%s.dat", prefix, identifier));
181: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, name, FILE_MODE_READ, &viewer));
182: PetscCall(MatLoad(A, viewer));
183: PetscCall(PetscViewerDestroy(&viewer));
184: PetscCall(PetscSNPrintf(name, sizeof(name), "%s%s_is_%d_%d.dat", prefix, identifier, rank, size));
185: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
186: PetscCall(ISLoad(is, viewer));
187: PetscCall(PetscViewerDestroy(&viewer));
188: PetscCall(PetscSNPrintf(name, sizeof(name), "%s%s_aux_%d_%d.dat", prefix, identifier, rank, size));
189: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, name, FILE_MODE_READ, &viewer));
190: PetscCall(MatLoad(aux, viewer));
191: PetscCall(PetscViewerDestroy(&viewer));
192: PetscFunctionReturn(PETSC_SUCCESS);
193: }
195: /*TEST
197: testset:
198: requires: datafilespath hpddm slepc double !complex !defined(PETSC_USE_64BIT_INDICES) defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
199: nsize: 4
200: args: -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -ksp_monitor -ksp_rtol 1e-4 -fieldsplit_ksp_max_it 100 -fieldsplit_pc_hpddm_levels_1_eps_nev 10 -fieldsplit_pc_hpddm_levels_1_st_share_sub_ksp -fieldsplit_pc_hpddm_has_neumann -fieldsplit_pc_hpddm_define_subdomains -fieldsplit_1_pc_hpddm_schur_precondition geneo -fieldsplit_pc_hpddm_coarse_pc_type redundant -fieldsplit_pc_hpddm_coarse_redundant_pc_type cholesky -fieldsplit_pc_hpddm_levels_1_sub_pc_type lu -fieldsplit_ksp_type fgmres -ksp_type fgmres -ksp_max_it 10 -fieldsplit_1_pc_hpddm_coarse_correction balanced -fieldsplit_1_pc_hpddm_levels_1_eps_gen_non_hermitian -fieldsplit_1_pc_hpddm_coarse_p 2
201: test:
202: requires: mumps
203: suffix: 1
204: args: -viewer -system {{elasticity stokes}separate output} -fieldsplit_1_pc_hpddm_ksp_pc_side left -fieldsplit_1_pc_hpddm_levels_1_sub_mat_mumps_icntl_26 1
205: filter: grep -v -e "action of " -e " " -e "block size" -e "total: nonzeros=" -e "using I-node" -e "aij" -e "transpose" -e "diagonal" -e "total number of" -e " rows="
206: test:
207: requires: mumps
208: suffix: 2
209: output_file: output/ex87_1_system-stokes.out
210: args: -viewer -system stokes -empty_A11 -transpose {{false true}shared output} -permute {{false true}shared output} -fieldsplit_1_pc_hpddm_ksp_pc_side right -fieldsplit_1_pc_hpddm_coarse_mat_type baij -fieldsplit_1_pc_hpddm_levels_1_sub_mat_mumps_icntl_26 1 -explicit {{false true}shared output}
211: filter: grep -v -e "action of " -e " " -e "block size" -e "total: nonzeros=" -e "using I-node" -e "aij" -e "transpose" -e "diagonal" -e "total number of" -e " rows=" | sed -e "s/ right preconditioning/ left preconditioning/g" -e "s/ using UNPRECONDITIONED/ using PRECONDITIONED/g"
212: test:
213: suffix: 1_petsc
214: args: -system {{elasticity stokes}separate output} -fieldsplit_1_pc_hpddm_ksp_pc_side left -fieldsplit_1_pc_hpddm_levels_1_sub_pc_factor_mat_solver_type petsc -fieldsplit_1_pc_hpddm_levels_1_eps_threshold 0.3 -permute
215: test:
216: suffix: 2_petsc
217: output_file: output/ex87_1_petsc_system-stokes.out
218: args: -system stokes -empty_A11 -transpose -fieldsplit_1_pc_hpddm_ksp_pc_side right -fieldsplit_1_pc_hpddm_levels_1_sub_pc_factor_mat_solver_type petsc -fieldsplit_1_pc_hpddm_coarse_mat_type baij -fieldsplit_1_pc_hpddm_levels_1_eps_threshold 0.3 -fieldsplit_1_pc_hpddm_levels_1_sub_pc_factor_shift_type inblocks -successive_solves
219: filter: sed -e "s/type: transpose/type: hermitiantranspose/g"
220: test:
221: suffix: threshold
222: output_file: output/ex87_1_petsc_system-elasticity.out
223: args: -fieldsplit_1_pc_hpddm_ksp_pc_side left -fieldsplit_1_pc_hpddm_levels_1_eps_threshold 0.2 -fieldsplit_1_pc_hpddm_coarse_mat_type {{baij sbaij}shared output} -successive_solves
225: test:
226: requires: datafilespath double !complex !defined(PETSC_USE_64BIT_INDICES) !hpddm !memkind
227: nsize: 4
228: suffix: selfp
229: output_file: output/ex41_1.out
230: filter: grep -v "CONVERGED_RTOL iterations"
231: args: -load_dir ${DATAFILESPATH}/matrices/hpddm/GENEO -system stokes -ksp_rtol 1e-4 -ksp_converged_reason -ksp_max_it 30 -pc_type fieldsplit -pc_fieldsplit_type schur -fieldsplit_ksp_type preonly -pc_fieldsplit_schur_precondition selfp -fieldsplit_pc_type bjacobi -fieldsplit_sub_pc_type lu -transpose {{false true}shared output} -fieldsplit_1_mat_schur_complement_ainv_type lump
233: TEST*/