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, buf_size;
130:   PetscMPIInt     rank, 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;
314:     PetscMPIInt  nreqs;

316:     partial_sum = 0;
317:     PetscCall(PetscMalloc1(pcis->n_neigh, &ptrs_buffer));
318:     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &send_reqs));
319:     PetscCall(PetscMalloc1(PetscMax(pcis->n_neigh - 1, 0), &recv_reqs));
320:     PetscCall(PetscMalloc1(pcis->n + 1, &all_factors));
321:     if (pcis->n_neigh > 0) ptrs_buffer[0] = 0;
322:     for (i = 1; i < pcis->n_neigh; i++) {
323:       partial_sum += pcis->n_shared[i];
324:       ptrs_buffer[i] = ptrs_buffer[i - 1] + pcis->n_shared[i];
325:     }
326:     PetscCall(PetscMalloc1(partial_sum, &send_buffer));
327:     PetscCall(PetscMalloc1(partial_sum, &recv_buffer));
328:     PetscCall(PetscMalloc1(partial_sum, &all_factors[0]));
329:     for (i = 0; i < pcis->n - 1; i++) {
330:       j                  = count[i];
331:       all_factors[i + 1] = all_factors[i] + j;
332:     }

334:     /* scatter B scaling to N vec */
335:     PetscCall(VecScatterBegin(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
336:     PetscCall(VecScatterEnd(pcis->N_to_B, pcis->D, pcis->vec1_N, INSERT_VALUES, SCATTER_REVERSE));
337:     /* communications */
338:     PetscCall(VecGetArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
339:     for (i = 1; i < pcis->n_neigh; i++) {
340:       for (j = 0; j < pcis->n_shared[i]; j++) send_buffer[ptrs_buffer[i - 1] + j] = array[pcis->shared[i][j]];
341:       buf_size = ptrs_buffer[i] - ptrs_buffer[i - 1];
342:       PetscCall(PetscMPIIntCast(pcis->neigh[i], &neigh));
343:       PetscCallMPI(MPIU_Isend(&send_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &send_reqs[i - 1]));
344:       PetscCallMPI(MPIU_Irecv(&recv_buffer[ptrs_buffer[i - 1]], buf_size, MPIU_SCALAR, neigh, 0, comm, &recv_reqs[i - 1]));
345:     }
346:     PetscCall(VecRestoreArrayRead(pcis->vec1_N, (const PetscScalar **)&array));
347:     PetscCall(PetscMPIIntCast(pcis->n_neigh - 1, &nreqs));
348:     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(nreqs, recv_reqs, MPI_STATUSES_IGNORE));
349:     /* put values in correct places */
350:     for (i = 1; i < pcis->n_neigh; i++) {
351:       for (j = 0; j < pcis->n_shared[i]; j++) {
352:         k              = pcis->shared[i][j];
353:         neigh_position = 0;
354:         while (neighbours_set[k][neigh_position] != pcis->neigh[i]) { neigh_position++; }
355:         all_factors[k][neigh_position] = recv_buffer[ptrs_buffer[i - 1] + j];
356:       }
357:     }
358:     if (pcis->n_neigh > 0) PetscCallMPI(MPI_Waitall(nreqs, send_reqs, MPI_STATUSES_IGNORE));
359:     PetscCall(PetscFree(send_reqs));
360:     PetscCall(PetscFree(recv_reqs));
361:     PetscCall(PetscFree(send_buffer));
362:     PetscCall(PetscFree(recv_buffer));
363:     PetscCall(PetscFree(ptrs_buffer));
364:   }

366:   /* Compute B and B_delta (local actions) */
367:   PetscCall(PetscMalloc1(pcis->n_neigh, &aux_sums));
368:   PetscCall(PetscMalloc1(n_local_lambda, &l2g_indices));
369:   PetscCall(PetscMalloc1(n_local_lambda, &vals_B_delta));
370:   PetscCall(PetscMalloc1(n_local_lambda, &cols_B_delta));
371:   if (!pcbddc->use_deluxe_scaling) {
372:     PetscCall(PetscMalloc1(n_local_lambda, &scaling_factors));
373:   } else {
374:     scaling_factors = NULL;
375:     all_factors     = NULL;
376:   }
377:   PetscCall(ISGetIndices(subset_n, &aux_global_numbering));
378:   partial_sum = 0;
379:   cum         = 0;
380:   for (i = 0; i < dual_size; i++) {
381:     n_global_lambda = aux_global_numbering[cum];
382:     j               = count[aux_local_numbering_1[i]];
383:     aux_sums[0]     = 0;
384:     for (s = 1; s < j; s++) aux_sums[s] = aux_sums[s - 1] + j - s + 1;
385:     if (all_factors) array = all_factors[aux_local_numbering_1[i]];
386:     n_neg_values = 0;
387:     while (n_neg_values < j && neighbours_set[aux_local_numbering_1[i]][n_neg_values] < rank) { n_neg_values++; }
388:     n_pos_values = j - n_neg_values;
389:     if (fully_redundant) {
390:       for (s = 0; s < n_neg_values; s++) {
391:         l2g_indices[partial_sum + s]  = aux_sums[s] + n_neg_values - s - 1 + n_global_lambda;
392:         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
393:         vals_B_delta[partial_sum + s] = -1.0;
394:         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s] = array[s];
395:       }
396:       for (s = 0; s < n_pos_values; s++) {
397:         l2g_indices[partial_sum + s + n_neg_values]  = aux_sums[n_neg_values] + s + n_global_lambda;
398:         cols_B_delta[partial_sum + s + n_neg_values] = dual_dofs_boundary_indices[i];
399:         vals_B_delta[partial_sum + s + n_neg_values] = 1.0;
400:         if (!pcbddc->use_deluxe_scaling) scaling_factors[partial_sum + s + n_neg_values] = array[s + n_neg_values];
401:       }
402:       partial_sum += j;
403:     } else {
404:       /* l2g_indices and default cols and vals of B_delta */
405:       for (s = 0; s < j; s++) {
406:         l2g_indices[partial_sum + s]  = n_global_lambda + s;
407:         cols_B_delta[partial_sum + s] = dual_dofs_boundary_indices[i];
408:         vals_B_delta[partial_sum + s] = 0.0;
409:       }
410:       /* B_delta */
411:       if (n_neg_values > 0) { /* there's a rank next to me to the left */
412:         vals_B_delta[partial_sum + n_neg_values - 1] = -1.0;
413:       }
414:       if (n_neg_values < j) { /* there's a rank next to me to the right */
415:         vals_B_delta[partial_sum + n_neg_values] = 1.0;
416:       }
417:       /* scaling as in Klawonn-Widlund 1999 */
418:       if (!pcbddc->use_deluxe_scaling) {
419:         for (s = 0; s < n_neg_values; s++) {
420:           scalar_value = 0.0;
421:           for (k = 0; k < s + 1; k++) scalar_value += array[k];
422:           scaling_factors[partial_sum + s] = -scalar_value;
423:         }
424:         for (s = 0; s < n_pos_values; s++) {
425:           scalar_value = 0.0;
426:           for (k = s + n_neg_values; k < j; k++) scalar_value += array[k];
427:           scaling_factors[partial_sum + s + n_neg_values] = scalar_value;
428:         }
429:       }
430:       partial_sum += j;
431:     }
432:     cum += aux_local_numbering_2[i];
433:   }
434:   PetscCall(ISRestoreIndices(subset_n, &aux_global_numbering));
435:   PetscCall(ISDestroy(&subset_mult));
436:   PetscCall(ISDestroy(&subset_n));
437:   PetscCall(PetscFree(aux_sums));
438:   PetscCall(PetscFree(aux_local_numbering_1));
439:   PetscCall(PetscFree(dual_dofs_boundary_indices));
440:   if (all_factors) {
441:     PetscCall(PetscFree(all_factors[0]));
442:     PetscCall(PetscFree(all_factors));
443:   }
444:   if (pcis->n) PetscCall(PetscFree(neighbours_set[0]));
445:   PetscCall(PetscFree2(count, neighbours_set));

447:   /* Create local part of B_delta */
448:   PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_delta));
449:   PetscCall(MatSetSizes(fetidpmat_ctx->B_delta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
450:   PetscCall(MatSetType(fetidpmat_ctx->B_delta, MATSEQAIJ));
451:   PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_delta, 1, NULL));
452:   PetscCall(MatSetOption(fetidpmat_ctx->B_delta, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
453:   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));
454:   PetscCall(PetscFree(vals_B_delta));
455:   PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));
456:   PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_delta, MAT_FINAL_ASSEMBLY));

458:   BD1 = NULL;
459:   BD2 = NULL;
460:   if (fully_redundant) {
461:     PetscCheck(!pcbddc->use_deluxe_scaling, comm, PETSC_ERR_SUP, "Deluxe FETIDP with fully-redundant multipliers to be implemented");
462:     PetscCall(MatCreate(PETSC_COMM_SELF, &ScalingMat));
463:     PetscCall(MatSetSizes(ScalingMat, n_local_lambda, n_local_lambda, n_local_lambda, n_local_lambda));
464:     PetscCall(MatSetType(ScalingMat, MATSEQAIJ));
465:     PetscCall(MatSeqAIJSetPreallocation(ScalingMat, 1, NULL));
466:     for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(ScalingMat, i, i, scaling_factors[i], INSERT_VALUES));
467:     PetscCall(MatAssemblyBegin(ScalingMat, MAT_FINAL_ASSEMBLY));
468:     PetscCall(MatAssemblyEnd(ScalingMat, MAT_FINAL_ASSEMBLY));
469:     PetscCall(MatMatMult(ScalingMat, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &fetidpmat_ctx->B_Ddelta));
470:     PetscCall(MatDestroy(&ScalingMat));
471:   } else {
472:     PetscCall(MatCreate(PETSC_COMM_SELF, &fetidpmat_ctx->B_Ddelta));
473:     PetscCall(MatSetSizes(fetidpmat_ctx->B_Ddelta, n_local_lambda, pcis->n_B, n_local_lambda, pcis->n_B));
474:     if (!pcbddc->use_deluxe_scaling || !pcbddc->sub_schurs) {
475:       PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSEQAIJ));
476:       PetscCall(MatSeqAIJSetPreallocation(fetidpmat_ctx->B_Ddelta, 1, NULL));
477:       for (i = 0; i < n_local_lambda; i++) PetscCall(MatSetValue(fetidpmat_ctx->B_Ddelta, i, cols_B_delta[i], scaling_factors[i], INSERT_VALUES));
478:       PetscCall(MatAssemblyBegin(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
479:       PetscCall(MatAssemblyEnd(fetidpmat_ctx->B_Ddelta, MAT_FINAL_ASSEMBLY));
480:     } else {
481:       /* scaling as in Klawonn-Widlund 1999 */
482:       PCBDDCDeluxeScaling deluxe_ctx = pcbddc->deluxe_ctx;
483:       PCBDDCSubSchurs     sub_schurs = pcbddc->sub_schurs;
484:       Mat                 T;
485:       PetscScalar        *W, lwork, *Bwork;
486:       const PetscInt     *idxs = NULL;
487:       PetscInt            cum, mss, *nnz;
488:       PetscBLASInt       *pivots, B_lwork, B_N, B_ierr;

490:       PetscCheck(pcbddc->deluxe_singlemat, comm, PETSC_ERR_USER, "Cannot compute B_Ddelta! rerun with -pc_bddc_deluxe_singlemat");
491:       mss = 0;
492:       PetscCall(PetscCalloc1(pcis->n_B, &nnz));
493:       if (sub_schurs->is_Ej_all) {
494:         PetscCall(ISGetIndices(sub_schurs->is_Ej_all, &idxs));
495:         for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
496:           PetscInt subset_size;

498:           PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
499:           for (j = 0; j < subset_size; j++) nnz[idxs[j + cum]] = subset_size;
500:           mss = PetscMax(mss, subset_size);
501:           cum += subset_size;
502:         }
503:       }
504:       PetscCall(MatCreate(PETSC_COMM_SELF, &T));
505:       PetscCall(MatSetSizes(T, pcis->n_B, pcis->n_B, pcis->n_B, pcis->n_B));
506:       PetscCall(MatSetType(T, MATSEQAIJ));
507:       PetscCall(MatSeqAIJSetPreallocation(T, 0, nnz));
508:       PetscCall(PetscFree(nnz));

510:       /* workspace allocation */
511:       B_lwork = 0;
512:       if (mss) {
513:         PetscScalar dummy = 1;

515:         B_lwork = -1;
516:         PetscCall(PetscBLASIntCast(mss, &B_N));
517:         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
518:         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, &dummy, &B_N, &B_N, &lwork, &B_lwork, &B_ierr));
519:         PetscCall(PetscFPTrapPop());
520:         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in query to GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr);
521:         PetscCall(PetscBLASIntCast((PetscInt)PetscRealPart(lwork), &B_lwork));
522:       }
523:       PetscCall(PetscMalloc3(mss * mss, &W, mss, &pivots, B_lwork, &Bwork));

525:       for (i = 0, cum = 0; i < sub_schurs->n_subs; i++) {
526:         const PetscScalar *M;
527:         PetscInt           subset_size;

529:         PetscCall(ISGetLocalSize(sub_schurs->is_subs[i], &subset_size));
530:         PetscCall(PetscBLASIntCast(subset_size, &B_N));
531:         PetscCall(MatDenseGetArrayRead(deluxe_ctx->seq_mat[i], &M));
532:         PetscCall(PetscArraycpy(W, M, subset_size * subset_size));
533:         PetscCall(MatDenseRestoreArrayRead(deluxe_ctx->seq_mat[i], &M));
534:         PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
535:         PetscCallBLAS("LAPACKgetrf", LAPACKgetrf_(&B_N, &B_N, W, &B_N, pivots, &B_ierr));
536:         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRF Lapack routine %" PetscBLASInt_FMT, B_ierr);
537:         PetscCallBLAS("LAPACKgetri", LAPACKgetri_(&B_N, W, &B_N, pivots, Bwork, &B_lwork, &B_ierr));
538:         PetscCheck(!B_ierr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in GETRI Lapack routine %" PetscBLASInt_FMT, B_ierr);
539:         PetscCall(PetscFPTrapPop());
540:         /* silent static analyzer */
541:         PetscCheck(idxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "IDXS not present");
542:         PetscCall(MatSetValues(T, subset_size, idxs + cum, subset_size, idxs + cum, W, INSERT_VALUES));
543:         cum += subset_size;
544:       }
545:       PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
546:       PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
547:       PetscCall(MatMatTransposeMult(T, fetidpmat_ctx->B_delta, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &BD1));
548:       PetscCall(MatMatMult(fetidpmat_ctx->B_delta, BD1, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &BD2));
549:       PetscCall(MatDestroy(&T));
550:       PetscCall(PetscFree3(W, pivots, Bwork));
551:       if (sub_schurs->is_Ej_all) PetscCall(ISRestoreIndices(sub_schurs->is_Ej_all, &idxs));
552:     }
553:   }
554:   PetscCall(PetscFree(scaling_factors));
555:   PetscCall(PetscFree(cols_B_delta));

557:   /* Layout of multipliers */
558:   PetscCall(PetscLayoutCreate(comm, &llay));
559:   PetscCall(PetscLayoutSetBlockSize(llay, 1));
560:   PetscCall(PetscLayoutSetSize(llay, fetidpmat_ctx->n_lambda));
561:   PetscCall(PetscLayoutSetUp(llay));
562:   PetscCall(PetscLayoutGetLocalSize(llay, &fetidpmat_ctx->n));

564:   /* Local work vector of multipliers */
565:   PetscCall(VecCreate(PETSC_COMM_SELF, &fetidpmat_ctx->lambda_local));
566:   PetscCall(VecSetSizes(fetidpmat_ctx->lambda_local, n_local_lambda, n_local_lambda));
567:   PetscCall(VecSetType(fetidpmat_ctx->lambda_local, VECSEQ));

569:   if (BD2) {
570:     ISLocalToGlobalMapping l2g;
571:     Mat                    T, TA, *pT;
572:     IS                     is;
573:     PetscInt               nl, N;
574:     BDdelta_DN             ctx;

576:     PetscCall(PetscLayoutGetLocalSize(llay, &nl));
577:     PetscCall(PetscLayoutGetSize(llay, &N));
578:     PetscCall(MatCreate(comm, &T));
579:     PetscCall(MatSetSizes(T, nl, nl, N, N));
580:     PetscCall(MatSetType(T, MATIS));
581:     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &l2g));
582:     PetscCall(MatSetLocalToGlobalMapping(T, l2g, l2g));
583:     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
584:     PetscCall(MatISSetLocalMat(T, BD2));
585:     PetscCall(MatAssemblyBegin(T, MAT_FINAL_ASSEMBLY));
586:     PetscCall(MatAssemblyEnd(T, MAT_FINAL_ASSEMBLY));
587:     PetscCall(MatDestroy(&BD2));
588:     PetscCall(MatConvert(T, MATAIJ, MAT_INITIAL_MATRIX, &TA));
589:     PetscCall(MatDestroy(&T));
590:     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_USE_POINTER, &is));
591:     PetscCall(MatCreateSubMatrices(TA, 1, &is, &is, MAT_INITIAL_MATRIX, &pT));
592:     PetscCall(MatDestroy(&TA));
593:     PetscCall(ISDestroy(&is));
594:     BD2 = pT[0];
595:     PetscCall(PetscFree(pT));

597:     /* B_Ddelta for non-redundant multipliers with deluxe scaling */
598:     PetscCall(PetscNew(&ctx));
599:     PetscCall(MatSetType(fetidpmat_ctx->B_Ddelta, MATSHELL));
600:     PetscCall(MatShellSetContext(fetidpmat_ctx->B_Ddelta, ctx));
601:     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT, (void (*)(void))MatMult_BDdelta_deluxe_nonred));
602:     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_MULT_TRANSPOSE, (void (*)(void))MatMultTranspose_BDdelta_deluxe_nonred));
603:     PetscCall(MatShellSetOperation(fetidpmat_ctx->B_Ddelta, MATOP_DESTROY, (void (*)(void))MatDestroy_BDdelta_deluxe_nonred));
604:     PetscCall(MatSetUp(fetidpmat_ctx->B_Ddelta));

606:     PetscCall(PetscObjectReference((PetscObject)BD1));
607:     ctx->BD = BD1;
608:     PetscCall(KSPCreate(PETSC_COMM_SELF, &ctx->kBD));
609:     PetscCall(KSPSetNestLevel(ctx->kBD, fetidpmat_ctx->pc->kspnestlevel));
610:     PetscCall(KSPSetOperators(ctx->kBD, BD2, BD2));
611:     PetscCall(VecDuplicate(fetidpmat_ctx->lambda_local, &ctx->work));
612:     fetidpmat_ctx->deluxe_nonred = PETSC_TRUE;
613:   }
614:   PetscCall(MatDestroy(&BD1));
615:   PetscCall(MatDestroy(&BD2));

617:   /* fetidpmat sizes */
618:   fetidpmat_ctx->n += nPgl;
619:   fetidpmat_ctx->N = fetidpmat_ctx->n_lambda + nPg;

621:   /* Global vector for FETI-DP linear system */
622:   PetscCall(VecCreate(comm, &fetidp_global));
623:   PetscCall(VecSetSizes(fetidp_global, fetidpmat_ctx->n, fetidpmat_ctx->N));
624:   PetscCall(VecSetType(fetidp_global, VECMPI));
625:   PetscCall(VecSetUp(fetidp_global));

627:   /* Decide layout for fetidp dofs: if it is a saddle point problem
628:      pressure is ordered first in the local part of the global vector
629:      of the FETI-DP linear system */
630:   if (nPg) {
631:     Vec             v;
632:     IS              IS_l2g_p, ais;
633:     PetscLayout     alay;
634:     const PetscInt *idxs, *pranges, *aranges, *lranges;
635:     PetscInt       *l2g_indices_p, rst;
636:     PetscMPIInt     rank;

638:     PetscCall(PetscMalloc1(nPl, &l2g_indices_p));
639:     PetscCall(VecGetLayout(fetidp_global, &alay));
640:     PetscCall(PetscLayoutGetRanges(alay, &aranges));
641:     PetscCall(PetscLayoutGetRanges(play, &pranges));
642:     PetscCall(PetscLayoutGetRanges(llay, &lranges));

644:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)fetidp_global), &rank));
645:     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), pranges[rank + 1] - pranges[rank], aranges[rank], 1, &fetidpmat_ctx->pressure));
646:     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->pressure, "F_P"));
647:     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)fetidp_global), lranges[rank + 1] - lranges[rank], aranges[rank] + pranges[rank + 1] - pranges[rank], 1, &fetidpmat_ctx->lagrange));
648:     PetscCall(PetscObjectSetName((PetscObject)fetidpmat_ctx->lagrange, "F_L"));
649:     PetscCall(ISLocalToGlobalMappingGetIndices(l2gmap_p, &idxs));
650:     /* shift local to global indices for pressure */
651:     for (i = 0; i < nPl; i++) {
652:       PetscMPIInt owner;

654:       PetscCall(PetscLayoutFindOwner(play, idxs[i], &owner));
655:       l2g_indices_p[i] = idxs[i] - pranges[owner] + aranges[owner];
656:     }
657:     PetscCall(ISLocalToGlobalMappingRestoreIndices(l2gmap_p, &idxs));
658:     PetscCall(ISCreateGeneral(comm, nPl, l2g_indices_p, PETSC_OWN_POINTER, &IS_l2g_p));

660:     /* local to global scatter for pressure */
661:     PetscCall(VecScatterCreate(fetidpmat_ctx->vP, NULL, fetidp_global, IS_l2g_p, &fetidpmat_ctx->l2g_p));
662:     PetscCall(ISDestroy(&IS_l2g_p));

664:     /* scatter for lagrange multipliers only */
665:     PetscCall(VecCreate(comm, &v));
666:     PetscCall(VecSetType(v, VECSTANDARD));
667:     PetscCall(VecSetLayout(v, llay));
668:     PetscCall(VecSetUp(v));
669:     PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_COPY_VALUES, &ais));
670:     PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, v, ais, &fetidpmat_ctx->l2g_lambda_only));
671:     PetscCall(ISDestroy(&ais));
672:     PetscCall(VecDestroy(&v));

674:     /* shift local to global indices for multipliers */
675:     for (i = 0; i < n_local_lambda; i++) {
676:       PetscInt    ps;
677:       PetscMPIInt owner;

679:       PetscCall(PetscLayoutFindOwner(llay, l2g_indices[i], &owner));
680:       ps             = pranges[owner + 1] - pranges[owner];
681:       l2g_indices[i] = l2g_indices[i] - lranges[owner] + aranges[owner] + ps;
682:     }

684:     /* scatter from alldofs to pressures global fetidp vector */
685:     PetscCall(PetscLayoutGetRange(alay, &rst, NULL));
686:     PetscCall(ISCreateStride(comm, nPgl, rst, 1, &ais));
687:     PetscCall(VecScatterCreate(pcis->vec1_global, pP, fetidp_global, ais, &fetidpmat_ctx->g2g_p));
688:     PetscCall(ISDestroy(&ais));
689:   }
690:   PetscCall(PetscLayoutDestroy(&llay));
691:   PetscCall(PetscLayoutDestroy(&play));
692:   PetscCall(ISCreateGeneral(comm, n_local_lambda, l2g_indices, PETSC_OWN_POINTER, &IS_l2g_lambda));

694:   /* scatter from local to global multipliers */
695:   PetscCall(VecScatterCreate(fetidpmat_ctx->lambda_local, NULL, fetidp_global, IS_l2g_lambda, &fetidpmat_ctx->l2g_lambda));
696:   PetscCall(ISDestroy(&IS_l2g_lambda));
697:   PetscCall(ISLocalToGlobalMappingDestroy(&l2gmap_p));
698:   PetscCall(VecDestroy(&fetidp_global));

700:   /* Create some work vectors needed by fetidp */
701:   PetscCall(VecDuplicate(pcis->vec1_B, &fetidpmat_ctx->temp_solution_B));
702:   PetscCall(VecDuplicate(pcis->vec1_D, &fetidpmat_ctx->temp_solution_D));
703:   PetscFunctionReturn(PETSC_SUCCESS);
704: }

706: PetscErrorCode PCBDDCSetupFETIDPPCContext(Mat fetimat, FETIDPPC_ctx fetidppc_ctx)
707: {
708:   FETIDPMat_ctx mat_ctx;
709:   PC_BDDC      *pcbddc = (PC_BDDC *)fetidppc_ctx->pc->data;
710:   PC_IS        *pcis   = (PC_IS *)fetidppc_ctx->pc->data;
711:   PetscBool     lumped = PETSC_FALSE;

713:   PetscFunctionBegin;
714:   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
715:   /* get references from objects created when setting up feti mat context */
716:   PetscCall(PetscObjectReference((PetscObject)mat_ctx->lambda_local));
717:   fetidppc_ctx->lambda_local = mat_ctx->lambda_local;
718:   PetscCall(PetscObjectReference((PetscObject)mat_ctx->B_Ddelta));
719:   fetidppc_ctx->B_Ddelta = mat_ctx->B_Ddelta;
720:   if (mat_ctx->deluxe_nonred) {
721:     PC            pc, mpc;
722:     BDdelta_DN    ctx;
723:     MatSolverType solver;
724:     const char   *prefix;

726:     PetscCall(MatShellGetContext(mat_ctx->B_Ddelta, &ctx));
727:     PetscCall(KSPSetType(ctx->kBD, KSPPREONLY));
728:     PetscCall(KSPGetPC(ctx->kBD, &mpc));
729:     PetscCall(KSPGetPC(pcbddc->ksp_D, &pc));
730:     PetscCall(PCSetType(mpc, PCLU));
731:     PetscCall(PCFactorGetMatSolverType(pc, &solver));
732:     if (solver) PetscCall(PCFactorSetMatSolverType(mpc, solver));
733:     PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
734:     PetscCall(KSPSetOptionsPrefix(ctx->kBD, prefix));
735:     PetscCall(KSPAppendOptionsPrefix(ctx->kBD, "bddelta_"));
736:     PetscCall(KSPSetFromOptions(ctx->kBD));
737:   }

739:   if (mat_ctx->l2g_lambda_only) {
740:     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda_only));
741:     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda_only;
742:   } else {
743:     PetscCall(PetscObjectReference((PetscObject)mat_ctx->l2g_lambda));
744:     fetidppc_ctx->l2g_lambda = mat_ctx->l2g_lambda;
745:   }
746:   /* Dirichlet preconditioner */
747:   PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_lumped", &lumped, NULL));
748:   if (!lumped) {
749:     IS        iV;
750:     PetscBool discrete_harmonic = PETSC_FALSE;

752:     PetscCall(PetscObjectQuery((PetscObject)fetidppc_ctx->pc, "__KSPFETIDP_iV", (PetscObject *)&iV));
753:     if (iV) PetscCall(PetscOptionsGetBool(NULL, ((PetscObject)fetimat)->prefix, "-pc_discrete_harmonic", &discrete_harmonic, NULL));
754:     if (discrete_harmonic) {
755:       KSP             sksp;
756:       PC              pc;
757:       PCBDDCSubSchurs sub_schurs = pcbddc->sub_schurs;
758:       Mat             A_II, A_IB, A_BI;
759:       IS              iP = NULL;
760:       PetscBool       isshell, reuse = PETSC_FALSE;
761:       KSPType         ksptype;
762:       const char     *prefix;

764:       /*
765:         We constructs a Schur complement for

767:         | A_II A_ID |
768:         | A_DI A_DD |

770:         instead of

772:         | A_II  B^t_II A_ID |
773:         | B_II -C_II   B_ID |
774:         | A_DI  B^t_ID A_DD |

776:       */
777:       if (sub_schurs && sub_schurs->reuse_solver) {
778:         PetscCall(PetscObjectQuery((PetscObject)sub_schurs->A, "__KSPFETIDP_iP", (PetscObject *)&iP));
779:         if (iP) reuse = PETSC_TRUE;
780:       }
781:       if (!reuse) {
782:         IS       aB;
783:         PetscInt nb;
784:         PetscCall(ISGetLocalSize(pcis->is_B_local, &nb));
785:         PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pcis->A_II), nb, 0, 1, &aB));
786:         PetscCall(MatCreateSubMatrix(pcis->A_II, iV, iV, MAT_INITIAL_MATRIX, &A_II));
787:         PetscCall(MatCreateSubMatrix(pcis->A_IB, iV, aB, MAT_INITIAL_MATRIX, &A_IB));
788:         PetscCall(MatCreateSubMatrix(pcis->A_BI, aB, iV, MAT_INITIAL_MATRIX, &A_BI));
789:         PetscCall(ISDestroy(&aB));
790:       } else {
791:         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_I_local, pcis->is_B_local, MAT_INITIAL_MATRIX, &A_IB));
792:         PetscCall(MatCreateSubMatrix(sub_schurs->A, pcis->is_B_local, pcis->is_I_local, MAT_INITIAL_MATRIX, &A_BI));
793:         PetscCall(PetscObjectReference((PetscObject)pcis->A_II));
794:         A_II = pcis->A_II;
795:       }
796:       PetscCall(MatCreateSchurComplement(A_II, A_II, A_IB, A_BI, pcis->A_BB, &fetidppc_ctx->S_j));

798:       /* propagate settings of solver */
799:       PetscCall(MatSchurComplementGetKSP(fetidppc_ctx->S_j, &sksp));
800:       PetscCall(KSPGetType(pcis->ksp_D, &ksptype));
801:       PetscCall(KSPSetType(sksp, ksptype));
802:       PetscCall(KSPGetPC(pcis->ksp_D, &pc));
803:       PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCSHELL, &isshell));
804:       if (!isshell) {
805:         MatSolverType solver;
806:         PCType        pctype;

808:         PetscCall(PCGetType(pc, &pctype));
809:         PetscCall(PCFactorGetMatSolverType(pc, &solver));
810:         PetscCall(KSPGetPC(sksp, &pc));
811:         PetscCall(PCSetType(pc, pctype));
812:         if (solver) PetscCall(PCFactorSetMatSolverType(pc, solver));
813:       } else {
814:         PetscCall(KSPGetPC(sksp, &pc));
815:         PetscCall(PCSetType(pc, PCLU));
816:       }
817:       PetscCall(MatDestroy(&A_II));
818:       PetscCall(MatDestroy(&A_IB));
819:       PetscCall(MatDestroy(&A_BI));
820:       PetscCall(MatGetOptionsPrefix(fetimat, &prefix));
821:       PetscCall(KSPSetOptionsPrefix(sksp, prefix));
822:       PetscCall(KSPAppendOptionsPrefix(sksp, "harmonic_"));
823:       PetscCall(KSPSetFromOptions(sksp));
824:       if (reuse) {
825:         PetscCall(KSPSetPC(sksp, sub_schurs->reuse_solver->interior_solver));
826:         PetscCall(PetscObjectIncrementTabLevel((PetscObject)sub_schurs->reuse_solver->interior_solver, (PetscObject)sksp, 0));
827:       }
828:     } else { /* default Dirichlet preconditioner is pde-harmonic */
829:       PetscCall(MatCreateSchurComplement(pcis->A_II, pcis->A_II, pcis->A_IB, pcis->A_BI, pcis->A_BB, &fetidppc_ctx->S_j));
830:       PetscCall(MatSchurComplementSetKSP(fetidppc_ctx->S_j, pcis->ksp_D));
831:     }
832:   } else {
833:     PetscCall(PetscObjectReference((PetscObject)pcis->A_BB));
834:     fetidppc_ctx->S_j = pcis->A_BB;
835:   }
836:   /* saddle-point */
837:   if (mat_ctx->xPg) {
838:     PetscCall(PetscObjectReference((PetscObject)mat_ctx->xPg));
839:     fetidppc_ctx->xPg = mat_ctx->xPg;
840:     PetscCall(PetscObjectReference((PetscObject)mat_ctx->yPg));
841:     fetidppc_ctx->yPg = mat_ctx->yPg;
842:   }
843:   PetscFunctionReturn(PETSC_SUCCESS);
844: }

846: static PetscErrorCode FETIDPMatMult_Kernel(Mat fetimat, Vec x, Vec y, PetscBool trans)
847: {
848:   FETIDPMat_ctx mat_ctx;
849:   PC_BDDC      *pcbddc;
850:   PC_IS        *pcis;

852:   PetscFunctionBegin;
853:   PetscCall(MatShellGetContext(fetimat, &mat_ctx));
854:   pcis   = (PC_IS *)mat_ctx->pc->data;
855:   pcbddc = (PC_BDDC *)mat_ctx->pc->data;
856:   /* Application of B_delta^T */
857:   PetscCall(VecSet(pcis->vec1_B, 0.));
858:   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
859:   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, x, mat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
860:   PetscCall(MatMultTranspose(mat_ctx->B_delta, mat_ctx->lambda_local, pcis->vec1_B));

862:   /* Add contribution from saddle point */
863:   if (mat_ctx->l2g_p) {
864:     PetscCall(VecScatterBegin(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
865:     PetscCall(VecScatterEnd(mat_ctx->l2g_p, x, mat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE));
866:     if (pcbddc->switch_static) {
867:       if (trans) {
868:         PetscCall(MatMultTranspose(mat_ctx->B_BI, mat_ctx->vP, pcis->vec1_D));
869:       } else {
870:         PetscCall(MatMult(mat_ctx->Bt_BI, mat_ctx->vP, pcis->vec1_D));
871:       }
872:     }
873:     if (trans) {
874:       PetscCall(MatMultTransposeAdd(mat_ctx->B_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
875:     } else {
876:       PetscCall(MatMultAdd(mat_ctx->Bt_BB, mat_ctx->vP, pcis->vec1_B, pcis->vec1_B));
877:     }
878:   } else {
879:     if (pcbddc->switch_static) PetscCall(VecSet(pcis->vec1_D, 0.0));
880:   }
881:   /* Application of \widetilde{S}^-1 */
882:   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
883:   PetscCall(PCBDDCApplyInterfacePreconditioner(mat_ctx->pc, trans));
884:   PetscCall(PetscArrayzero(pcbddc->benign_p0, pcbddc->benign_n));
885:   PetscCall(VecSet(y, 0.0));
886:   /* Application of B_delta */
887:   PetscCall(MatMult(mat_ctx->B_delta, pcis->vec1_B, mat_ctx->lambda_local));
888:   /* Contribution from boundary pressures */
889:   if (mat_ctx->C) {
890:     const PetscScalar *lx;
891:     PetscScalar       *ly;

893:     /* pressure ordered first in the local part of x and y */
894:     PetscCall(VecGetArrayRead(x, &lx));
895:     PetscCall(VecGetArray(y, &ly));
896:     PetscCall(VecPlaceArray(mat_ctx->xPg, lx));
897:     PetscCall(VecPlaceArray(mat_ctx->yPg, ly));
898:     if (trans) {
899:       PetscCall(MatMultTranspose(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
900:     } else {
901:       PetscCall(MatMult(mat_ctx->C, mat_ctx->xPg, mat_ctx->yPg));
902:     }
903:     PetscCall(VecResetArray(mat_ctx->xPg));
904:     PetscCall(VecResetArray(mat_ctx->yPg));
905:     PetscCall(VecRestoreArrayRead(x, &lx));
906:     PetscCall(VecRestoreArray(y, &ly));
907:   }
908:   /* Add contribution from saddle point */
909:   if (mat_ctx->l2g_p) {
910:     if (trans) {
911:       PetscCall(MatMultTranspose(mat_ctx->Bt_BB, pcis->vec1_B, mat_ctx->vP));
912:     } else {
913:       PetscCall(MatMult(mat_ctx->B_BB, pcis->vec1_B, mat_ctx->vP));
914:     }
915:     if (pcbddc->switch_static) {
916:       if (trans) {
917:         PetscCall(MatMultTransposeAdd(mat_ctx->Bt_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
918:       } else {
919:         PetscCall(MatMultAdd(mat_ctx->B_BI, pcis->vec1_D, mat_ctx->vP, mat_ctx->vP));
920:       }
921:     }
922:     PetscCall(VecScatterBegin(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
923:     PetscCall(VecScatterEnd(mat_ctx->l2g_p, mat_ctx->vP, y, ADD_VALUES, SCATTER_FORWARD));
924:   }
925:   PetscCall(VecScatterBegin(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
926:   PetscCall(VecScatterEnd(mat_ctx->l2g_lambda, mat_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
927:   PetscFunctionReturn(PETSC_SUCCESS);
928: }

930: PetscErrorCode FETIDPMatMult(Mat fetimat, Vec x, Vec y)
931: {
932:   PetscFunctionBegin;
933:   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_FALSE));
934:   PetscFunctionReturn(PETSC_SUCCESS);
935: }

937: PetscErrorCode FETIDPMatMultTranspose(Mat fetimat, Vec x, Vec y)
938: {
939:   PetscFunctionBegin;
940:   PetscCall(FETIDPMatMult_Kernel(fetimat, x, y, PETSC_TRUE));
941:   PetscFunctionReturn(PETSC_SUCCESS);
942: }

944: static PetscErrorCode FETIDPPCApply_Kernel(PC fetipc, Vec x, Vec y, PetscBool trans)
945: {
946:   FETIDPPC_ctx pc_ctx;
947:   PC_IS       *pcis;

949:   PetscFunctionBegin;
950:   PetscCall(PCShellGetContext(fetipc, &pc_ctx));
951:   pcis = (PC_IS *)pc_ctx->pc->data;
952:   /* Application of B_Ddelta^T */
953:   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
954:   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, x, pc_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE));
955:   PetscCall(VecSet(pcis->vec2_B, 0.0));
956:   PetscCall(MatMultTranspose(pc_ctx->B_Ddelta, pc_ctx->lambda_local, pcis->vec2_B));
957:   /* Application of local Schur complement */
958:   if (trans) {
959:     PetscCall(MatMultTranspose(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
960:   } else {
961:     PetscCall(MatMult(pc_ctx->S_j, pcis->vec2_B, pcis->vec1_B));
962:   }
963:   /* Application of B_Ddelta */
964:   PetscCall(MatMult(pc_ctx->B_Ddelta, pcis->vec1_B, pc_ctx->lambda_local));
965:   PetscCall(VecSet(y, 0.0));
966:   PetscCall(VecScatterBegin(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
967:   PetscCall(VecScatterEnd(pc_ctx->l2g_lambda, pc_ctx->lambda_local, y, ADD_VALUES, SCATTER_FORWARD));
968:   PetscFunctionReturn(PETSC_SUCCESS);
969: }

971: PetscErrorCode FETIDPPCApply(PC pc, Vec x, Vec y)
972: {
973:   PetscFunctionBegin;
974:   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_FALSE));
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }

978: PetscErrorCode FETIDPPCApplyTranspose(PC pc, Vec x, Vec y)
979: {
980:   PetscFunctionBegin;
981:   PetscCall(FETIDPPCApply_Kernel(pc, x, y, PETSC_TRUE));
982:   PetscFunctionReturn(PETSC_SUCCESS);
983: }

985: PetscErrorCode FETIDPPCView(PC pc, PetscViewer viewer)
986: {
987:   FETIDPPC_ctx pc_ctx;
988:   PetscBool    iascii;
989:   PetscViewer  sviewer;

991:   PetscFunctionBegin;
992:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
993:   if (iascii) {
994:     PetscMPIInt rank;
995:     PetscBool   isschur, isshell;

997:     PetscCall(PCShellGetContext(pc, &pc_ctx));
998:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
999:     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->S_j, MATSCHURCOMPLEMENT, &isschur));
1000:     if (isschur) {
1001:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Dirichlet preconditioner (just from rank 0)\n"));
1002:     } else {
1003:       PetscCall(PetscViewerASCIIPrintf(viewer, "  Lumped preconditioner (just from rank 0)\n"));
1004:     }
1005:     PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1006:     if (rank == 0) {
1007:       PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
1008:       PetscCall(PetscViewerASCIIPushTab(sviewer));
1009:       PetscCall(MatView(pc_ctx->S_j, sviewer));
1010:       PetscCall(PetscViewerASCIIPopTab(sviewer));
1011:       PetscCall(PetscViewerPopFormat(sviewer));
1012:     }
1013:     PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1014:     PetscCall(PetscObjectTypeCompare((PetscObject)pc_ctx->B_Ddelta, MATSHELL, &isshell));
1015:     if (isshell) {
1016:       BDdelta_DN ctx;
1017:       PetscCall(PetscViewerASCIIPrintf(viewer, "  FETI-DP BDdelta: DB^t * (B D^-1 B^t)^-1 for deluxe scaling (just from rank 0)\n"));
1018:       PetscCall(MatShellGetContext(pc_ctx->B_Ddelta, &ctx));
1019:       PetscCall(PetscViewerGetSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1020:       if (rank == 0) {
1021:         PetscInt tl;

1023:         PetscCall(PetscViewerASCIIGetTab(sviewer, &tl));
1024:         PetscCall(PetscObjectSetTabLevel((PetscObject)ctx->kBD, tl));
1025:         PetscCall(KSPView(ctx->kBD, sviewer));
1026:         PetscCall(PetscViewerPushFormat(sviewer, PETSC_VIEWER_ASCII_INFO));
1027:         PetscCall(MatView(ctx->BD, sviewer));
1028:         PetscCall(PetscViewerPopFormat(sviewer));
1029:       }
1030:       PetscCall(PetscViewerRestoreSubViewer(viewer, PetscObjectComm((PetscObject)pc_ctx->S_j), &sviewer));
1031:     }
1032:     PetscCall(PetscViewerFlush(viewer));
1033:   }
1034:   PetscFunctionReturn(PETSC_SUCCESS);
1035: }