Actual source code: fetidp.c

  1: #include <petsc/private/kspimpl.h>
  2: #include <petsc/private/pcbddcimpl.h>
  3: #include <petsc/private/pcbddcprivateimpl.h>
  4: #include <petscdm.h>

  6: static PetscBool  cited       = PETSC_FALSE;
  7: static PetscBool  cited2      = PETSC_FALSE;
  8: static const char citation[]  = "@article{ZampiniPCBDDC,\n"
  9:                                 "author = {Stefano Zampini},\n"
 10:                                 "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n"
 11:                                 "journal = {SIAM Journal on Scientific Computing},\n"
 12:                                 "volume = {38},\n"
 13:                                 "number = {5},\n"
 14:                                 "pages = {S282-S306},\n"
 15:                                 "year = {2016},\n"
 16:                                 "doi = {10.1137/15M1025785},\n"
 17:                                 "URL = {http://dx.doi.org/10.1137/15M1025785},\n"
 18:                                 "eprint = {http://dx.doi.org/10.1137/15M1025785}\n"
 19:                                 "}\n"
 20:                                 "@article{ZampiniDualPrimal,\n"
 21:                                 "author = {Stefano Zampini},\n"
 22:                                 "title = {{D}ual-{P}rimal methods for the cardiac {B}idomain model},\n"
 23:                                 "volume = {24},\n"
 24:                                 "number = {04},\n"
 25:                                 "pages = {667-696},\n"
 26:                                 "year = {2014},\n"
 27:                                 "doi = {10.1142/S0218202513500632},\n"
 28:                                 "URL = {https://www.worldscientific.com/doi/abs/10.1142/S0218202513500632},\n"
 29:                                 "eprint = {https://www.worldscientific.com/doi/pdf/10.1142/S0218202513500632}\n"
 30:                                 "}\n";
 31: static const char citation2[] = "@article{li2013nonoverlapping,\n"
 32:                                 "title={A nonoverlapping domain decomposition method for incompressible Stokes equations with continuous pressures},\n"
 33:                                 "author={Li, Jing and Tu, Xuemin},\n"
 34:                                 "journal={SIAM Journal on Numerical Analysis},\n"
 35:                                 "volume={51},\n"
 36:                                 "number={2},\n"
 37:                                 "pages={1235--1253},\n"
 38:                                 "year={2013},\n"
 39:                                 "publisher={Society for Industrial and Applied Mathematics}\n"
 40:                                 "}\n";

 42: /*
 43:     This file implements the FETI-DP method in PETSc as part of KSP.
 44: */
 45: typedef struct {
 46:   KSP parentksp;
 47: } KSP_FETIDPMon;

 49: typedef struct {
 50:   KSP              innerksp;        /* the KSP for the Lagrange multipliers */
 51:   PC               innerbddc;       /* the inner BDDC object */
 52:   PetscBool        fully_redundant; /* true for using a fully redundant set of multipliers */
 53:   PetscBool        userbddc;        /* true if the user provided the PCBDDC object */
 54:   PetscBool        saddlepoint;     /* support for saddle point problems */
 55:   IS               pP;              /* index set for pressure variables */
 56:   Vec              rhs_flip;        /* see KSPFETIDPSetUpOperators */
 57:   KSP_FETIDPMon   *monctx;          /* monitor context, used to pass user defined monitors
 58:                                         in the physical space */
 59:   PetscObjectState matstate;        /* these are needed just in the saddle point case */
 60:   PetscObjectState matnnzstate;     /* where we are going to use MatZeroRows on pmat */
 61:   PetscBool        statechanged;
 62:   PetscBool        check;
 63: } KSP_FETIDP;

 65: static PetscErrorCode KSPFETIDPSetPressureOperator_FETIDP(KSP ksp, Mat P)
 66: {
 67:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

 69:   if (P) fetidp->saddlepoint = PETSC_TRUE;
 70:   PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_PPmat", (PetscObject)P);
 71:   return 0;
 72: }

 74: /*@
 75:  KSPFETIDPSetPressureOperator - Sets the operator used to setup the pressure preconditioner for saddle point FETI-DP.

 77:    Collective on ksp

 79:    Input Parameters:
 80: +  ksp - the FETI-DP Krylov solver
 81: -  P - the linear operator to be preconditioned, usually the mass matrix.

 83:    Level: advanced

 85:    Notes:
 86:     The operator can be either passed in a) monolithic global ordering, b) pressure-only global ordering
 87:           or c) interface pressure ordering (if -ksp_fetidp_pressure_all false).
 88:           In cases b) and c), the pressure ordering of dofs needs to satisfy
 89:              pid_1 < pid_2  iff  gid_1 < gid_2
 90:           where pid_1 and pid_2 are two different pressure dof numbers and gid_1 and gid_2 the corresponding
 91:           id in the monolithic global ordering.

 93: .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPGetInnerBDDC`, `KSPFETIDPGetInnerKSP`, `KSPSetOperators`
 94: @*/
 95: PetscErrorCode KSPFETIDPSetPressureOperator(KSP ksp, Mat P)
 96: {
 99:   PetscTryMethod(ksp, "KSPFETIDPSetPressureOperator_C", (KSP, Mat), (ksp, P));
100:   return 0;
101: }

103: static PetscErrorCode KSPFETIDPGetInnerKSP_FETIDP(KSP ksp, KSP *innerksp)
104: {
105:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

107:   *innerksp = fetidp->innerksp;
108:   return 0;
109: }

111: /*@
112:  KSPFETIDPGetInnerKSP - Gets the KSP object for the Lagrange multipliers

114:    Input Parameters:
115: +  ksp - the FETI-DP KSP
116: -  innerksp - the KSP for the multipliers

118:    Level: advanced

120:    Notes:

122: .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPSetInnerBDDC`, `KSPFETIDPGetInnerBDDC`
123: @*/
124: PetscErrorCode KSPFETIDPGetInnerKSP(KSP ksp, KSP *innerksp)
125: {
128:   PetscUseMethod(ksp, "KSPFETIDPGetInnerKSP_C", (KSP, KSP *), (ksp, innerksp));
129:   return 0;
130: }

132: static PetscErrorCode KSPFETIDPGetInnerBDDC_FETIDP(KSP ksp, PC *pc)
133: {
134:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

136:   *pc = fetidp->innerbddc;
137:   return 0;
138: }

140: /*@
141:  KSPFETIDPGetInnerBDDC - Gets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers

143:    Input Parameters:
144: +  ksp - the FETI-DP Krylov solver
145: -  pc - the BDDC preconditioner

147:    Level: advanced

149:    Notes:

151: .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPSetInnerBDDC`, `KSPFETIDPGetInnerKSP`
152: @*/
153: PetscErrorCode KSPFETIDPGetInnerBDDC(KSP ksp, PC *pc)
154: {
157:   PetscUseMethod(ksp, "KSPFETIDPGetInnerBDDC_C", (KSP, PC *), (ksp, pc));
158:   return 0;
159: }

161: static PetscErrorCode KSPFETIDPSetInnerBDDC_FETIDP(KSP ksp, PC pc)
162: {
163:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

165:   PetscObjectReference((PetscObject)pc);
166:   PCDestroy(&fetidp->innerbddc);
167:   fetidp->innerbddc = pc;
168:   fetidp->userbddc  = PETSC_TRUE;
169:   return 0;
170: }

172: /*@
173:  KSPFETIDPSetInnerBDDC - Sets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers

175:    Collective on ksp

177:    Input Parameters:
178: +  ksp - the FETI-DP Krylov solver
179: -  pc - the BDDC preconditioner

181:    Level: advanced

183:    Notes:

185: .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPGetInnerBDDC`, `KSPFETIDPGetInnerKSP`
186: @*/
187: PetscErrorCode KSPFETIDPSetInnerBDDC(KSP ksp, PC pc)
188: {
189:   PetscBool isbddc;

193:   PetscObjectTypeCompare((PetscObject)pc, PCBDDC, &isbddc);
195:   PetscTryMethod(ksp, "KSPFETIDPSetInnerBDDC_C", (KSP, PC), (ksp, pc));
196:   return 0;
197: }

199: static PetscErrorCode KSPBuildSolution_FETIDP(KSP ksp, Vec v, Vec *V)
200: {
201:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
202:   Mat         F;
203:   Vec         Xl;

205:   KSPGetOperators(fetidp->innerksp, &F, NULL);
206:   KSPBuildSolution(fetidp->innerksp, NULL, &Xl);
207:   if (v) {
208:     PCBDDCMatFETIDPGetSolution(F, Xl, v);
209:     *V = v;
210:   } else {
211:     PCBDDCMatFETIDPGetSolution(F, Xl, *V);
212:   }
213:   return 0;
214: }

216: static PetscErrorCode KSPMonitor_FETIDP(KSP ksp, PetscInt it, PetscReal rnorm, void *ctx)
217: {
218:   KSP_FETIDPMon *monctx = (KSP_FETIDPMon *)ctx;

220:   KSPMonitor(monctx->parentksp, it, rnorm);
221:   return 0;
222: }

224: static PetscErrorCode KSPComputeEigenvalues_FETIDP(KSP ksp, PetscInt nmax, PetscReal *r, PetscReal *c, PetscInt *neig)
225: {
226:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

228:   KSPComputeEigenvalues(fetidp->innerksp, nmax, r, c, neig);
229:   return 0;
230: }

232: static PetscErrorCode KSPComputeExtremeSingularValues_FETIDP(KSP ksp, PetscReal *emax, PetscReal *emin)
233: {
234:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

236:   KSPComputeExtremeSingularValues(fetidp->innerksp, emax, emin);
237:   return 0;
238: }

240: static PetscErrorCode KSPFETIDPCheckOperators(KSP ksp, PetscViewer viewer)
241: {
242:   KSP_FETIDP     *fetidp = (KSP_FETIDP *)ksp->data;
243:   PC_BDDC        *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
244:   PC_IS          *pcis   = (PC_IS *)fetidp->innerbddc->data;
245:   Mat_IS         *matis  = (Mat_IS *)fetidp->innerbddc->pmat->data;
246:   Mat             F;
247:   FETIDPMat_ctx   fetidpmat_ctx;
248:   Vec             test_vec, test_vec_p = NULL, fetidp_global;
249:   IS              dirdofs, isvert;
250:   MPI_Comm        comm = PetscObjectComm((PetscObject)ksp);
251:   PetscScalar     sval, *array;
252:   PetscReal       val, rval;
253:   const PetscInt *vertex_indices;
254:   PetscInt        i, n_vertices;
255:   PetscBool       isascii;

258:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);
260:   PetscViewerASCIIPrintf(viewer, "----------FETI-DP MAT  --------------\n");
261:   PetscViewerASCIIAddTab(viewer, 2);
262:   KSPGetOperators(fetidp->innerksp, &F, NULL);
263:   PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO);
264:   MatView(F, viewer);
265:   PetscViewerPopFormat(viewer);
266:   PetscViewerASCIISubtractTab(viewer, 2);
267:   MatShellGetContext(F, &fetidpmat_ctx);
268:   PetscViewerASCIIPrintf(viewer, "----------FETI-DP TESTS--------------\n");
269:   PetscViewerASCIIPrintf(viewer, "All tests should return zero!\n");
270:   PetscViewerASCIIPrintf(viewer, "FETIDP MAT context in the ");
271:   if (fetidp->fully_redundant) {
272:     PetscViewerASCIIPrintf(viewer, "fully redundant case for lagrange multipliers.\n");
273:   } else {
274:     PetscViewerASCIIPrintf(viewer, "Non-fully redundant case for lagrange multiplier.\n");
275:   }
276:   PetscViewerFlush(viewer);

278:   /* Get Vertices used to define the BDDC */
279:   PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph, NULL, NULL, NULL, NULL, &isvert);
280:   ISGetLocalSize(isvert, &n_vertices);
281:   ISGetIndices(isvert, &vertex_indices);

283:   /******************************************************************/
284:   /* TEST A/B: Test numbering of global fetidp dofs                 */
285:   /******************************************************************/
286:   MatCreateVecs(F, &fetidp_global, NULL);
287:   VecDuplicate(fetidpmat_ctx->lambda_local, &test_vec);
288:   VecSet(fetidp_global, 1.0);
289:   VecSet(test_vec, 1.);
290:   VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
291:   VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
292:   if (fetidpmat_ctx->l2g_p) {
293:     VecDuplicate(fetidpmat_ctx->vP, &test_vec_p);
294:     VecSet(test_vec_p, 1.);
295:     VecScatterBegin(fetidpmat_ctx->l2g_p, fetidp_global, fetidpmat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE);
296:     VecScatterEnd(fetidpmat_ctx->l2g_p, fetidp_global, fetidpmat_ctx->vP, INSERT_VALUES, SCATTER_REVERSE);
297:   }
298:   VecAXPY(test_vec, -1.0, fetidpmat_ctx->lambda_local);
299:   VecNorm(test_vec, NORM_INFINITY, &val);
300:   VecDestroy(&test_vec);
301:   MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm);
302:   PetscViewerASCIIPrintf(viewer, "A: CHECK glob to loc: % 1.14e\n", (double)rval);

304:   if (fetidpmat_ctx->l2g_p) {
305:     VecAXPY(test_vec_p, -1.0, fetidpmat_ctx->vP);
306:     VecNorm(test_vec_p, NORM_INFINITY, &val);
307:     MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm);
308:     PetscViewerASCIIPrintf(viewer, "A: CHECK glob to loc (p): % 1.14e\n", (double)rval);
309:   }

311:   if (fetidp->fully_redundant) {
312:     VecSet(fetidp_global, 0.0);
313:     VecSet(fetidpmat_ctx->lambda_local, 0.5);
314:     VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
315:     VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
316:     VecSum(fetidp_global, &sval);
317:     val = PetscRealPart(sval) - fetidpmat_ctx->n_lambda;
318:     MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm);
319:     PetscViewerASCIIPrintf(viewer, "B: CHECK loc to glob: % 1.14e\n", (double)rval);
320:   }

322:   if (fetidpmat_ctx->l2g_p) {
323:     VecSet(pcis->vec1_N, 1.0);
324:     VecSet(pcis->vec1_global, 0.0);
325:     VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE);
326:     VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE);

328:     VecSet(fetidp_global, 0.0);
329:     VecSet(fetidpmat_ctx->vP, -1.0);
330:     VecScatterBegin(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
331:     VecScatterEnd(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
332:     VecScatterBegin(fetidpmat_ctx->g2g_p, fetidp_global, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE);
333:     VecScatterEnd(fetidpmat_ctx->g2g_p, fetidp_global, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE);
334:     VecScatterBegin(fetidpmat_ctx->g2g_p, pcis->vec1_global, fetidp_global, INSERT_VALUES, SCATTER_FORWARD);
335:     VecScatterEnd(fetidpmat_ctx->g2g_p, pcis->vec1_global, fetidp_global, INSERT_VALUES, SCATTER_FORWARD);
336:     VecSum(fetidp_global, &sval);
337:     val = PetscRealPart(sval);
338:     MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm);
339:     PetscViewerASCIIPrintf(viewer, "B: CHECK loc to glob (p): % 1.14e\n", (double)rval);
340:   }

342:   /******************************************************************/
343:   /* TEST C: It should hold B_delta*w=0, w\in\widehat{W}            */
344:   /* This is the meaning of the B matrix                            */
345:   /******************************************************************/

347:   VecSetRandom(pcis->vec1_N, NULL);
348:   VecSet(pcis->vec1_global, 0.0);
349:   VecScatterBegin(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE);
350:   VecScatterEnd(matis->rctx, pcis->vec1_N, pcis->vec1_global, ADD_VALUES, SCATTER_REVERSE);
351:   VecScatterBegin(matis->rctx, pcis->vec1_global, pcis->vec1_N, INSERT_VALUES, SCATTER_FORWARD);
352:   VecScatterEnd(matis->rctx, pcis->vec1_global, pcis->vec1_N, INSERT_VALUES, SCATTER_FORWARD);
353:   VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD);
354:   VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD);
355:   /* Action of B_delta */
356:   MatMult(fetidpmat_ctx->B_delta, pcis->vec1_B, fetidpmat_ctx->lambda_local);
357:   VecSet(fetidp_global, 0.0);
358:   VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
359:   VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
360:   VecNorm(fetidp_global, NORM_INFINITY, &val);
361:   PetscViewerASCIIPrintf(viewer, "C: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n", (double)val);

363:   /******************************************************************/
364:   /* TEST D: It should hold E_Dw = w - P_Dw w\in\widetilde{W}       */
365:   /* E_D = R_D^TR                                                   */
366:   /* P_D = B_{D,delta}^T B_{delta}                                  */
367:   /* eq.44 Mandel Tezaur and Dohrmann 2005                          */
368:   /******************************************************************/

370:   /* compute a random vector in \widetilde{W} */
371:   VecSetRandom(pcis->vec1_N, NULL);
372:   /* set zero at vertices and essential dofs */
373:   VecGetArray(pcis->vec1_N, &array);
374:   for (i = 0; i < n_vertices; i++) array[vertex_indices[i]] = 0.0;
375:   PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph, &dirdofs);
376:   if (dirdofs) {
377:     const PetscInt *idxs;
378:     PetscInt        ndir;

380:     ISGetLocalSize(dirdofs, &ndir);
381:     ISGetIndices(dirdofs, &idxs);
382:     for (i = 0; i < ndir; i++) array[idxs[i]] = 0.0;
383:     ISRestoreIndices(dirdofs, &idxs);
384:   }
385:   VecRestoreArray(pcis->vec1_N, &array);
386:   /* store w for final comparison */
387:   VecDuplicate(pcis->vec1_B, &test_vec);
388:   VecScatterBegin(pcis->N_to_B, pcis->vec1_N, test_vec, INSERT_VALUES, SCATTER_FORWARD);
389:   VecScatterEnd(pcis->N_to_B, pcis->vec1_N, test_vec, INSERT_VALUES, SCATTER_FORWARD);

391:   /* Jump operator P_D : results stored in pcis->vec1_B */
392:   /* Action of B_delta */
393:   MatMult(fetidpmat_ctx->B_delta, test_vec, fetidpmat_ctx->lambda_local);
394:   VecSet(fetidp_global, 0.0);
395:   VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
396:   VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
397:   /* Action of B_Ddelta^T */
398:   VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
399:   VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
400:   MatMultTranspose(fetidpmat_ctx->B_Ddelta, fetidpmat_ctx->lambda_local, pcis->vec1_B);

402:   /* Average operator E_D : results stored in pcis->vec2_B */
403:   PCBDDCScalingExtension(fetidpmat_ctx->pc, test_vec, pcis->vec1_global);
404:   VecScatterBegin(pcis->global_to_B, pcis->vec1_global, pcis->vec2_B, INSERT_VALUES, SCATTER_FORWARD);
405:   VecScatterEnd(pcis->global_to_B, pcis->vec1_global, pcis->vec2_B, INSERT_VALUES, SCATTER_FORWARD);

407:   /* test E_D=I-P_D */
408:   VecAXPY(pcis->vec1_B, 1.0, pcis->vec2_B);
409:   VecAXPY(pcis->vec1_B, -1.0, test_vec);
410:   VecNorm(pcis->vec1_B, NORM_INFINITY, &val);
411:   VecDestroy(&test_vec);
412:   MPI_Reduce(&val, &rval, 1, MPIU_REAL, MPIU_MAX, 0, comm);
413:   PetscViewerASCIIPrintf(viewer, "%d: CHECK infty norm of E_D + P_D - I: %1.14e\n", PetscGlobalRank, (double)val);

415:   /******************************************************************/
416:   /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W}           */
417:   /* eq.48 Mandel Tezaur and Dohrmann 2005                          */
418:   /******************************************************************/

420:   VecSetRandom(pcis->vec1_N, NULL);
421:   /* set zero at vertices and essential dofs */
422:   VecGetArray(pcis->vec1_N, &array);
423:   for (i = 0; i < n_vertices; i++) array[vertex_indices[i]] = 0.0;
424:   if (dirdofs) {
425:     const PetscInt *idxs;
426:     PetscInt        ndir;

428:     ISGetLocalSize(dirdofs, &ndir);
429:     ISGetIndices(dirdofs, &idxs);
430:     for (i = 0; i < ndir; i++) array[idxs[i]] = 0.0;
431:     ISRestoreIndices(dirdofs, &idxs);
432:   }
433:   VecRestoreArray(pcis->vec1_N, &array);

435:   /* Jump operator P_D : results stored in pcis->vec1_B */

437:   VecScatterBegin(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD);
438:   VecScatterEnd(pcis->N_to_B, pcis->vec1_N, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD);
439:   /* Action of B_delta */
440:   MatMult(fetidpmat_ctx->B_delta, pcis->vec1_B, fetidpmat_ctx->lambda_local);
441:   VecSet(fetidp_global, 0.0);
442:   VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
443:   VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, fetidp_global, ADD_VALUES, SCATTER_FORWARD);
444:   /* Action of B_Ddelta^T */
445:   VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
446:   VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
447:   MatMultTranspose(fetidpmat_ctx->B_Ddelta, fetidpmat_ctx->lambda_local, pcis->vec1_B);
448:   /* scaling */
449:   PCBDDCScalingExtension(fetidpmat_ctx->pc, pcis->vec1_B, pcis->vec1_global);
450:   VecNorm(pcis->vec1_global, NORM_INFINITY, &val);
451:   PetscViewerASCIIPrintf(viewer, "E: CHECK infty norm of R^T_D P_D: % 1.14e\n", (double)val);

453:   if (!fetidp->fully_redundant) {
454:     /******************************************************************/
455:     /* TEST F: It should holds B_{delta}B^T_{D,delta}=I               */
456:     /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005               */
457:     /******************************************************************/
458:     VecDuplicate(fetidp_global, &test_vec);
459:     VecSetRandom(fetidp_global, NULL);
460:     if (fetidpmat_ctx->l2g_p) {
461:       VecSet(fetidpmat_ctx->vP, 0.);
462:       VecScatterBegin(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, INSERT_VALUES, SCATTER_FORWARD);
463:       VecScatterEnd(fetidpmat_ctx->l2g_p, fetidpmat_ctx->vP, fetidp_global, INSERT_VALUES, SCATTER_FORWARD);
464:     }
465:     /* Action of B_Ddelta^T */
466:     VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
467:     VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidp_global, fetidpmat_ctx->lambda_local, INSERT_VALUES, SCATTER_REVERSE);
468:     MatMultTranspose(fetidpmat_ctx->B_Ddelta, fetidpmat_ctx->lambda_local, pcis->vec1_B);
469:     /* Action of B_delta */
470:     MatMult(fetidpmat_ctx->B_delta, pcis->vec1_B, fetidpmat_ctx->lambda_local);
471:     VecSet(test_vec, 0.0);
472:     VecScatterBegin(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, test_vec, ADD_VALUES, SCATTER_FORWARD);
473:     VecScatterEnd(fetidpmat_ctx->l2g_lambda, fetidpmat_ctx->lambda_local, test_vec, ADD_VALUES, SCATTER_FORWARD);
474:     VecAXPY(fetidp_global, -1., test_vec);
475:     VecNorm(fetidp_global, NORM_INFINITY, &val);
476:     PetscViewerASCIIPrintf(viewer, "E: CHECK infty norm of P^T_D - I: % 1.14e\n", (double)val);
477:     VecDestroy(&test_vec);
478:   }
479:   PetscViewerASCIIPrintf(viewer, "-------------------------------------\n");
480:   PetscViewerFlush(viewer);
481:   VecDestroy(&test_vec_p);
482:   ISDestroy(&dirdofs);
483:   VecDestroy(&fetidp_global);
484:   ISRestoreIndices(isvert, &vertex_indices);
485:   PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph, NULL, NULL, NULL, NULL, &isvert);
486:   return 0;
487: }

489: static PetscErrorCode KSPFETIDPSetUpOperators(KSP ksp)
490: {
491:   KSP_FETIDP      *fetidp = (KSP_FETIDP *)ksp->data;
492:   PC_BDDC         *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
493:   Mat              A, Ap;
494:   PetscInt         fid = -1;
495:   PetscMPIInt      size;
496:   PetscBool        ismatis, pisz = PETSC_FALSE, allp = PETSC_FALSE, schp = PETSC_FALSE;
497:   PetscBool        flip = PETSC_FALSE; /* Usually, Stokes is written (B = -\int_\Omega \nabla \cdot u q)
498:                            | A B'| | v | = | f |
499:                            | B 0 | | p | = | g |
500:                             If -ksp_fetidp_saddlepoint_flip is true, the code assumes it is written as
501:                            | A B'| | v | = | f |
502:                            |-B 0 | | p | = |-g |
503:                          */
504:   PetscObjectState matstate, matnnzstate;

506:   PetscOptionsBegin(PetscObjectComm((PetscObject)ksp), ((PetscObject)ksp)->prefix, "FETI-DP options", "PC");
507:   PetscOptionsInt("-ksp_fetidp_pressure_field", "Field id for pressures for saddle-point problems", NULL, fid, &fid, NULL);
508:   PetscOptionsBool("-ksp_fetidp_pressure_all", "Use the whole pressure set instead of just that at the interface", NULL, allp, &allp, NULL);
509:   PetscOptionsBool("-ksp_fetidp_saddlepoint_flip", "Flip the sign of the pressure-velocity (lower-left) block", NULL, flip, &flip, NULL);
510:   PetscOptionsBool("-ksp_fetidp_pressure_schur", "Use a BDDC solver for pressure", NULL, schp, &schp, NULL);
511:   PetscOptionsEnd();

513:   MPI_Comm_size(PetscObjectComm((PetscObject)ksp), &size);
514:   fetidp->saddlepoint = (fid >= 0 ? PETSC_TRUE : fetidp->saddlepoint);
515:   if (size == 1) fetidp->saddlepoint = PETSC_FALSE;

517:   KSPGetOperators(ksp, &A, &Ap);
518:   PetscObjectTypeCompare((PetscObject)A, MATIS, &ismatis);

521:   /* Quiet return if the matrix states are unchanged.
522:      Needed only for the saddle point case since it uses MatZeroRows
523:      on a matrix that may not have changed */
524:   PetscObjectStateGet((PetscObject)A, &matstate);
525:   MatGetNonzeroState(A, &matnnzstate);
526:   if (matstate == fetidp->matstate && matnnzstate == fetidp->matnnzstate) return 0;
527:   fetidp->matstate     = matstate;
528:   fetidp->matnnzstate  = matnnzstate;
529:   fetidp->statechanged = fetidp->saddlepoint;

531:   /* see if we have some fields attached */
532:   if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
533:     DM             dm;
534:     PetscContainer c;

536:     KSPGetDM(ksp, &dm);
537:     PetscObjectQuery((PetscObject)A, "_convert_nest_lfields", (PetscObject *)&c);
538:     if (dm) {
539:       IS      *fields;
540:       PetscInt nf, i;

542:       DMCreateFieldDecomposition(dm, &nf, NULL, &fields, NULL);
543:       PCBDDCSetDofsSplitting(fetidp->innerbddc, nf, fields);
544:       for (i = 0; i < nf; i++) ISDestroy(&fields[i]);
545:       PetscFree(fields);
546:     } else if (c) {
547:       MatISLocalFields lf;

549:       PetscContainerGetPointer(c, (void **)&lf);
550:       PCBDDCSetDofsSplittingLocal(fetidp->innerbddc, lf->nr, lf->rf);
551:     }
552:   }

554:   if (!fetidp->saddlepoint) {
555:     PCSetOperators(fetidp->innerbddc, A, A);
556:   } else {
557:     Mat          nA, lA, PPmat;
558:     MatNullSpace nnsp;
559:     IS           pP;
560:     PetscInt     totP;

562:     MatISGetLocalMat(A, &lA);
563:     PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lA", (PetscObject)lA);

565:     pP = fetidp->pP;
566:     if (!pP) { /* first time, need to compute pressure dofs */
567:       PC_IS                 *pcis  = (PC_IS *)fetidp->innerbddc->data;
568:       Mat_IS                *matis = (Mat_IS *)(A->data);
569:       ISLocalToGlobalMapping l2g;
570:       IS                     lP = NULL, II, pII, lPall, Pall, is1, is2;
571:       const PetscInt        *idxs;
572:       PetscInt               nl, ni, *widxs;
573:       PetscInt               i, j, n_neigh, *neigh, *n_shared, **shared, *count;
574:       PetscInt               rst, ren, n;
575:       PetscBool              ploc;

577:       MatGetLocalSize(A, &nl, NULL);
578:       MatGetOwnershipRange(A, &rst, &ren);
579:       MatGetLocalSize(lA, &n, NULL);
580:       MatISGetLocalToGlobalMapping(A, &l2g, NULL);

582:       if (!pcis->is_I_local) { /* need to compute interior dofs */
583:         PetscCalloc1(n, &count);
584:         ISLocalToGlobalMappingGetInfo(l2g, &n_neigh, &neigh, &n_shared, &shared);
585:         for (i = 1; i < n_neigh; i++)
586:           for (j = 0; j < n_shared[i]; j++) count[shared[i][j]] += 1;
587:         for (i = 0, j = 0; i < n; i++)
588:           if (!count[i]) count[j++] = i;
589:         ISLocalToGlobalMappingRestoreInfo(l2g, &n_neigh, &neigh, &n_shared, &shared);
590:         ISCreateGeneral(PETSC_COMM_SELF, j, count, PETSC_OWN_POINTER, &II);
591:       } else {
592:         PetscObjectReference((PetscObject)pcis->is_I_local);
593:         II = pcis->is_I_local;
594:       }

596:       /* interior dofs in layout */
597:       PetscArrayzero(matis->sf_leafdata, n);
598:       PetscArrayzero(matis->sf_rootdata, nl);
599:       ISGetLocalSize(II, &ni);
600:       ISGetIndices(II, &idxs);
601:       for (i = 0; i < ni; i++) matis->sf_leafdata[idxs[i]] = 1;
602:       ISRestoreIndices(II, &idxs);
603:       PetscSFReduceBegin(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE);
604:       PetscSFReduceEnd(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE);
605:       PetscMalloc1(PetscMax(nl, n), &widxs);
606:       for (i = 0, ni = 0; i < nl; i++)
607:         if (matis->sf_rootdata[i]) widxs[ni++] = i + rst;
608:       ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &pII);

610:       /* pressure dofs */
611:       Pall  = NULL;
612:       lPall = NULL;
613:       ploc  = PETSC_FALSE;
614:       if (fid < 0) { /* zero pressure block */
615:         PetscInt np;

617:         MatFindZeroDiagonals(A, &Pall);
618:         ISGetSize(Pall, &np);
619:         if (!np) { /* zero-block not found, defaults to last field (if set) */
620:           fid = pcbddc->n_ISForDofsLocal ? pcbddc->n_ISForDofsLocal - 1 : pcbddc->n_ISForDofs - 1;
621:           ISDestroy(&Pall);
622:         } else if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
623:           PCBDDCSetDofsSplitting(fetidp->innerbddc, 1, &Pall);
624:         }
625:       }
626:       if (!Pall) { /* look for registered fields */
627:         if (pcbddc->n_ISForDofsLocal) {
628:           PetscInt np;

631:           /* need a sequential IS */
632:           ISGetLocalSize(pcbddc->ISForDofsLocal[fid], &np);
633:           ISGetIndices(pcbddc->ISForDofsLocal[fid], &idxs);
634:           ISCreateGeneral(PETSC_COMM_SELF, np, idxs, PETSC_COPY_VALUES, &lPall);
635:           ISRestoreIndices(pcbddc->ISForDofsLocal[fid], &idxs);
636:           ploc = PETSC_TRUE;
637:         } else if (pcbddc->n_ISForDofs) {
639:           PetscObjectReference((PetscObject)pcbddc->ISForDofs[fid]);
640:           Pall = pcbddc->ISForDofs[fid];
641:         } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Cannot detect pressure field! Use KSPFETIDPGetInnerBDDC() + PCBDDCSetDofsSplitting or PCBDDCSetDofsSplittingLocal");
642:       }

644:       /* if the user requested the entire pressure,
645:          remove the interior pressure dofs from II (or pII) */
646:       if (allp) {
647:         if (ploc) {
648:           IS nII;
649:           ISDifference(II, lPall, &nII);
650:           ISDestroy(&II);
651:           II = nII;
652:         } else {
653:           IS nII;
654:           ISDifference(pII, Pall, &nII);
655:           ISDestroy(&pII);
656:           pII = nII;
657:         }
658:       }
659:       if (ploc) {
660:         ISDifference(lPall, II, &lP);
661:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lP", (PetscObject)lP);
662:       } else {
663:         ISDifference(Pall, pII, &pP);
664:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pP", (PetscObject)pP);
665:         /* need all local pressure dofs */
666:         PetscArrayzero(matis->sf_leafdata, n);
667:         PetscArrayzero(matis->sf_rootdata, nl);
668:         ISGetLocalSize(Pall, &ni);
669:         ISGetIndices(Pall, &idxs);
670:         for (i = 0; i < ni; i++) matis->sf_rootdata[idxs[i] - rst] = 1;
671:         ISRestoreIndices(Pall, &idxs);
672:         PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE);
673:         PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE);
674:         for (i = 0, ni = 0; i < n; i++)
675:           if (matis->sf_leafdata[i]) widxs[ni++] = i;
676:         ISCreateGeneral(PETSC_COMM_SELF, ni, widxs, PETSC_COPY_VALUES, &lPall);
677:       }

679:       if (!Pall) {
680:         PetscArrayzero(matis->sf_leafdata, n);
681:         PetscArrayzero(matis->sf_rootdata, nl);
682:         ISGetLocalSize(lPall, &ni);
683:         ISGetIndices(lPall, &idxs);
684:         for (i = 0; i < ni; i++) matis->sf_leafdata[idxs[i]] = 1;
685:         ISRestoreIndices(lPall, &idxs);
686:         PetscSFReduceBegin(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE);
687:         PetscSFReduceEnd(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE);
688:         for (i = 0, ni = 0; i < nl; i++)
689:           if (matis->sf_rootdata[i]) widxs[ni++] = i + rst;
690:         ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &Pall);
691:       }
692:       PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject)Pall);

694:       if (flip) {
695:         PetscInt npl;
696:         ISGetLocalSize(Pall, &npl);
697:         ISGetIndices(Pall, &idxs);
698:         MatCreateVecs(A, NULL, &fetidp->rhs_flip);
699:         VecSet(fetidp->rhs_flip, 1.);
700:         VecSetOption(fetidp->rhs_flip, VEC_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE);
701:         for (i = 0; i < npl; i++) VecSetValue(fetidp->rhs_flip, idxs[i], -1., INSERT_VALUES);
702:         VecAssemblyBegin(fetidp->rhs_flip);
703:         VecAssemblyEnd(fetidp->rhs_flip);
704:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_flip", (PetscObject)fetidp->rhs_flip);
705:         ISRestoreIndices(Pall, &idxs);
706:       }
707:       ISDestroy(&Pall);
708:       ISDestroy(&pII);

710:       /* local selected pressures in subdomain-wise and global ordering */
711:       PetscArrayzero(matis->sf_leafdata, n);
712:       PetscArrayzero(matis->sf_rootdata, nl);
713:       if (!ploc) {
714:         PetscInt *widxs2;

717:         ISGetLocalSize(pP, &ni);
718:         ISGetIndices(pP, &idxs);
719:         for (i = 0; i < ni; i++) matis->sf_rootdata[idxs[i] - rst] = 1;
720:         ISRestoreIndices(pP, &idxs);
721:         PetscSFBcastBegin(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE);
722:         PetscSFBcastEnd(matis->sf, MPIU_INT, matis->sf_rootdata, matis->sf_leafdata, MPI_REPLACE);
723:         for (i = 0, ni = 0; i < n; i++)
724:           if (matis->sf_leafdata[i]) widxs[ni++] = i;
725:         PetscMalloc1(ni, &widxs2);
726:         ISLocalToGlobalMappingApply(l2g, ni, widxs, widxs2);
727:         ISCreateGeneral(PETSC_COMM_SELF, ni, widxs, PETSC_COPY_VALUES, &lP);
728:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lP", (PetscObject)lP);
729:         ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs2, PETSC_OWN_POINTER, &is1);
730:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_gP", (PetscObject)is1);
731:         ISDestroy(&is1);
732:       } else {
734:         ISGetLocalSize(lP, &ni);
735:         ISGetIndices(lP, &idxs);
736:         for (i = 0; i < ni; i++)
737:           if (idxs[i] >= 0 && idxs[i] < n) matis->sf_leafdata[idxs[i]] = 1;
738:         ISRestoreIndices(lP, &idxs);
739:         PetscSFReduceBegin(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE);
740:         ISLocalToGlobalMappingApply(l2g, ni, idxs, widxs);
741:         ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &is1);
742:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_gP", (PetscObject)is1);
743:         ISDestroy(&is1);
744:         PetscSFReduceEnd(matis->sf, MPIU_INT, matis->sf_leafdata, matis->sf_rootdata, MPI_REPLACE);
745:         for (i = 0, ni = 0; i < nl; i++)
746:           if (matis->sf_rootdata[i]) widxs[ni++] = i + rst;
747:         ISCreateGeneral(PetscObjectComm((PetscObject)ksp), ni, widxs, PETSC_COPY_VALUES, &pP);
748:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pP", (PetscObject)pP);
749:       }
750:       PetscFree(widxs);

752:       /* If there's any "interior pressure",
753:          we may want to use a discrete harmonic solver instead
754:          of a Stokes harmonic for the Dirichlet preconditioner
755:          Need to extract the interior velocity dofs in interior dofs ordering (iV)
756:          and interior pressure dofs in local ordering (iP) */
757:       if (!allp) {
758:         ISLocalToGlobalMapping l2g_t;

760:         ISDifference(lPall, lP, &is1);
761:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_iP", (PetscObject)is1);
762:         ISDifference(II, is1, &is2);
763:         ISDestroy(&is1);
764:         ISLocalToGlobalMappingCreateIS(II, &l2g_t);
765:         ISGlobalToLocalMappingApplyIS(l2g_t, IS_GTOLM_DROP, is2, &is1);
766:         ISGetLocalSize(is1, &i);
767:         ISGetLocalSize(is2, &j);
769:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_iV", (PetscObject)is1);
770:         ISLocalToGlobalMappingDestroy(&l2g_t);
771:         ISDestroy(&is1);
772:         ISDestroy(&is2);
773:       }
774:       ISDestroy(&II);

776:       /* exclude selected pressures from the inner BDDC */
777:       if (pcbddc->DirichletBoundariesLocal) {
778:         IS       list[2], plP, isout;
779:         PetscInt np;

781:         /* need a parallel IS */
782:         ISGetLocalSize(lP, &np);
783:         ISGetIndices(lP, &idxs);
784:         ISCreateGeneral(PetscObjectComm((PetscObject)ksp), np, idxs, PETSC_USE_POINTER, &plP);
785:         list[0] = plP;
786:         list[1] = pcbddc->DirichletBoundariesLocal;
787:         ISConcatenate(PetscObjectComm((PetscObject)ksp), 2, list, &isout);
788:         ISSortRemoveDups(isout);
789:         ISDestroy(&plP);
790:         ISRestoreIndices(lP, &idxs);
791:         PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc, isout);
792:         ISDestroy(&isout);
793:       } else if (pcbddc->DirichletBoundaries) {
794:         IS list[2], isout;

796:         list[0] = pP;
797:         list[1] = pcbddc->DirichletBoundaries;
798:         ISConcatenate(PetscObjectComm((PetscObject)ksp), 2, list, &isout);
799:         ISSortRemoveDups(isout);
800:         PCBDDCSetDirichletBoundaries(fetidp->innerbddc, isout);
801:         ISDestroy(&isout);
802:       } else {
803:         IS       plP;
804:         PetscInt np;

806:         /* need a parallel IS */
807:         ISGetLocalSize(lP, &np);
808:         ISGetIndices(lP, &idxs);
809:         ISCreateGeneral(PetscObjectComm((PetscObject)ksp), np, idxs, PETSC_COPY_VALUES, &plP);
810:         PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc, plP);
811:         ISDestroy(&plP);
812:         ISRestoreIndices(lP, &idxs);
813:       }

815:       /* save CSR information for the pressure BDDC solver (if any) */
816:       if (schp) {
817:         PetscInt np, nt;

819:         MatGetSize(matis->A, &nt, NULL);
820:         ISGetLocalSize(lP, &np);
821:         if (np) {
822:           PetscInt *xadj = pcbddc->mat_graph->xadj;
823:           PetscInt *adjn = pcbddc->mat_graph->adjncy;
824:           PetscInt  nv   = pcbddc->mat_graph->nvtxs_csr;

826:           if (nv && nv == nt) {
827:             ISLocalToGlobalMapping pmap;
828:             PetscInt              *schp_csr, *schp_xadj, *schp_adjn, p;
829:             PetscContainer         c;

831:             ISLocalToGlobalMappingCreateIS(lPall, &pmap);
832:             ISGetIndices(lPall, &idxs);
833:             for (p = 0, nv = 0; p < np; p++) {
834:               PetscInt x, n = idxs[p];

836:               ISGlobalToLocalMappingApply(pmap, IS_GTOLM_DROP, xadj[n + 1] - xadj[n], adjn + xadj[n], &x, NULL);
837:               nv += x;
838:             }
839:             PetscMalloc1(np + 1 + nv, &schp_csr);
840:             schp_xadj = schp_csr;
841:             schp_adjn = schp_csr + np + 1;
842:             for (p = 0, schp_xadj[0] = 0; p < np; p++) {
843:               PetscInt x, n = idxs[p];

845:               ISGlobalToLocalMappingApply(pmap, IS_GTOLM_DROP, xadj[n + 1] - xadj[n], adjn + xadj[n], &x, schp_adjn + schp_xadj[p]);
846:               schp_xadj[p + 1] = schp_xadj[p] + x;
847:             }
848:             ISRestoreIndices(lPall, &idxs);
849:             ISLocalToGlobalMappingDestroy(&pmap);
850:             PetscContainerCreate(PETSC_COMM_SELF, &c);
851:             PetscContainerSetPointer(c, schp_csr);
852:             PetscContainerSetUserDestroy(c, PetscContainerUserDestroyDefault);
853:             PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pCSR", (PetscObject)c);
854:             PetscContainerDestroy(&c);
855:           }
856:         }
857:       }
858:       ISDestroy(&lPall);
859:       ISDestroy(&lP);
860:       fetidp->pP = pP;
861:     }

863:     /* total number of selected pressure dofs */
864:     ISGetSize(fetidp->pP, &totP);

866:     /* Set operator for inner BDDC */
867:     if (totP || fetidp->rhs_flip) {
868:       MatDuplicate(A, MAT_COPY_VALUES, &nA);
869:     } else {
870:       PetscObjectReference((PetscObject)A);
871:       nA = A;
872:     }
873:     if (fetidp->rhs_flip) {
874:       MatDiagonalScale(nA, fetidp->rhs_flip, NULL);
875:       if (totP) {
876:         Mat lA2;

878:         MatISGetLocalMat(nA, &lA);
879:         MatDuplicate(lA, MAT_COPY_VALUES, &lA2);
880:         PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lA", (PetscObject)lA2);
881:         MatDestroy(&lA2);
882:       }
883:     }

885:     if (totP) {
886:       MatSetOption(nA, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE);
887:       MatZeroRowsColumnsIS(nA, fetidp->pP, 1., NULL, NULL);
888:     } else {
889:       PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_lA", NULL);
890:     }
891:     MatGetNearNullSpace(Ap, &nnsp);
892:     if (!nnsp) MatGetNullSpace(Ap, &nnsp);
893:     if (!nnsp) MatGetNearNullSpace(A, &nnsp);
894:     if (!nnsp) MatGetNullSpace(A, &nnsp);
895:     MatSetNearNullSpace(nA, nnsp);
896:     PCSetOperators(fetidp->innerbddc, nA, nA);
897:     MatDestroy(&nA);

899:     /* non-zero rhs on interior dofs when applying the preconditioner */
900:     if (totP) pcbddc->switch_static = PETSC_TRUE;

902:     /* if there are no interface pressures, set inner bddc flag for benign saddle point */
903:     if (!totP) {
904:       pcbddc->benign_saddle_point = PETSC_TRUE;
905:       pcbddc->compute_nonetflux   = PETSC_TRUE;
906:     }

908:     /* Operators for pressure preconditioner */
909:     if (totP) {
910:       /* Extract pressure block if needed */
911:       if (!pisz) {
912:         Mat C;
913:         IS  nzrows = NULL;

915:         MatCreateSubMatrix(A, fetidp->pP, fetidp->pP, MAT_INITIAL_MATRIX, &C);
916:         MatFindNonzeroRows(C, &nzrows);
917:         if (nzrows) {
918:           PetscInt i;

920:           ISGetSize(nzrows, &i);
921:           ISDestroy(&nzrows);
922:           if (!i) pisz = PETSC_TRUE;
923:         }
924:         if (!pisz) {
925:           MatScale(C, -1.); /* i.e. Almost Incompressible Elasticity, Stokes discretized with Q1xQ1_stabilized */
926:           PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_C", (PetscObject)C);
927:         }
928:         MatDestroy(&C);
929:       }
930:       /* Divergence mat */
931:       if (!pcbddc->divudotp) {
932:         Mat       B;
933:         IS        P;
934:         IS        l2l = NULL;
935:         PetscBool save;

937:         PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject *)&P);
938:         if (!pisz) {
939:           IS       F, V;
940:           PetscInt m, M;

942:           MatGetOwnershipRange(A, &m, &M);
943:           ISCreateStride(PetscObjectComm((PetscObject)A), M - m, m, 1, &F);
944:           ISComplement(P, m, M, &V);
945:           MatCreateSubMatrix(A, P, V, MAT_INITIAL_MATRIX, &B);
946:           {
947:             Mat_IS *Bmatis = (Mat_IS *)B->data;
948:             PetscObjectReference((PetscObject)Bmatis->getsub_cis);
949:             l2l = Bmatis->getsub_cis;
950:           }
951:           ISDestroy(&V);
952:           ISDestroy(&F);
953:         } else {
954:           MatCreateSubMatrix(A, P, NULL, MAT_INITIAL_MATRIX, &B);
955:         }
956:         save = pcbddc->compute_nonetflux; /* SetDivergenceMat activates nonetflux computation */
957:         PCBDDCSetDivergenceMat(fetidp->innerbddc, B, PETSC_FALSE, l2l);
958:         pcbddc->compute_nonetflux = save;
959:         MatDestroy(&B);
960:         ISDestroy(&l2l);
961:       }
962:       if (A != Ap) { /* user has provided a different Pmat, this always superseeds the setter (TODO: is it OK?) */
963:         /* use monolithic operator, we restrict later */
964:         KSPFETIDPSetPressureOperator(ksp, Ap);
965:       }
966:       PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_PPmat", (PetscObject *)&PPmat);

968:       /* PPmat not present, use some default choice */
969:       if (!PPmat) {
970:         Mat C;

972:         PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_C", (PetscObject *)&C);
973:         if (!schp && C) { /* non-zero pressure block, most likely Almost Incompressible Elasticity */
974:           KSPFETIDPSetPressureOperator(ksp, C);
975:         } else if (!pisz && schp) { /* we need the whole pressure mass matrix to define the interface BDDC */
976:           IS P;

978:           PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject *)&P);
979:           MatCreateSubMatrix(A, P, P, MAT_INITIAL_MATRIX, &C);
980:           MatScale(C, -1.);
981:           KSPFETIDPSetPressureOperator(ksp, C);
982:           MatDestroy(&C);
983:         } else { /* identity (need to be scaled properly by the user using e.g. a Richardson method */
984:           PetscInt nl;

986:           ISGetLocalSize(fetidp->pP, &nl);
987:           MatCreate(PetscObjectComm((PetscObject)ksp), &C);
988:           MatSetSizes(C, nl, nl, totP, totP);
989:           MatSetType(C, MATAIJ);
990:           MatMPIAIJSetPreallocation(C, 1, NULL, 0, NULL);
991:           MatSeqAIJSetPreallocation(C, 1, NULL);
992:           MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY);
993:           MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY);
994:           MatShift(C, 1.);
995:           KSPFETIDPSetPressureOperator(ksp, C);
996:           MatDestroy(&C);
997:         }
998:       }

1000:       /* Preconditioned operator for the pressure block */
1001:       PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_PPmat", (PetscObject *)&PPmat);
1002:       if (PPmat) {
1003:         Mat      C;
1004:         IS       Pall;
1005:         PetscInt AM, PAM, PAN, pam, pan, am, an, pl, pIl, pAg, pIg;

1007:         PetscObjectQuery((PetscObject)fetidp->innerbddc, "__KSPFETIDP_aP", (PetscObject *)&Pall);
1008:         MatGetSize(A, &AM, NULL);
1009:         MatGetSize(PPmat, &PAM, &PAN);
1010:         ISGetSize(Pall, &pAg);
1011:         ISGetSize(fetidp->pP, &pIg);
1012:         MatGetLocalSize(PPmat, &pam, &pan);
1013:         MatGetLocalSize(A, &am, &an);
1014:         ISGetLocalSize(Pall, &pIl);
1015:         ISGetLocalSize(fetidp->pP, &pl);
1020:         if (PAM == AM) { /* monolithic ordering, restrict to pressure */
1021:           if (schp) {
1022:             MatCreateSubMatrix(PPmat, Pall, Pall, MAT_INITIAL_MATRIX, &C);
1023:           } else {
1024:             MatCreateSubMatrix(PPmat, fetidp->pP, fetidp->pP, MAT_INITIAL_MATRIX, &C);
1025:           }
1026:         } else if (pAg == PAM) { /* global ordering for pressure only */
1027:           if (!allp && !schp) {  /* solving for interface pressure only */
1028:             IS restr;

1030:             ISRenumber(fetidp->pP, NULL, NULL, &restr);
1031:             MatCreateSubMatrix(PPmat, restr, restr, MAT_INITIAL_MATRIX, &C);
1032:             ISDestroy(&restr);
1033:           } else {
1034:             PetscObjectReference((PetscObject)PPmat);
1035:             C = PPmat;
1036:           }
1037:         } else if (pIg == PAM) { /* global ordering for selected pressure only */
1039:           PetscObjectReference((PetscObject)PPmat);
1040:           C = PPmat;
1041:         } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_USER, "Unable to use the pressure matrix");

1043:         KSPFETIDPSetPressureOperator(ksp, C);
1044:         MatDestroy(&C);
1045:       } else SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "Missing Pmat for pressure block");
1046:     } else { /* totP == 0 */
1047:       PetscObjectCompose((PetscObject)fetidp->innerbddc, "__KSPFETIDP_pP", NULL);
1048:     }
1049:   }
1050:   return 0;
1051: }

1053: static PetscErrorCode KSPSetUp_FETIDP(KSP ksp)
1054: {
1055:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1056:   PC_BDDC    *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
1057:   PetscBool   flg;

1059:   KSPFETIDPSetUpOperators(ksp);
1060:   /* set up BDDC */
1061:   PCSetErrorIfFailure(fetidp->innerbddc, ksp->errorifnotconverged);
1062:   PCSetUp(fetidp->innerbddc);
1063:   /* FETI-DP as it is implemented needs an exact coarse solver */
1064:   if (pcbddc->coarse_ksp) {
1065:     KSPSetTolerances(pcbddc->coarse_ksp, PETSC_SMALL, PETSC_SMALL, PETSC_DEFAULT, 1000);
1066:     KSPSetNormType(pcbddc->coarse_ksp, KSP_NORM_DEFAULT);
1067:   }
1068:   /* FETI-DP as it is implemented needs exact local Neumann solvers */
1069:   KSPSetTolerances(pcbddc->ksp_R, PETSC_SMALL, PETSC_SMALL, PETSC_DEFAULT, 1000);
1070:   KSPSetNormType(pcbddc->ksp_R, KSP_NORM_DEFAULT);

1072:   /* setup FETI-DP operators
1073:      If fetidp->statechanged is true, we need to update the operators
1074:      needed in the saddle-point case. This should be replaced
1075:      by a better logic when the FETI-DP matrix and preconditioner will
1076:      have their own classes */
1077:   if (pcbddc->new_primal_space || fetidp->statechanged) {
1078:     Mat F; /* the FETI-DP matrix */
1079:     PC  D; /* the FETI-DP preconditioner */
1080:     KSPReset(fetidp->innerksp);
1081:     PCBDDCCreateFETIDPOperators(fetidp->innerbddc, fetidp->fully_redundant, ((PetscObject)ksp)->prefix, &F, &D);
1082:     KSPSetOperators(fetidp->innerksp, F, F);
1083:     KSPSetTolerances(fetidp->innerksp, ksp->rtol, ksp->abstol, ksp->divtol, ksp->max_it);
1084:     KSPSetPC(fetidp->innerksp, D);
1085:     PetscObjectIncrementTabLevel((PetscObject)D, (PetscObject)fetidp->innerksp, 0);
1086:     KSPSetFromOptions(fetidp->innerksp);
1087:     MatCreateVecs(F, &(fetidp->innerksp)->vec_rhs, &(fetidp->innerksp)->vec_sol);
1088:     MatDestroy(&F);
1089:     PCDestroy(&D);
1090:     if (fetidp->check) {
1091:       PetscViewer viewer;

1093:       if (!pcbddc->dbg_viewer) {
1094:         viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
1095:       } else {
1096:         viewer = pcbddc->dbg_viewer;
1097:       }
1098:       KSPFETIDPCheckOperators(ksp, viewer);
1099:     }
1100:   }
1101:   fetidp->statechanged     = PETSC_FALSE;
1102:   pcbddc->new_primal_space = PETSC_FALSE;

1104:   /* propagate settings to the inner solve */
1105:   KSPGetComputeSingularValues(ksp, &flg);
1106:   KSPSetComputeSingularValues(fetidp->innerksp, flg);
1107:   if (ksp->res_hist) KSPSetResidualHistory(fetidp->innerksp, ksp->res_hist, ksp->res_hist_max, ksp->res_hist_reset);
1108:   KSPSetErrorIfNotConverged(fetidp->innerksp, ksp->errorifnotconverged);
1109:   KSPSetUp(fetidp->innerksp);
1110:   return 0;
1111: }

1113: static PetscErrorCode KSPSolve_FETIDP(KSP ksp)
1114: {
1115:   Mat                F, A;
1116:   MatNullSpace       nsp;
1117:   Vec                X, B, Xl, Bl;
1118:   KSP_FETIDP        *fetidp = (KSP_FETIDP *)ksp->data;
1119:   PC_BDDC           *pcbddc = (PC_BDDC *)fetidp->innerbddc->data;
1120:   KSPConvergedReason reason;
1121:   PC                 pc;
1122:   PCFailedReason     pcreason;
1123:   PetscInt           hist_len;

1125:   PetscCitationsRegister(citation, &cited);
1126:   if (fetidp->saddlepoint) PetscCitationsRegister(citation2, &cited2);
1127:   KSPGetOperators(ksp, &A, NULL);
1128:   KSPGetRhs(ksp, &B);
1129:   KSPGetSolution(ksp, &X);
1130:   KSPGetOperators(fetidp->innerksp, &F, NULL);
1131:   KSPGetRhs(fetidp->innerksp, &Bl);
1132:   KSPGetSolution(fetidp->innerksp, &Xl);
1133:   PCBDDCMatFETIDPGetRHS(F, B, Bl);
1134:   if (ksp->transpose_solve) {
1135:     KSPSolveTranspose(fetidp->innerksp, Bl, Xl);
1136:   } else {
1137:     KSPSolve(fetidp->innerksp, Bl, Xl);
1138:   }
1139:   KSPGetConvergedReason(fetidp->innerksp, &reason);
1140:   KSPGetPC(fetidp->innerksp, &pc);
1141:   PCGetFailedReason(pc, &pcreason);
1142:   if ((reason < 0 && reason != KSP_DIVERGED_ITS) || pcreason) {
1143:     PetscInt its;
1144:     KSPGetIterationNumber(fetidp->innerksp, &its);
1145:     ksp->reason = KSP_DIVERGED_PC_FAILED;
1146:     VecSetInf(Xl);
1147:     PetscInfo(ksp, "Inner KSP solve failed: %s %s at iteration %" PetscInt_FMT, KSPConvergedReasons[reason], PCFailedReasons[pcreason], its);
1148:   }
1149:   PCBDDCMatFETIDPGetSolution(F, Xl, X);
1150:   MatGetNullSpace(A, &nsp);
1151:   if (nsp) MatNullSpaceRemove(nsp, X);
1152:   /* update ksp with stats from inner ksp */
1153:   KSPGetConvergedReason(fetidp->innerksp, &ksp->reason);
1154:   KSPGetIterationNumber(fetidp->innerksp, &ksp->its);
1155:   ksp->totalits += ksp->its;
1156:   KSPGetResidualHistory(fetidp->innerksp, NULL, &hist_len);
1157:   ksp->res_hist_len = (size_t)hist_len;
1158:   /* restore defaults for inner BDDC (Pre/PostSolve flags) */
1159:   pcbddc->temp_solution_used        = PETSC_FALSE;
1160:   pcbddc->rhs_change                = PETSC_FALSE;
1161:   pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1162:   return 0;
1163: }

1165: static PetscErrorCode KSPReset_FETIDP(KSP ksp)
1166: {
1167:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1168:   PC_BDDC    *pcbddc;

1170:   ISDestroy(&fetidp->pP);
1171:   VecDestroy(&fetidp->rhs_flip);
1172:   /* avoid PCReset that does not take into account ref counting */
1173:   PCDestroy(&fetidp->innerbddc);
1174:   PCCreate(PetscObjectComm((PetscObject)ksp), &fetidp->innerbddc);
1175:   PCSetType(fetidp->innerbddc, PCBDDC);
1176:   pcbddc                   = (PC_BDDC *)fetidp->innerbddc->data;
1177:   pcbddc->symmetric_primal = PETSC_FALSE;
1178:   KSPDestroy(&fetidp->innerksp);
1179:   fetidp->saddlepoint  = PETSC_FALSE;
1180:   fetidp->matstate     = -1;
1181:   fetidp->matnnzstate  = -1;
1182:   fetidp->statechanged = PETSC_TRUE;
1183:   return 0;
1184: }

1186: static PetscErrorCode KSPDestroy_FETIDP(KSP ksp)
1187: {
1188:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

1190:   KSPReset_FETIDP(ksp);
1191:   PCDestroy(&fetidp->innerbddc);
1192:   KSPDestroy(&fetidp->innerksp);
1193:   PetscFree(fetidp->monctx);
1194:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetInnerBDDC_C", NULL);
1195:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerBDDC_C", NULL);
1196:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerKSP_C", NULL);
1197:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetPressureOperator_C", NULL);
1198:   PetscFree(ksp->data);
1199:   return 0;
1200: }

1202: static PetscErrorCode KSPView_FETIDP(KSP ksp, PetscViewer viewer)
1203: {
1204:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;
1205:   PetscBool   iascii;

1207:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
1208:   if (iascii) {
1209:     PetscViewerASCIIPrintf(viewer, "  fully redundant: %d\n", fetidp->fully_redundant);
1210:     PetscViewerASCIIPrintf(viewer, "  saddle point:    %d\n", fetidp->saddlepoint);
1211:     PetscViewerASCIIPrintf(viewer, "Inner KSP solver details\n");
1212:   }
1213:   PetscViewerASCIIPushTab(viewer);
1214:   KSPView(fetidp->innerksp, viewer);
1215:   PetscViewerASCIIPopTab(viewer);
1216:   if (iascii) PetscViewerASCIIPrintf(viewer, "Inner BDDC solver details\n");
1217:   PetscViewerASCIIPushTab(viewer);
1218:   PCView(fetidp->innerbddc, viewer);
1219:   PetscViewerASCIIPopTab(viewer);
1220:   return 0;
1221: }

1223: static PetscErrorCode KSPSetFromOptions_FETIDP(KSP ksp, PetscOptionItems *PetscOptionsObject)
1224: {
1225:   KSP_FETIDP *fetidp = (KSP_FETIDP *)ksp->data;

1227:   /* set options prefixes for the inner objects, since the parent prefix will be valid at this point */
1228:   PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerksp, ((PetscObject)ksp)->prefix);
1229:   PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerksp, "fetidp_");
1230:   if (!fetidp->userbddc) {
1231:     PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerbddc, ((PetscObject)ksp)->prefix);
1232:     PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerbddc, "fetidp_bddc_");
1233:   }
1234:   PetscOptionsHeadBegin(PetscOptionsObject, "KSP FETIDP options");
1235:   PetscOptionsBool("-ksp_fetidp_fullyredundant", "Use fully redundant multipliers", "none", fetidp->fully_redundant, &fetidp->fully_redundant, NULL);
1236:   PetscOptionsBool("-ksp_fetidp_saddlepoint", "Activates support for saddle-point problems", NULL, fetidp->saddlepoint, &fetidp->saddlepoint, NULL);
1237:   PetscOptionsBool("-ksp_fetidp_check", "Activates verbose debugging output FETI-DP operators", NULL, fetidp->check, &fetidp->check, NULL);
1238:   PetscOptionsHeadEnd();
1239:   PCSetFromOptions(fetidp->innerbddc);
1240:   return 0;
1241: }

1243: /*MC
1244:      KSPFETIDP - The FETI-DP method

1246:    This class implements the FETI-DP method [1].
1247:    The matrix for the KSP must be of type MATIS.
1248:    The FETI-DP linear system (automatically generated constructing an internal PCBDDC object) is solved using an internal KSP object.

1250:    Options Database Keys:
1251: +   -ksp_fetidp_fullyredundant <false>   - use a fully redundant set of Lagrange multipliers
1252: .   -ksp_fetidp_saddlepoint <false>      - activates support for saddle point problems, see [2]
1253: .   -ksp_fetidp_saddlepoint_flip <false> - usually, an incompressible Stokes problem is written as
1254:                                            | A B^T | | v | = | f |
1255:                                            | B 0   | | p | = | g |
1256:                                            with B representing -\int_\Omega \nabla \cdot u q.
1257:                                            If -ksp_fetidp_saddlepoint_flip is true, the code assumes that the user provides it as
1258:                                            | A B^T | | v | = | f |
1259:                                            |-B 0   | | p | = |-g |
1260: .   -ksp_fetidp_pressure_field <-1>      - activates support for saddle point problems, and identifies the pressure field id.
1261:                                            If this information is not provided, the pressure field is detected by using MatFindZeroDiagonals().
1262: -   -ksp_fetidp_pressure_all <false>     - if false, uses the interface pressures, as described in [2]. If true, uses the entire pressure field.

1264:    Level: Advanced

1266:    Notes:
1267:     Options for the inner KSP and for the customization of the PCBDDC object can be specified at command line by using the prefixes -fetidp_ and -fetidp_bddc_. E.g.,
1268: .vb
1269:       -fetidp_ksp_type gmres -fetidp_bddc_pc_bddc_symmetric false
1270: .ve
1271:    will use GMRES for the solution of the linear system on the Lagrange multipliers, generated using a non-symmetric PCBDDC.

1273:    For saddle point problems with continuous pressures, the preconditioned operator for the pressure solver can be specified with KSPFETIDPSetPressureOperator().
1274:    Alternatively, the pressure operator is extracted from the precondioned matrix (if it is different from the linear solver matrix).
1275:    If none of the above, an identity matrix will be created; the user then needs to scale it through a Richardson solver.
1276:    Options for the pressure solver can be prefixed with -fetidp_fielsplit_p_, E.g.
1277: .vb
1278:       -fetidp_fielsplit_p_ksp_type preonly -fetidp_fielsplit_p_pc_type lu -fetidp_fielsplit_p_pc_factor_mat_solver_type mumps
1279: .ve
1280:    In order to use the deluxe version of FETI-DP, you must customize the inner BDDC operator with -fetidp_bddc_pc_bddc_use_deluxe_scaling -fetidp_bddc_pc_bddc_deluxe_singlemat and use
1281:    non-redundant multipliers, i.e. -ksp_fetidp_fullyredundant false. Options for the scaling solver are prefixed by -fetidp_bddelta_, E.g.
1282: .vb
1283:       -fetidp_bddelta_pc_factor_mat_solver_type mumps -fetidp_bddelta_pc_type lu
1284: .ve

1286:    Some of the basic options such as the maximum number of iterations and tolerances are automatically passed from this KSP to the inner KSP that actually performs the iterations.

1288:    The converged reason and number of iterations computed are passed from the inner KSP to this KSP at the end of the solution.

1290:    Developer Notes:
1291:     Even though this method does not directly use any norms, the user is allowed to set the KSPNormType to any value.
1292:     This is so users do not have to change KSPNormType options when they switch from other KSP methods to this one.

1294:    References:
1295: +  * - C. Farhat, M. Lesoinne, P. LeTallec, K. Pierson, and D. Rixen, FETI-DP: a dual-primal unified FETI method. I. A faster alternative to the two-level FETI method, Internat. J. Numer. Methods Engrg., 50 (2001), pp. 1523--1544
1296: -  * - X. Tu, J. Li, A FETI-DP type domain decomposition algorithm for three-dimensional incompressible Stokes equations, SIAM J. Numer. Anal., 53 (2015), pp. 720-742

1298: .seealso: `MATIS`, `PCBDDC`, `KSPFETIDPSetInnerBDDC()`, `KSPFETIDPGetInnerBDDC()`, `KSPFETIDPGetInnerKSP()`
1299: M*/
1300: PETSC_EXTERN PetscErrorCode KSPCreate_FETIDP(KSP ksp)
1301: {
1302:   KSP_FETIDP    *fetidp;
1303:   KSP_FETIDPMon *monctx;
1304:   PC_BDDC       *pcbddc;
1305:   PC             pc;

1307:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_LEFT, 3);
1308:   KSPSetSupportedNorm(ksp, KSP_NORM_NONE, PC_RIGHT, 2);
1309:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_LEFT, 2);
1310:   KSPSetSupportedNorm(ksp, KSP_NORM_PRECONDITIONED, PC_RIGHT, 2);
1311:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_LEFT, 2);
1312:   KSPSetSupportedNorm(ksp, KSP_NORM_UNPRECONDITIONED, PC_RIGHT, 2);
1313:   KSPSetSupportedNorm(ksp, KSP_NORM_NATURAL, PC_LEFT, 2);

1315:   PetscNew(&fetidp);
1316:   fetidp->matstate     = -1;
1317:   fetidp->matnnzstate  = -1;
1318:   fetidp->statechanged = PETSC_TRUE;

1320:   ksp->data                              = (void *)fetidp;
1321:   ksp->ops->setup                        = KSPSetUp_FETIDP;
1322:   ksp->ops->solve                        = KSPSolve_FETIDP;
1323:   ksp->ops->destroy                      = KSPDestroy_FETIDP;
1324:   ksp->ops->computeeigenvalues           = KSPComputeEigenvalues_FETIDP;
1325:   ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_FETIDP;
1326:   ksp->ops->view                         = KSPView_FETIDP;
1327:   ksp->ops->setfromoptions               = KSPSetFromOptions_FETIDP;
1328:   ksp->ops->buildsolution                = KSPBuildSolution_FETIDP;
1329:   ksp->ops->buildresidual                = KSPBuildResidualDefault;
1330:   KSPGetPC(ksp, &pc);
1331:   PCSetType(pc, PCNONE);
1332:   /* create the inner KSP for the Lagrange multipliers */
1333:   KSPCreate(PetscObjectComm((PetscObject)ksp), &fetidp->innerksp);
1334:   KSPGetPC(fetidp->innerksp, &pc);
1335:   PCSetType(pc, PCNONE);
1336:   /* monitor */
1337:   PetscNew(&monctx);
1338:   monctx->parentksp = ksp;
1339:   fetidp->monctx    = monctx;
1340:   KSPMonitorSet(fetidp->innerksp, KSPMonitor_FETIDP, fetidp->monctx, NULL);
1341:   /* create the inner BDDC */
1342:   PCCreate(PetscObjectComm((PetscObject)ksp), &fetidp->innerbddc);
1343:   PCSetType(fetidp->innerbddc, PCBDDC);
1344:   /* make sure we always obtain a consistent FETI-DP matrix
1345:      for symmetric problems, the user can always customize it through the command line */
1346:   pcbddc                   = (PC_BDDC *)fetidp->innerbddc->data;
1347:   pcbddc->symmetric_primal = PETSC_FALSE;
1348:   /* composed functions */
1349:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetInnerBDDC_C", KSPFETIDPSetInnerBDDC_FETIDP);
1350:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerBDDC_C", KSPFETIDPGetInnerBDDC_FETIDP);
1351:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPGetInnerKSP_C", KSPFETIDPGetInnerKSP_FETIDP);
1352:   PetscObjectComposeFunction((PetscObject)ksp, "KSPFETIDPSetPressureOperator_C", KSPFETIDPSetPressureOperator_FETIDP);
1353:   /* need to call KSPSetUp_FETIDP even with KSP_SETUP_NEWMATRIX */
1354:   ksp->setupnewmatrix = PETSC_TRUE;
1355:   return 0;
1356: }