Actual source code: bddcfetidp.c
1: #include <petsc/private/pcbddcimpl.h>
2: #include <petsc/private/pcbddcprivateimpl.h>
3: #include <petscblaslapack.h>
5: static PetscErrorCode MatMult_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
6: {
7: BDdelta_DN ctx;
9: PetscFunctionBegin;
10: PetscCall(MatShellGetContext(A, &ctx));
11: PetscCall(MatMultTranspose(ctx->BD, x, ctx->work));
12: PetscCall(KSPSolveTranspose(ctx->kBD, ctx->work, y));
13: /* No PC so cannot propagate up failure in KSPSolveTranspose() */
14: PetscFunctionReturn(PETSC_SUCCESS);
15: }
17: static PetscErrorCode MatMultTranspose_BDdelta_deluxe_nonred(Mat A, Vec x, Vec y)
18: {
19: BDdelta_DN ctx;
21: PetscFunctionBegin;
22: PetscCall(MatShellGetContext(A, &ctx));
23: PetscCall(KSPSolve(ctx->kBD, x, ctx->work));
24: /* No PC so cannot propagate up failure in KSPSolve() */
25: PetscCall(MatMult(ctx->BD, ctx->work, y));
26: PetscFunctionReturn(PETSC_SUCCESS);
27: }
29: static PetscErrorCode MatDestroy_BDdelta_deluxe_nonred(Mat A)
30: {
31: BDdelta_DN ctx;
33: PetscFunctionBegin;
34: PetscCall(MatShellGetContext(A, &ctx));
35: PetscCall(MatDestroy(&ctx->BD));
36: PetscCall(KSPDestroy(&ctx->kBD));
37: PetscCall(VecDestroy(&ctx->work));
38: PetscCall(PetscFree(ctx));
39: PetscFunctionReturn(PETSC_SUCCESS);
40: }
42: PetscErrorCode PCBDDCCreateFETIDPMatContext(PC pc, FETIDPMat_ctx *fetidpmat_ctx)
43: {
44: FETIDPMat_ctx newctx;
46: PetscFunctionBegin;
47: PetscCall(PetscNew(&newctx));
48: /* increase the reference count for BDDC preconditioner */
49: PetscCall(PetscObjectReference((PetscObject)pc));
50: newctx->pc = pc;
51: *fetidpmat_ctx = newctx;
52: PetscFunctionReturn(PETSC_SUCCESS);
53: }
55: PetscErrorCode PCBDDCCreateFETIDPPCContext(PC pc, FETIDPPC_ctx *fetidppc_ctx)
56: {
57: FETIDPPC_ctx newctx;
59: PetscFunctionBegin;
60: PetscCall(PetscNew(&newctx));
61: /* increase the reference count for BDDC preconditioner */
62: PetscCall(PetscObjectReference((PetscObject)pc));
63: newctx->pc = pc;
64: *fetidppc_ctx = newctx;
65: PetscFunctionReturn(PETSC_SUCCESS);
66: }
68: PetscErrorCode PCBDDCDestroyFETIDPMat(Mat A)
69: {
70: FETIDPMat_ctx mat_ctx;
72: PetscFunctionBegin;
73: PetscCall(MatShellGetContext(A, &mat_ctx));
74: PetscCall(VecDestroy(&mat_ctx->lambda_local));
75: PetscCall(VecDestroy(&mat_ctx->temp_solution_D));
76: PetscCall(VecDestroy(&mat_ctx->temp_solution_B));
77: PetscCall(MatDestroy(&mat_ctx->B_delta));
78: PetscCall(MatDestroy(&mat_ctx->B_Ddelta));
79: PetscCall(MatDestroy(&mat_ctx->B_BB));
80: PetscCall(MatDestroy(&mat_ctx->B_BI));
81: PetscCall(MatDestroy(&mat_ctx->Bt_BB));
82: PetscCall(MatDestroy(&mat_ctx->Bt_BI));
83: PetscCall(MatDestroy(&mat_ctx->C));
84: PetscCall(VecDestroy(&mat_ctx->rhs_flip));
85: PetscCall(VecDestroy(&mat_ctx->vP));
86: PetscCall(VecDestroy(&mat_ctx->xPg));
87: PetscCall(VecDestroy(&mat_ctx->yPg));
88: PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda));
89: PetscCall(VecScatterDestroy(&mat_ctx->l2g_lambda_only));
90: PetscCall(VecScatterDestroy(&mat_ctx->l2g_p));
91: PetscCall(VecScatterDestroy(&mat_ctx->g2g_p));
92: PetscCall(PCDestroy(&mat_ctx->pc)); /* decrease PCBDDC reference count */
93: PetscCall(ISDestroy(&mat_ctx->pressure));
94: PetscCall(ISDestroy(&mat_ctx->lagrange));
95: PetscCall(PetscFree(mat_ctx));
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: PetscErrorCode PCBDDCDestroyFETIDPPC(PC pc)
100: {
101: FETIDPPC_ctx pc_ctx;
103: PetscFunctionBegin;
104: PetscCall(PCShellGetContext(pc, &pc_ctx));
105: PetscCall(VecDestroy(&pc_ctx->lambda_local));
106: PetscCall(MatDestroy(&pc_ctx->B_Ddelta));
107: PetscCall(VecScatterDestroy(&pc_ctx->l2g_lambda));
108: PetscCall(MatDestroy(&pc_ctx->S_j));
109: PetscCall(PCDestroy(&pc_ctx->pc)); /* decrease PCBDDC reference count */
110: PetscCall(VecDestroy(&pc_ctx->xPg));
111: PetscCall(VecDestroy(&pc_ctx->yPg));
112: PetscCall(PetscFree(pc_ctx));
113: PetscFunctionReturn(PETSC_SUCCESS);
114: }
116: PetscErrorCode PCBDDCSetupFETIDPMatContext(FETIDPMat_ctx fetidpmat_ctx)
117: {
118: PC_IS *pcis = (PC_IS *)fetidpmat_ctx->pc->data;
119: PC_BDDC *pcbddc = (PC_BDDC *)fetidpmat_ctx->pc->data;
120: PCBDDCGraph mat_graph = pcbddc->mat_graph;
121: Mat_IS *matis = (Mat_IS *)fetidpmat_ctx->pc->pmat->data;
122: MPI_Comm comm;
123: Mat ScalingMat, BD1, BD2;
124: Vec fetidp_global;
125: IS IS_l2g_lambda;
126: IS subset, subset_mult, subset_n, isvert;
127: PetscBool skip_node, fully_redundant;
128: PetscInt i, j, k, s, n_boundary_dofs, n_global_lambda, n_vertices, partial_sum;
129: PetscInt cum, n_local_lambda, n_lambda_for_dof, dual_size, n_neg_values, n_pos_values;
130: PetscMPIInt rank, size, buf_size, neigh;
131: PetscScalar scalar_value;
132: const PetscInt *vertex_indices;
133: PetscInt *dual_dofs_boundary_indices, *aux_local_numbering_1;
134: const PetscInt *aux_global_numbering;
135: PetscInt *aux_sums, *cols_B_delta, *l2g_indices;
136: PetscScalar *array, *scaling_factors, *vals_B_delta;
137: PetscScalar **all_factors;
138: PetscInt *aux_local_numbering_2;
139: PetscInt *count, **neighbours_set;
140: PetscLayout llay;
142: /* saddlepoint */
143: ISLocalToGlobalMapping l2gmap_p;
144: PetscLayout play;
145: IS gP, pP;
146: PetscInt nPl, nPg, nPgl;
148: PetscFunctionBegin;
149: PetscCall(PetscObjectGetComm((PetscObject)fetidpmat_ctx->pc, &comm));
150: PetscCallMPI(MPI_Comm_rank(comm, &rank));
151: PetscCallMPI(MPI_Comm_size(comm, &size));
153: /* saddlepoint */
154: nPl = 0;
155: nPg = 0;
156: nPgl = 0;
157: gP = NULL;
158: pP = NULL;
159: l2gmap_p = NULL;
160: play = NULL;
161: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_pP", (PetscObject *)&pP));
162: if (pP) { /* saddle point */
163: /* subdomain pressures in global numbering */
164: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_gP", (PetscObject *)&gP));
165: PetscCheck(gP, PETSC_COMM_SELF, PETSC_ERR_PLIB, "gP not present");
166: PetscCall(ISGetLocalSize(gP, &nPl));
167: PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->vP));
168: PetscCall(VecSetSizes(fetidpmat_ctx->vP, nPl, nPl));
169: PetscCall(VecSetType(fetidpmat_ctx->vP, VECSTANDARD));
170: PetscCall(VecSetUp(fetidpmat_ctx->vP));
172: /* pressure matrix */
173: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_C", (PetscObject *)&fetidpmat_ctx->C));
174: if (!fetidpmat_ctx->C) { /* null pressure block, compute layout and global numbering for pressures */
175: IS Pg;
177: PetscCall(ISRenumber(gP, NULL, &nPg, &Pg));
178: PetscCall(ISLocalToGlobalMappingCreateIS(Pg, &l2gmap_p));
179: PetscCall(ISDestroy(&Pg));
180: PetscCall(PetscLayoutCreate(comm, &play));
181: PetscCall(PetscLayoutSetBlockSize(play, 1));
182: PetscCall(PetscLayoutSetSize(play, nPg));
183: PetscCall(ISGetLocalSize(pP, &nPgl));
184: PetscCall(PetscLayoutSetLocalSize(play, nPgl));
185: PetscCall(PetscLayoutSetUp(play));
186: } else {
187: PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->C));
188: PetscCall(MatISGetLocalToGlobalMapping(fetidpmat_ctx->C, &l2gmap_p, NULL));
189: PetscCall(PetscObjectReference((PetscObject)l2gmap_p));
190: PetscCall(MatGetSize(fetidpmat_ctx->C, &nPg, NULL));
191: PetscCall(MatGetLocalSize(fetidpmat_ctx->C, NULL, &nPgl));
192: PetscCall(MatGetLayouts(fetidpmat_ctx->C, NULL, &llay));
193: PetscCall(PetscLayoutReference(llay, &play));
194: }
195: PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->xPg));
196: PetscCall(VecCreateMPIWithArray(comm, 1, nPgl, nPg, NULL, &fetidpmat_ctx->yPg));
198: /* import matrices for pressures coupling */
199: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BI", (PetscObject *)&fetidpmat_ctx->B_BI));
200: PetscCheck(fetidpmat_ctx->B_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BI not present");
201: PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BI));
203: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_B_BB", (PetscObject *)&fetidpmat_ctx->B_BB));
204: PetscCheck(fetidpmat_ctx->B_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B_BB not present");
205: PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->B_BB));
207: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BI", (PetscObject *)&fetidpmat_ctx->Bt_BI));
208: PetscCheck(fetidpmat_ctx->Bt_BI, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BI not present");
209: PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BI));
211: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_Bt_BB", (PetscObject *)&fetidpmat_ctx->Bt_BB));
212: PetscCheck(fetidpmat_ctx->Bt_BB, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bt_BB not present");
213: PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->Bt_BB));
215: PetscCall(PetscObjectQuery((PetscObject)fetidpmat_ctx->pc, "__KSPFETIDP_flip", (PetscObject *)&fetidpmat_ctx->rhs_flip));
216: if (fetidpmat_ctx->rhs_flip) PetscCall(PetscObjectReference((PetscObject)fetidpmat_ctx->rhs_flip));
217: }
219: /* Default type of lagrange multipliers is non-redundant */
220: fully_redundant = fetidpmat_ctx->fully_redundant;
222: /* Evaluate local and global number of lagrange multipliers */
223: PetscCall(VecSet(pcis->vec1_N, 0.0));
224: n_local_lambda = 0;
225: partial_sum = 0;
226: n_boundary_dofs = 0;
228: /* Get Vertices used to define the BDDC */
229: PetscCall(PCBDDCGraphGetCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert));
230: PetscCall(ISGetLocalSize(isvert, &n_vertices));
231: PetscCall(ISGetIndices(isvert, &vertex_indices));
233: dual_size = pcis->n_B - n_vertices;
234: PetscCall(PetscMalloc1(dual_size, &dual_dofs_boundary_indices));
235: PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_1));
236: PetscCall(PetscMalloc1(dual_size, &aux_local_numbering_2));
238: /* the code below does not support multiple subdomains per process
239: error out in this case
240: TODO: I guess I can use PetscSFGetMultiSF and the code will be easier and more general */
241: PetscCall(PetscMalloc2(pcis->n, &count, pcis->n, &neighbours_set));
242: for (i = 0, j = 0; i < pcis->n; i++) j += mat_graph->nodes[i].count;
243: if (pcis->n) PetscCall(PetscMalloc1(j, &neighbours_set[0]));
244: for (i = 0; i < pcis->n; i++) {
245: PCBDDCGraphNode *node = &mat_graph->nodes[i];
247: count[i] = 0;
248: for (j = 0; j < node->count; j++) {
249: if (node->neighbours_set[j] == rank) continue;
250: neighbours_set[i][count[i]++] = node->neighbours_set[j];
251: }
252: PetscCheck(count[i] == node->count - 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported");
253: s = count[i];
254: PetscCall(PetscSortRemoveDupsInt(count + i, neighbours_set[i]));
255: PetscCheck(s == count[i], PETSC_COMM_SELF, PETSC_ERR_SUP, "Multiple subdomains per process not supported");
256: if (i != pcis->n - 1) neighbours_set[i + 1] = neighbours_set[i] + count[i];
257: }
259: PetscCall(VecGetArray(pcis->vec1_N, &array));
260: for (i = 0, s = 0; i < pcis->n; i++) {
261: j = count[i]; /* RECALL: count[i] does not count myself */
262: if (j > 0) n_boundary_dofs++;
263: skip_node = PETSC_FALSE;
264: if (s < n_vertices && vertex_indices[s] == i) { /* it works for a sorted set of vertices */
265: skip_node = PETSC_TRUE;
266: s++;
267: }
268: if (j < 1) skip_node = PETSC_TRUE;
269: if (mat_graph->nodes[i].special_dof == PCBDDCGRAPH_DIRICHLET_MARK) skip_node = PETSC_TRUE;
270: if (!skip_node) {
271: if (fully_redundant) {
272: /* fully redundant set of lagrange multipliers */
273: n_lambda_for_dof = (j * (j + 1)) / 2;
274: } else {
275: n_lambda_for_dof = j;
276: }
277: n_local_lambda += j;
278: /* needed to evaluate global number of lagrange multipliers */
279: array[i] = (1.0 * n_lambda_for_dof) / (j + 1.0); /* already scaled for the next global sum */
280: /* store some data needed */
281: dual_dofs_boundary_indices[partial_sum] = n_boundary_dofs - 1;
282: aux_local_numbering_1[partial_sum] = i;
283: aux_local_numbering_2[partial_sum] = n_lambda_for_dof;
284: partial_sum++;
285: }
286: }
287: PetscCall(VecRestoreArray(pcis->vec1_N, &array));
288: PetscCall(ISRestoreIndices(isvert, &vertex_indices));
289: PetscCall(PCBDDCGraphRestoreCandidatesIS(mat_graph, NULL, NULL, NULL, NULL, &isvert));
290: dual_size = partial_sum;
292: /* compute global ordering of lagrange multipliers and associate l2g map */
293: PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_1, PETSC_COPY_VALUES, &subset_n));
294: PetscCall(ISLocalToGlobalMappingApplyIS(pcis->mapping, subset_n, &subset));
295: PetscCall(ISDestroy(&subset_n));
296: PetscCall(ISCreateGeneral(comm, partial_sum, aux_local_numbering_2, PETSC_OWN_POINTER, &subset_mult));
297: PetscCall(ISRenumber(subset, subset_mult, &fetidpmat_ctx->n_lambda, &subset_n));
298: PetscCall(ISDestroy(&subset));
300: if (PetscDefined(USE_DEBUG)) {
301: PetscCall(VecSet(pcis->vec1_global, 0.0));
302: PetscCall(VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
303: PetscCall(VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE));
304: PetscCall(VecSum(pcis->vec1_global, &scalar_value));
305: i = (PetscInt)PetscRealPart(scalar_value);
306: PetscCheck(i == fetidpmat_ctx->n_lambda, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Global number of multipliers mismatch! (%" PetscInt_FMT " != %" PetscInt_FMT ")", fetidpmat_ctx->n_lambda, i);
307: }
309: /* init data for scaling factors exchange */
310: if (!pcbddc->use_deluxe_scaling) {
311: PetscInt *ptrs_buffer, neigh_position;
312: PetscScalar *send_buffer, *recv_buffer;
313: MPI_Request *send_reqs, *recv_reqs;
315: partial_sum = 0;
316: PetscCall(PetscMalloc1(pcis->n_neigh, &ptrs_buffer));
317: PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &send_reqs));
318: PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &recv_reqs));
319: PetscCall(PetscMalloc1(pcis->n + 1, &all_factors));
320: if (pcis->n_neigh > 0) ptrs_buffer[0] = 0;
321: for (i = 1; i < pcis->n_neigh; i++) {
322: partial_sum += pcis->n_shared[i];
323: ptrs_buffer[i] = ptrs_buffer[i - 1] + pcis->n_shared[i];
324: }
325: PetscCall(PetscMalloc1(partial_sum, &send_buffer));
326: PetscCall(PetscMalloc1(partial_sum, &recv_buffer));
327: PetscCall(PetscMalloc1(partial_sum, &all_factors[0]));
328: for (i = 0; i < pcis->n - 1; i++) {
329: j = count[i];
330: all_factors[i + 1] = all_factors[i] + j;
331: }
333: /* scatter B scaling to N vec */
334: PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
335: PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
336: /* communications */
337: PetscCall(VecGetArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
338: for (i = 1; i < pcis->n_neigh; i++) {
339: for (j = 0; j < pcis->n_shared[i]; j++) send_buffer[ptrs_buffer[i - 1] + j] = array[pcis->shared[i][j]];
340: PetscCall(PetscMPIIntCast(ptrs_buffer[i] - ptrs_buffer[i - 1], &buf_size));
341: PetscCall(PetscMPIIntCast(pcis->neigh[i], &neigh));
342: PetscCallMPI(MPI_Isend(&send_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &send_reqs[i - 1]));
343: PetscCallMPI(MPI_Irecv(&recv_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &recv_reqs[i - 1]));
344: }
345: PetscCall(VecRestoreArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
346: if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(pcis->n_neigh - 1, recv_reqs, MPI_STATUSES_IGNORE));
347: /* put values in correct places */
348: for (i = 1; i < pcis->n_neigh; i++) {
349: for (j = 0; j < pcis->n_shared[i]; j++) {
350: k = pcis->shared[i][j];
351: neigh_position = 0;
352: while (neighbours_set[k][neigh_position] != pcis->neigh[i]) { neigh_position++; }
353: all_factors[k][neigh_position] = recv_buffer[ptrs_buffer[i - 1] + j];
354: }
355: }
356: if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(pcis->n_neigh - 1, send_reqs, MPI_STATUSES_IGNORE));
357: PetscCall(PetscFree(send_reqs));
358: PetscCall(PetscFree(recv_reqs));
359: PetscCall(PetscFree(send_buffer));
360: PetscCall(PetscFree(recv_buffer));
361: PetscCall(PetscFree(ptrs_buffer));
362: }
364: /* Compute B and B_delta (local actions) */
365: PetscCall(PetscMalloc1(pcis->n_neigh, &aux_sums));
366: PetscCall(PetscMalloc1(n_local_lambda, &l2g_indices));
367: PetscCall(PetscMalloc1(n_local_lambda, &vals_B_delta));
368: PetscCall(PetscMalloc1(n_local_lambda, &cols_B_delta));
369: if (!pcbddc->use_deluxe_scaling) {
370: PetscCall(PetscMalloc1(n_local_lambda, &scaling_factors));
371: } else {
372: scaling_factors = NULL;
373: all_factors = NULL;
374: }
375: PetscCall(ISGetIndices(subset_n, &aux_global_numbering));
376: partial_sum = 0;
377: cum = 0;
378: for (i = 0; i < dual_size; i++) {
379: n_global_lambda = aux_global_numbering[cum];
380: j = count[aux_local_numbering_1[i]];
381: aux_sums[0] = 0;
382: for (s = 1; s < j; s++) aux_sums[s] = aux_sums[s - 1] + j - s + 1;
383: if (all_factors) array = all_factors[aux_local_numbering_1[i]];
384: n_neg_values = 0;
385: while (n_neg_values < j && neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) { n_neg_values++; }
386: n_pos_values = j - n_neg_values;
387: if (fully_redundant) {
388: for (s = 0; s < n_neg_values; s++) {
389: l2g_indices[partial_sum + s] = aux_sums[s] + n_neg_values - s - 1 + n_global_lambda;
390: cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
391: vals_B_delta[partial_sum + s] = -1.0;
392: if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s] = array[s];
393: }
394: for (s = 0; s < n_pos_values; s++) {
395: l2g_indices[partial_sum + s + n_neg_values] = aux_sums[n_neg_values] + s + n_global_lambda;
396: cols_B_delta[partial_sum + s + n_neg_values] = dual_dofs_boundary_indices[i];
397: vals_B_delta[partial_sum + s + n_neg_values] = 1.0;
398: if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s + n_neg_values] = array[s + n_neg_values];
399: }
400: partial_sum += j;
401: } else {
402: /* l2g_indices and default cols and vals of B_delta */
403: for (s = 0; s < j; s++) {
404: l2g_indices[partial_sum + s] = n_global_lambda + s;
405: cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
406: vals_B_delta[partial_sum + s] = 0.0;
407: }
408: /* B_delta */
409: if (n_neg_values > 0) { /* there's a rank next to me to the left */
410: vals_B_delta[partial_sum + n_neg_values - 1] = -1.0;
411: }
412: if (n_neg_values < j) { /* there's a rank next to me to the right */
413: vals_B_delta[partial_sum + n_neg_values] = 1.0;
414: }
415: /* scaling as in Klawonn-Widlund 1999 */
416: if (!pcbddc->use_deluxe_scaling) {
417: for (s = 0; s < n_neg_values; s++) {
418: scalar_value = 0.0;
419: for (k = 0; k < s + 1; k++) scalar_value += array[k];
420: scaling_factors[partial_sum + s] = -scalar_value;
421: }
422: for (s = 0; s < n_pos_values; s++) {
423: scalar_value = 0.0;
424: for (k = s + n_neg_values; k < j; k++) scalar_value += array[k];
425: scaling_factors[partial_sum + s + n_neg_values] = scalar_value;
426: }
427: }
428: partial_sum += j;
429: }
430: cum += aux_local_numbering_2[i];
431: }
432: PetscCall(ISRestoreIndices(subset_n, &aux_global_numbering));
433: PetscCall(ISDestroy(&subset_mult));
434: PetscCall(ISDestroy(&subset_n));
435: PetscCall(PetscFree(aux_sums));
436: PetscCall(PetscFree(aux_local_numbering_1));
437: PetscCall(PetscFree(dual_dofs_boundary_indices));
438: if (all_factors) {
439: PetscCall(PetscFree(all_factors[0]));
440: PetscCall(PetscFree(all_factors));
441: }
442: if (pcis->n) PetscCall(PetscFree(neighbours_set[0]));
443: PetscCall(PetscFree2(count, neighbours_set));
445: /* Create local part of B_delta */
446: PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_delta));
447: PetscCall(MatSetSizes(fetidpmat_ctx->B_delta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
448: PetscCall(MatSetType(fetidpmat_ctx->B_delta, MATSEQAIJ));
449: PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta, 1, NULL));
450: PetscCall(MatSetOption(fetidpmat_ctx->B_delta, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
451: for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_delta, i, cols_B_delta[i], vals_B_delta[i], INSERT_VALUES));
452: PetscCall(PetscFree(vals_B_delta));
453: PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
454: PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
456: BD1 = NULL;
457: BD2 = NULL;
458: if (fully_redundant) {
459: PetscCheck(!pcbddc->use_deluxe_scaling, comm, PETSC_ERR_SUP, "Deluxe FETIDP with fully-redundant multipliers to be implemented");
460: PetscCall(MatCreate(PETSC_COMM_SELF, &ScalingMat));
461: PetscCall(MatSetSizes(ScalingMat, n_local_lambda, n_local_lambda, n_local_lambda, n_local_lambda));
462: PetscCall(MatSetType(ScalingMat, MATSEQAIJ));
463: PetscCall(MatSeqAIJSetPreallocation(ScalingMat, 1, NULL));
464: for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(ScalingMat, i, i, scaling_factors[i], INSERT_VALUES));
465: PetscCall(MatAssemblyBegin(ScalingMat, MAT_FINAL_ASSEMBLY));
466: PetscCall(MatAssemblyEnd(ScalingMat, MAT_FINAL_ASSEMBLY));
467: PetscCall(MatMatMult(ScalingMat, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &fetidpmat_ctx->B_Ddelta));
468: PetscCall(MatDestroy(&ScalingMat));
469: } else {
470: PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_Ddelta));
471: PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
472: if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
473: PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSEQAIJ));
474: PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta, 1, NULL));
475: for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta, i, cols_B_delta[i], scaling_factors[i], INSERT_VALUES));
476: PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
477: PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
478: } else {
479: /* scaling as in Klawonn-Widlund 1999 */
480: PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
481: PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
482: Mat T;
483: PetscScalar *W, lwork, *Bwork;
484: const PetscInt *idxs = NULL;
485: PetscInt cum, mss, *nnz;
486: PetscBLASInt *pivots, B_lwork, B_N, B_ierr;
488: PetscCheck(pcbddc->deluxe_singlemat, comm, PETSC_ERR_USER, "Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
489: mss = 0;
490: PetscCall(PetscCalloc1(pcis->n_B, &nnz));
491: if (sub_schurs->is_Ej_all) {
492: PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs));
493: for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
494: PetscInt subset_size;
496: PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
497: for (j = 0; j < subset_size; j++) nnz[idxs[j + cum]] = subset_size;
498: mss = PetscMax(mss, subset_size);
499: cum += subset_size;
500: }
501: }
502: PetscCall(MatCreate(PETSC_COMM_SELF, &T));
503: PetscCall(MatSetSizes(T, pcis->n_B, pcis->n_B, pcis->n_B, pcis->n_B));
504: PetscCall(MatSetType(T, MATSEQAIJ));
505: PetscCall(MatSeqAIJSetPreallocation(T, 0, nnz));
506: PetscCall(PetscFree(nnz));
508: /* workspace allocation */
509: B_lwork = 0;
510: if (mss) {
511: PetscScalar dummy = 1;
513: B_lwork = -1;
514: PetscCall(PetscBLASIntCast(mss, &B_N));
515: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
516: PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummy, &B_N, &B_N, &lwork, &B_lwork, &B_ierr));
517: PetscCall(PetscFPTrapPop());
518: PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %d", (int)B_ierr);
519: PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork));
520: }
521: PetscCall(PetscMalloc3(mss * mss, &W, mss, &pivots, B_lwork, &Bwork));
523: for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
524: const PetscScalar *M;
525: PetscInt subset_size;
527: PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
528: PetscCall(PetscBLASIntCast(subset_size, &B_N));
529: PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i], &M));
530: PetscCall(PetscArraycpy(W, M, subset_size * subset_size));
531: PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i], &M));
532: PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
533: PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, W, &B_N, pivots, &B_ierr));
534: PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %d", (int)B_ierr);
535: PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, W, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
536: PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %d", (int)B_ierr);
537: PetscCall(PetscFPTrapPop());
538: /* silent static analyzer */
539: PetscCheck(idxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IDXS not present");
540: PetscCall(MatSetValues(T, subset_size, idxs + cum, subset_size, idxs + cum, W, INSERT_VALUES));
541: cum += subset_size;
542: }
543: PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
544: PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
545: PetscCall(MatMatTransposeMult(T, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &BD1));
546: PetscCall(MatMatMult(fetidpmat_ctx->B_delta, BD1, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &BD2));
547: PetscCall(MatDestroy(&T));
548: PetscCall(PetscFree3(W, pivots, Bwork));
549: if (sub_schurs->is_Ej_all) PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs));
550: }
551: }
552: PetscCall(PetscFree(scaling_factors));
553: PetscCall(PetscFree(cols_B_delta));
555: /* Layout of multipliers */
556: PetscCall(PetscLayoutCreate(comm, &llay));
557: PetscCall(PetscLayoutSetBlockSize(llay, 1));
558: PetscCall(PetscLayoutSetSize(llay, fetidpmat_ctx->n_lambda));
559: PetscCall(PetscLayoutSetUp(llay));
560: PetscCall(PetscLayoutGetLocalSize(llay, &fetidpmat_ctx->n));
562: /* Local work vector of multipliers */
563: PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->lambda_local));
564: PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local, n_local_lambda, n_local_lambda));
565: PetscCall(VecSetType(fetidpmat_ctx->lambda_local, VECSEQ));
567: if (BD2) {
568: ISLocalToGlobalMapping l2g;
569: Mat T, TA, *pT;
570: IS is;
571: PetscInt nl, N;
572: BDdelta_DN ctx;
574: PetscCall(PetscLayoutGetLocalSize(llay, &nl));
575: PetscCall(PetscLayoutGetSize(llay, &N));
576: PetscCall(MatCreate(comm, &T));
577: PetscCall(MatSetSizes(T, nl, nl, N, N));
578: PetscCall(MatSetType(T, MATIS));
579: PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &l2g));
580: PetscCall(MatSetLocalToGlobalMapping(T, l2g, l2g));
581: PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
582: PetscCall(MatISSetLocalMat(T, BD2));
583: PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
584: PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
585: PetscCall(MatDestroy(&BD2));
586: PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &TA));
587: PetscCall(MatDestroy(&T));
588: PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_USE_POINTER, &is));
589: PetscCall(MatCreateSubMatrices(TA, 1, &is, &is, MAT_INITIAL_MATRIX, &pT));
590: PetscCall(MatDestroy(&TA));
591: PetscCall(ISDestroy(&is));
592: BD2 = pT[0];
593: PetscCall(PetscFree(pT));
595: /* B_Ddelta for non-redundant multipliers with deluxe scaling */
596: PetscCall(PetscNew(&ctx));
597: PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSHELL));
598: PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta, ctx));
599: PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT, (void (*)(void))MatMult_BDdelta_deluxe_nonred));
600: PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred));
601: PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_DESTROY, (void (*)(void))MatDestroy_BDdelta_deluxe_nonred));
602: PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta));
604: PetscCall(PetscObjectReference((PetscObject)BD1));
605: ctx->BD = BD1;
606: PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->kBD));
607: PetscCall(KSPSetNestLevel(ctx->kBD, fetidpmat_ctx->pc->kspnestlevel));
608: PetscCall(KSPSetOperators(ctx->kBD, BD2, BD2));
609: PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &ctx->work));
610: fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
611: }
612: PetscCall(MatDestroy(&BD1));
613: PetscCall(MatDestroy(&BD2));
615: /* fetidpmat sizes */
616: fetidpmat_ctx->n += nPgl;
617: fetidpmat_ctx->N = fetidpmat_ctx->n_lambda + nPg;
619: /* Global vector for FETI-DP linear system */
620: PetscCall(VecCreate(comm, &fetidp_global));
621: PetscCall(VecSetSizes(fetidp_global, fetidpmat_ctx->n, fetidpmat_ctx->N));
622: PetscCall(VecSetType(fetidp_global, VECMPI));
623: PetscCall(VecSetUp(fetidp_global));
625: /* Decide layout for fetidp dofs: if it is a saddle point problem
626: pressure is ordered first in the local part of the global vector
627: of the FETI-DP linear system */
628: if (nPg) {
629: Vec v;
630: IS IS_l2g_p, ais;
631: PetscLayout alay;
632: const PetscInt *idxs, *pranges, *aranges, *lranges;
633: PetscInt *l2g_indices_p, rst;
634: PetscMPIInt rank;
636: PetscCall(PetscMalloc1(nPl, &l2g_indices_p));
637: PetscCall(VecGetLayout(fetidp_global, &alay));
638: PetscCall(PetscLayoutGetRanges(alay, &aranges));
639: PetscCall(PetscLayoutGetRanges(play, &pranges));
640: PetscCall(PetscLayoutGetRanges(llay, &lranges));
642: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global), &rank));
643: PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), pranges[rank + 1] - pranges[rank], aranges[rank], 1, &fetidpmat_ctx->pressure));
644: PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure, "F_P"));
645: PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), lranges[rank + 1] - lranges[rank], aranges[rank] + pranges[rank + 1] - pranges[rank], 1, &fetidpmat_ctx->lagrange));
646: PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange, "F_L"));
647: PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p, &idxs));
648: /* shift local to global indices for pressure */
649: for (i = 0; i < nPl; i++) {
650: PetscMPIInt owner;
652: PetscCall(PetscLayoutFindOwner(play, idxs[i], &owner));
653: l2g_indices_p[i] = idxs[i] - pranges[owner] + aranges[owner];
654: }
655: PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p, &idxs));
656: PetscCall(ISCreateGeneral(comm, nPl, l2g_indices_p, PETSC_OWN_POINTER, &IS_l2g_p));
658: /* local to global scatter for pressure */
659: PetscCall(VecScatterCreate(fetidpmat_ctx->vP, NULL, fetidp_global, IS_l2g_p, &fetidpmat_ctx->l2g_p));
660: PetscCall(ISDestroy(&IS_l2g_p));
662: /* scatter for lagrange multipliers only */
663: PetscCall(VecCreate(comm, &v));
664: PetscCall(VecSetType(v, VECSTANDARD));
665: PetscCall(VecSetLayout(v, llay));
666: PetscCall(VecSetUp(v));
667: PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &ais));
668: PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, v, ais, &fetidpmat_ctx->l2g_lambda_only));
669: PetscCall(ISDestroy(&ais));
670: PetscCall(VecDestroy(&v));
672: /* shift local to global indices for multipliers */
673: for (i = 0; i < n_local_lambda; i++) {
674: PetscInt ps;
675: PetscMPIInt owner;
677: PetscCall(PetscLayoutFindOwner(llay, l2g_indices[i], &owner));
678: ps = pranges[owner + 1] - pranges[owner];
679: l2g_indices[i] = l2g_indices[i] - lranges[owner] + aranges[owner] + ps;
680: }
682: /* scatter from alldofs to pressures global fetidp vector */
683: PetscCall(PetscLayoutGetRange(alay, &rst, NULL));
684: PetscCall(ISCreateStride(comm, nPgl, rst, 1, &ais));
685: PetscCall(VecScatterCreate(pcis->vec1_global, pP, fetidp_global, ais, &fetidpmat_ctx->g2g_p));
686: PetscCall(ISDestroy(&ais));
687: }
688: PetscCall(PetscLayoutDestroy(&llay));
689: PetscCall(PetscLayoutDestroy(&play));
690: PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_OWN_POINTER, &IS_l2g_lambda));
692: /* scatter from local to global multipliers */
693: PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, fetidp_global, IS_l2g_lambda, &fetidpmat_ctx->l2g_lambda));
694: PetscCall(ISDestroy(&IS_l2g_lambda));
695: PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p));
696: PetscCall(VecDestroy(&fetidp_global));
698: /* Create some work vectors needed by fetidp */
699: PetscCall(VecDuplicate(pcis->vec1_B, &fetidpmat_ctx->temp_solution_B));
700: PetscCall(VecDuplicate(pcis->vec1_D, &fetidpmat_ctx->temp_solution_D));
701: PetscFunctionReturn(PETSC_SUCCESS);
702: }
704: PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
705: {
706: FETIDPMat_ctx mat_ctx;
707: PC_BDDC *pcbddc = (PC_BDDC *)fetidppc_ctx->pc->data;
708: PC_IS *pcis = (PC_IS *)fetidppc_ctx->pc->data;
709: PetscBool lumped = PETSC_FALSE;
711: PetscFunctionBegin;
712: PetscCall(MatShellGetContext(fetimat, &mat_ctx));
713: /* get references from objects created when setting up feti mat context */
714: PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local));
715: fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
716: PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta));
717: fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
718: if (mat_ctx->deluxe_nonred) {
719: PC pc, mpc;
720: BDdelta_DN ctx;
721: MatSolverType solver;
722: const char *prefix;
724: PetscCall(MatShellGetContext(mat_ctx->B_Ddelta, &ctx));
725: PetscCall(KSPSetType(ctx->kBD, KSPPREONLY));
726: PetscCall(KSPGetPC(ctx->kBD, &mpc));
727: PetscCall(KSPGetPC(pcbddc->ksp_D, &pc));
728: PetscCall(PCSetType(mpc, PCLU));
729: PetscCall(PCFactorGetMatSolverType(pc, (MatSolverType *)&solver));
730: if (solver) PetscCall(PCFactorSetMatSolverType(mpc, solver));
731: PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
732: PetscCall(KSPSetOptionsPrefix(ctx->kBD, prefix));
733: PetscCall(KSPAppendOptionsPrefix(ctx->kBD, "bddelta_"));
734: PetscCall(KSPSetFromOptions(ctx->kBD));
735: }
737: if (mat_ctx->l2g_lambda_only) {
738: PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only));
739: fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
740: } else {
741: PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda));
742: fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
743: }
744: /* Dirichlet preconditioner */
745: PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_lumped", &lumped, NULL));
746: if (!lumped) {
747: IS iV;
748: PetscBool discrete_harmonic = PETSC_FALSE;
750: PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc, "__KSPFETIDP_iV", (PetscObject *)&iV));
751: if (iV) PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_discrete_harmonic", &discrete_harmonic, NULL));
752: if (discrete_harmonic) {
753: KSP sksp;
754: PC pc;
755: PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
756: Mat A_II, A_IB, A_BI;
757: IS iP = NULL;
758: PetscBool isshell, reuse = PETSC_FALSE;
759: KSPType ksptype;
760: const char *prefix;
762: /*
763: We constructs a Schur complement for
765: | A_II A_ID |
766: | A_DI A_DD |
768: instead of
770: | A_II B^t_II A_ID |
771: | B_II -C_II B_ID |
772: | A_DI B^t_ID A_DD |
774: */
775: if (sub_schurs && sub_schurs->reuse_solver) {
776: PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A, "__KSPFETIDP_iP", (PetscObject *)&iP));
777: if (iP) reuse = PETSC_TRUE;
778: }
779: if (!reuse) {
780: IS aB;
781: PetscInt nb;
782: PetscCall(ISGetLocalSize(pcis->is_B_local, &nb));
783: PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II), nb, 0, 1, &aB));
784: PetscCall(MatCreateSubMatrix(pcis->A_II, iV, iV, MAT_INITIAL_MATRIX, &A_II));
785: PetscCall(MatCreateSubMatrix(pcis->A_IB, iV, aB, MAT_INITIAL_MATRIX, &A_IB));
786: PetscCall(MatCreateSubMatrix(pcis->A_BI, aB, iV, MAT_INITIAL_MATRIX, &A_BI));
787: PetscCall(ISDestroy(&aB));
788: } else {
789: PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &A_IB));
790: PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &A_BI));
791: PetscCall(PetscObjectReference((PetscObject)pcis->A_II));
792: A_II = pcis->A_II;
793: }
794: PetscCall(MatCreateSchurComplement(A_II, A_II, A_IB, A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
796: /* propagate settings of solver */
797: PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j, &sksp));
798: PetscCall(KSPGetType(pcis->ksp_D, &ksptype));
799: PetscCall(KSPSetType(sksp, ksptype));
800: PetscCall(KSPGetPC(pcis->ksp_D, &pc));
801: PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &isshell));
802: if (!isshell) {
803: MatSolverType solver;
804: PCType pctype;
806: PetscCall(PCGetType(pc, &pctype));
807: PetscCall(PCFactorGetMatSolverType(pc, (MatSolverType *)&solver));
808: PetscCall(KSPGetPC(sksp, &pc));
809: PetscCall(PCSetType(pc, pctype));
810: if (solver) PetscCall(PCFactorSetMatSolverType(pc, solver));
811: } else {
812: PetscCall(KSPGetPC(sksp, &pc));
813: PetscCall(PCSetType(pc, PCLU));
814: }
815: PetscCall(MatDestroy(&A_II));
816: PetscCall(MatDestroy(&A_IB));
817: PetscCall(MatDestroy(&A_BI));
818: PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
819: PetscCall(KSPSetOptionsPrefix(sksp, prefix));
820: PetscCall(KSPAppendOptionsPrefix(sksp, "harmonic_"));
821: PetscCall(KSPSetFromOptions(sksp));
822: if (reuse) {
823: PetscCall(KSPSetPC(sksp, sub_schurs->reuse_solver->interior_solver));
824: PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver, (PetscObject)sksp, 0));
825: }
826: } else { /* default Dirichlet preconditioner is pde-harmonic */
827: PetscCall(MatCreateSchurComplement(pcis->A_II, pcis->A_II, pcis->A_IB, pcis->A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
828: PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j, pcis->ksp_D));
829: }
830: } else {
831: PetscCall(PetscObjectReference((PetscObject)pcis->A_BB));
832: fetidppc_ctx->S_j = pcis->A_BB;
833: }
834: /* saddle-point */
835: if (mat_ctx->xPg) {
836: PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg));
837: fetidppc_ctx->xPg = mat_ctx->xPg;
838: PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg));
839: fetidppc_ctx->yPg = mat_ctx->yPg;
840: }
841: PetscFunctionReturn(PETSC_SUCCESS);
842: }
844: static PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
845: {
846: FETIDPMat_ctx mat_ctx;
847: PC_BDDC *pcbddc;
848: PC_IS *pcis;
850: PetscFunctionBegin;
851: PetscCall(MatShellGetContext(fetimat, &mat_ctx));
852: pcis = (PC_IS *)mat_ctx->pc->data;
853: pcbddc = (PC_BDDC *)mat_ctx->pc->data;
854: /* Application of B_delta^T */
855: PetscCall(VecSet(pcis->vec1_B, 0.));
856: PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
857: PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
858: PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B));
860: /* Add contribution from saddle point */
861: if (mat_ctx->l2g_p) {
862: PetscCall(VecScatterBegin(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
863: PetscCall(VecScatterEnd(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
864: if (pcbddc->switch_static) {
865: if (trans) {
866: PetscCall(MatMultTranspose(mat_ctx->B_BI, mat_ctx->vP, pcis->vec1_D));
867: } else {
868: PetscCall(MatMult(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D));
869: }
870: }
871: if (trans) {
872: PetscCall(MatMultTransposeAdd(mat_ctx->B_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
873: } else {
874: PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
875: }
876: } else {
877: if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D, 0.0));
878: }
879: /* Application of \widetilde{S}^-1 */
880: PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
881: PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, trans));
882: PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
883: PetscCall(VecSet(y, 0.0));
884: /* Application of B_delta */
885: PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local));
886: /* Contribution from boundary pressures */
887: if (mat_ctx->C) {
888: const PetscScalar *lx;
889: PetscScalar *ly;
891: /* pressure ordered first in the local part of x and y */
892: PetscCall(VecGetArrayRead(x, &lx));
893: PetscCall(VecGetArray(y, &ly));
894: PetscCall(VecPlaceArray(mat_ctx->xPg, lx));
895: PetscCall(VecPlaceArray(mat_ctx->yPg, ly));
896: if (trans) {
897: PetscCall(MatMultTranspose(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
898: } else {
899: PetscCall(MatMult(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
900: }
901: PetscCall(VecResetArray(mat_ctx->xPg));
902: PetscCall(VecResetArray(mat_ctx->yPg));
903: PetscCall(VecRestoreArrayRead(x, &lx));
904: PetscCall(VecRestoreArray(y, &ly));
905: }
906: /* Add contribution from saddle point */
907: if (mat_ctx->l2g_p) {
908: if (trans) {
909: PetscCall(MatMultTranspose(mat_ctx->Bt_BB, pcis->vec1_B, mat_ctx->vP));
910: } else {
911: PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP));
912: }
913: if (pcbddc->switch_static) {
914: if (trans) {
915: PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
916: } else {
917: PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
918: }
919: }
920: PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
921: PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
922: }
923: PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
924: PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
925: PetscFunctionReturn(PETSC_SUCCESS);
926: }
928: PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
929: {
930: PetscFunctionBegin;
931: PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_FALSE));
932: PetscFunctionReturn(PETSC_SUCCESS);
933: }
935: PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
936: {
937: PetscFunctionBegin;
938: PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_TRUE));
939: PetscFunctionReturn(PETSC_SUCCESS);
940: }
942: static PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
943: {
944: FETIDPPC_ctx pc_ctx;
945: PC_IS *pcis;
947: PetscFunctionBegin;
948: PetscCall(PCShellGetContext(fetipc, &pc_ctx));
949: pcis = (PC_IS *)pc_ctx->pc->data;
950: /* Application of B_Ddelta^T */
951: PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
952: PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
953: PetscCall(VecSet(pcis->vec2_B, 0.0));
954: PetscCall(MatMultTranspose(pc_ctx->B_Ddelta, pc_ctx->lambda_local, pcis->vec2_B));
955: /* Application of local Schur complement */
956: if (trans) {
957: PetscCall(MatMultTranspose(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
958: } else {
959: PetscCall(MatMult(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
960: }
961: /* Application of B_Ddelta */
962: PetscCall(MatMult(pc_ctx->B_Ddelta, pcis->vec1_B, pc_ctx->lambda_local));
963: PetscCall(VecSet(y, 0.0));
964: PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
965: PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
966: PetscFunctionReturn(PETSC_SUCCESS);
967: }
969: PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
970: {
971: PetscFunctionBegin;
972: PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_FALSE));
973: PetscFunctionReturn(PETSC_SUCCESS);
974: }
976: PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
977: {
978: PetscFunctionBegin;
979: PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_TRUE));
980: PetscFunctionReturn(PETSC_SUCCESS);
981: }
983: PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
984: {
985: FETIDPPC_ctx pc_ctx;
986: PetscBool iascii;
987: PetscViewer sviewer;
989: PetscFunctionBegin;
990: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
991: if (iascii) {
992: PetscMPIInt rank;
993: PetscBool isschur, isshell;
995: PetscCall(PCShellGetContext(pc, &pc_ctx));
996: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
997: PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j, MATSCHURCOMPLEMENT, &isschur));
998: if (isschur) {
999: PetscCall(PetscViewerASCIIPrintf(viewer, " Dirichlet preconditioner (just from rank 0)\n"));
1000: } else {
1001: PetscCall(PetscViewerASCIIPrintf(viewer, " Lumped preconditioner (just from rank 0)\n"));
1002: }
1003: PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1004: if (rank == 0) {
1005: PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
1006: PetscCall(PetscViewerASCIIPushTab(sviewer));
1007: PetscCall(MatView(pc_ctx->S_j, sviewer));
1008: PetscCall(PetscViewerASCIIPopTab(sviewer));
1009: PetscCall(PetscViewerPopFormat(sviewer));
1010: }
1011: PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1012: PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta, MATSHELL, &isshell));
1013: if (isshell) {
1014: BDdelta_DN ctx;
1015: PetscCall(PetscViewerASCIIPrintf(viewer, " FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n"));
1016: PetscCall(MatShellGetContext(pc_ctx->B_Ddelta, &ctx));
1017: PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1018: if (rank == 0) {
1019: PetscInt tl;
1021: PetscCall(PetscViewerASCIIGetTab(sviewer, &tl));
1022: PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD, tl));
1023: PetscCall(KSPView(ctx->kBD, sviewer));
1024: PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
1025: PetscCall(MatView(ctx->BD, sviewer));
1026: PetscCall(PetscViewerPopFormat(sviewer));
1027: }
1028: PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1029: }
1030: PetscCall(PetscViewerFlush(viewer));
1031: }
1032: PetscFunctionReturn(PETSC_SUCCESS);
1033: }