Actual source code: fetidp.c
petsc-3.10.5 2019-03-28
1: #include <petsc/private/kspimpl.h> /*I <petscksp.h> I*/
2: #include <../src/ksp/pc/impls/bddc/bddc.h>
3: #include <../src/ksp/pc/impls/bddc/bddcprivate.h>
4: #include <petscdm.h>
6: static PetscBool cited = PETSC_FALSE;
7: static PetscBool cited2 = PETSC_FALSE;
8: static const char citation[] =
9: "@article{ZampiniPCBDDC,\n"
10: "author = {Stefano Zampini},\n"
11: "title = {{PCBDDC}: A Class of Robust Dual-Primal Methods in {PETS}c},\n"
12: "journal = {SIAM Journal on Scientific Computing},\n"
13: "volume = {38},\n"
14: "number = {5},\n"
15: "pages = {S282-S306},\n"
16: "year = {2016},\n"
17: "doi = {10.1137/15M1025785},\n"
18: "URL = {http://dx.doi.org/10.1137/15M1025785},\n"
19: "eprint = {http://dx.doi.org/10.1137/15M1025785}\n"
20: "}\n"
21: "@article{ZampiniDualPrimal,\n"
22: "author = {Stefano Zampini},\n"
23: "title = {{D}ual-{P}rimal methods for the cardiac {B}idomain model},\n"
24: "volume = {24},\n"
25: "number = {04},\n"
26: "pages = {667-696},\n"
27: "year = {2014},\n"
28: "doi = {10.1142/S0218202513500632},\n"
29: "URL = {http://www.worldscientific.com/doi/abs/10.1142/S0218202513500632},\n"
30: "eprint = {http://www.worldscientific.com/doi/pdf/10.1142/S0218202513500632}\n"
31: "}\n";
32: static const char citation2[] =
33: "@article{li2013nonoverlapping,\n"
34: "title={A nonoverlapping domain decomposition method for incompressible Stokes equations with continuous pressures},\n"
35: "author={Li, Jing and Tu, Xuemin},\n"
36: "journal={SIAM Journal on Numerical Analysis},\n"
37: "volume={51},\n"
38: "number={2},\n"
39: "pages={1235--1253},\n"
40: "year={2013},\n"
41: "publisher={Society for Industrial and Applied Mathematics}\n"
42: "}\n";
44: /*
45: This file implements the FETI-DP method in PETSc as part of KSP.
46: */
47: typedef struct {
48: KSP parentksp;
49: } KSP_FETIDPMon;
51: typedef struct {
52: KSP innerksp; /* the KSP for the Lagrange multipliers */
53: PC innerbddc; /* the inner BDDC object */
54: PetscBool fully_redundant; /* true for using a fully redundant set of multipliers */
55: PetscBool userbddc; /* true if the user provided the PCBDDC object */
56: PetscBool saddlepoint; /* support for saddle point problems */
57: IS pP; /* index set for pressure variables */
58: Vec rhs_flip; /* see KSPFETIDPSetUpOperators */
59: KSP_FETIDPMon *monctx; /* monitor context, used to pass user defined monitors
60: in the physical space */
61: PetscObjectState matstate; /* these are needed just in the saddle point case */
62: PetscObjectState matnnzstate; /* where we are going to use MatZeroRows on pmat */
63: PetscBool statechanged;
64: PetscBool check;
65: } KSP_FETIDP;
67: static PetscErrorCode KSPFETIDPSetPressureOperator_FETIDP(KSP ksp, Mat P)
68: {
69: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
73: if (P) fetidp->saddlepoint = PETSC_TRUE;
74: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)P);
75: return(0);
76: }
78: /*@
79: KSPFETIDPSetPressureOperator - Sets the operator used to setup the pressure preconditioner for saddle point FETI-DP.
81: Collective on KSP
83: Input Parameters:
84: + ksp - the FETI-DP Krylov solver
85: - P - the linear operator to be preconditioned, usually the mass matrix.
87: Level: advanced
89: Notes:
90: The operator can be either passed in a) monolithic global ordering, b) pressure-only global ordering
91: or c) interface pressure ordering (if -ksp_fetidp_pressure_all false).
92: In cases b) and c), the pressure ordering of dofs needs to satisfy
93: pid_1 < pid_2 iff gid_1 < gid_2
94: where pid_1 and pid_2 are two different pressure dof numbers and gid_1 and gid_2 the corresponding
95: id in the monolithic global ordering.
97: .seealso: MATIS, PCBDDC, KSPFETIDPGetInnerBDDC, KSPFETIDPGetInnerKSP, KSPSetOperators
98: @*/
99: PetscErrorCode KSPFETIDPSetPressureOperator(KSP ksp, Mat P)
100: {
106: PetscTryMethod(ksp,"KSPFETIDPSetPressureOperator_C",(KSP,Mat),(ksp,P));
107: return(0);
108: }
110: static PetscErrorCode KSPFETIDPGetInnerKSP_FETIDP(KSP ksp, KSP* innerksp)
111: {
112: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
115: *innerksp = fetidp->innerksp;
116: return(0);
117: }
119: /*@
120: KSPFETIDPGetInnerKSP - Gets the KSP object for the Lagrange multipliers
122: Input Parameters:
123: + ksp - the FETI-DP KSP
124: - innerksp - the KSP for the multipliers
126: Level: advanced
128: Notes:
130: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC, KSPFETIDPGetInnerBDDC
131: @*/
132: PetscErrorCode KSPFETIDPGetInnerKSP(KSP ksp, KSP* innerksp)
133: {
139: PetscUseMethod(ksp,"KSPFETIDPGetInnerKSP_C",(KSP,KSP*),(ksp,innerksp));
140: return(0);
141: }
143: static PetscErrorCode KSPFETIDPGetInnerBDDC_FETIDP(KSP ksp, PC* pc)
144: {
145: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
148: *pc = fetidp->innerbddc;
149: return(0);
150: }
152: /*@
153: KSPFETIDPGetInnerBDDC - Gets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers
155: Input Parameters:
156: + ksp - the FETI-DP Krylov solver
157: - pc - the BDDC preconditioner
159: Level: advanced
161: Notes:
163: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC, KSPFETIDPGetInnerKSP
164: @*/
165: PetscErrorCode KSPFETIDPGetInnerBDDC(KSP ksp, PC* pc)
166: {
172: PetscUseMethod(ksp,"KSPFETIDPGetInnerBDDC_C",(KSP,PC*),(ksp,pc));
173: return(0);
174: }
176: static PetscErrorCode KSPFETIDPSetInnerBDDC_FETIDP(KSP ksp, PC pc)
177: {
178: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
182: PetscObjectReference((PetscObject)pc);
183: PCDestroy(&fetidp->innerbddc);
184: fetidp->innerbddc = pc;
185: fetidp->userbddc = PETSC_TRUE;
186: return(0);
187: }
189: /*@
190: KSPFETIDPSetInnerBDDC - Sets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers
192: Collective on KSP
194: Input Parameters:
195: + ksp - the FETI-DP Krylov solver
196: - pc - the BDDC preconditioner
198: Level: advanced
200: Notes:
202: .seealso: MATIS, PCBDDC, KSPFETIDPGetInnerBDDC, KSPFETIDPGetInnerKSP
203: @*/
204: PetscErrorCode KSPFETIDPSetInnerBDDC(KSP ksp, PC pc)
205: {
206: PetscBool isbddc;
212: PetscObjectTypeCompare((PetscObject)pc,PCBDDC,&isbddc);
213: if (!isbddc) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"KSPFETIDPSetInnerBDDC need a PCBDDC preconditioner");
214: PetscTryMethod(ksp,"KSPFETIDPSetInnerBDDC_C",(KSP,PC),(ksp,pc));
215: return(0);
216: }
218: static PetscErrorCode KSPBuildSolution_FETIDP(KSP ksp,Vec v,Vec *V)
219: {
220: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
221: Mat F;
222: Vec Xl;
226: KSPGetOperators(fetidp->innerksp,&F,NULL);
227: KSPBuildSolution(fetidp->innerksp,NULL,&Xl);
228: if (v) {
229: PCBDDCMatFETIDPGetSolution(F,Xl,v);
230: *V = v;
231: } else {
232: PCBDDCMatFETIDPGetSolution(F,Xl,*V);
233: }
234: return(0);
235: }
237: static PetscErrorCode KSPMonitor_FETIDP(KSP ksp,PetscInt it,PetscReal rnorm,void* ctx)
238: {
239: KSP_FETIDPMon *monctx = (KSP_FETIDPMon*)ctx;
243: KSPMonitor(monctx->parentksp,it,rnorm);
244: return(0);
245: }
247: static PetscErrorCode KSPComputeEigenvalues_FETIDP(KSP ksp,PetscInt nmax,PetscReal *r,PetscReal *c,PetscInt *neig)
248: {
249: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
253: KSPComputeEigenvalues(fetidp->innerksp,nmax,r,c,neig);
254: return(0);
255: }
257: static PetscErrorCode KSPComputeExtremeSingularValues_FETIDP(KSP ksp,PetscReal *emax,PetscReal *emin)
258: {
259: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
263: KSPComputeExtremeSingularValues(fetidp->innerksp,emax,emin);
264: return(0);
265: }
267: static PetscErrorCode KSPFETIDPCheckOperators(KSP ksp, PetscViewer viewer)
268: {
269: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
270: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
271: PC_IS *pcis = (PC_IS*)fetidp->innerbddc->data;
272: Mat_IS *matis = (Mat_IS*)fetidp->innerbddc->pmat->data;
273: Mat F;
274: FETIDPMat_ctx fetidpmat_ctx;
275: Vec test_vec,test_vec_p = NULL,fetidp_global;
276: IS dirdofs,isvert;
277: MPI_Comm comm = PetscObjectComm((PetscObject)ksp);
278: PetscScalar sval,*array;
279: PetscReal val,rval;
280: const PetscInt *vertex_indices;
281: PetscInt i,n_vertices;
282: PetscBool isascii;
287: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
288: if (!isascii) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported viewer");
289: PetscViewerASCIIPrintf(viewer,"----------FETI-DP MAT --------------\n");
290: PetscViewerASCIIAddTab(viewer,2);
291: KSPGetOperators(fetidp->innerksp,&F,NULL);
292: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
293: MatView(F,viewer);
294: PetscViewerPopFormat(viewer);
295: PetscViewerASCIISubtractTab(viewer,2);
296: MatShellGetContext(F,(void**)&fetidpmat_ctx);
297: PetscViewerASCIIPrintf(viewer,"----------FETI-DP TESTS--------------\n");
298: PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");
299: PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");
300: if (fetidp->fully_redundant) {
301: PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");
302: } else {
303: PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");
304: }
305: PetscViewerFlush(viewer);
307: /* Get Vertices used to define the BDDC */
308: PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);
309: ISGetLocalSize(isvert,&n_vertices);
310: ISGetIndices(isvert,&vertex_indices);
312: /******************************************************************/
313: /* TEST A/B: Test numbering of global fetidp dofs */
314: /******************************************************************/
315: MatCreateVecs(F,&fetidp_global,NULL);
316: VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);
317: VecSet(fetidp_global,1.0);
318: VecSet(test_vec,1.);
319: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
320: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
321: if (fetidpmat_ctx->l2g_p) {
322: VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);
323: VecSet(test_vec_p,1.);
324: VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
325: VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
326: }
327: VecAXPY(test_vec,-1.0,fetidpmat_ctx->lambda_local);
328: VecNorm(test_vec,NORM_INFINITY,&val);
329: VecDestroy(&test_vec);
330: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
331: PetscViewerASCIIPrintf(viewer,"A: CHECK glob to loc: % 1.14e\n",rval);
333: if (fetidpmat_ctx->l2g_p) {
334: VecAXPY(test_vec_p,-1.0,fetidpmat_ctx->vP);
335: VecNorm(test_vec_p,NORM_INFINITY,&val);
336: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
337: PetscViewerASCIIPrintf(viewer,"A: CHECK glob to loc (p): % 1.14e\n",rval);
338: }
340: if (fetidp->fully_redundant) {
341: VecSet(fetidp_global,0.0);
342: VecSet(fetidpmat_ctx->lambda_local,0.5);
343: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
344: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
345: VecSum(fetidp_global,&sval);
346: val = PetscRealPart(sval)-fetidpmat_ctx->n_lambda;
347: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
348: PetscViewerASCIIPrintf(viewer,"B: CHECK loc to glob: % 1.14e\n",rval);
349: }
351: if (fetidpmat_ctx->l2g_p) {
352: VecSet(pcis->vec1_N,1.0);
353: VecSet(pcis->vec1_global,0.0);
354: VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
355: VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
357: VecSet(fetidp_global,0.0);
358: VecSet(fetidpmat_ctx->vP,-1.0);
359: VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
360: VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
361: VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
362: VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
363: VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
364: VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
365: VecSum(fetidp_global,&sval);
366: val = PetscRealPart(sval);
367: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
368: PetscViewerASCIIPrintf(viewer,"B: CHECK loc to glob (p): % 1.14e\n",rval);
369: }
371: /******************************************************************/
372: /* TEST C: It should hold B_delta*w=0, w\in\widehat{W} */
373: /* This is the meaning of the B matrix */
374: /******************************************************************/
376: VecSetRandom(pcis->vec1_N,NULL);
377: VecSet(pcis->vec1_global,0.0);
378: VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
379: VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
380: VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
381: VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
382: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
383: VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
384: /* Action of B_delta */
385: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
386: VecSet(fetidp_global,0.0);
387: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
388: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
389: VecNorm(fetidp_global,NORM_INFINITY,&val);
390: PetscViewerASCIIPrintf(viewer,"C: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",val);
392: /******************************************************************/
393: /* TEST D: It should hold E_Dw = w - P_Dw w\in\widetilde{W} */
394: /* E_D = R_D^TR */
395: /* P_D = B_{D,delta}^T B_{delta} */
396: /* eq.44 Mandel Tezaur and Dohrmann 2005 */
397: /******************************************************************/
399: /* compute a random vector in \widetilde{W} */
400: VecSetRandom(pcis->vec1_N,NULL);
401: /* set zero at vertices and essential dofs */
402: VecGetArray(pcis->vec1_N,&array);
403: for (i=0;i<n_vertices;i++) array[vertex_indices[i]] = 0.0;
404: PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirdofs);
405: if (dirdofs) {
406: const PetscInt *idxs;
407: PetscInt ndir;
409: ISGetLocalSize(dirdofs,&ndir);
410: ISGetIndices(dirdofs,&idxs);
411: for (i=0;i<ndir;i++) array[idxs[i]] = 0.0;
412: ISRestoreIndices(dirdofs,&idxs);
413: }
414: VecRestoreArray(pcis->vec1_N,&array);
415: /* store w for final comparison */
416: VecDuplicate(pcis->vec1_B,&test_vec);
417: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
418: VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
420: /* Jump operator P_D : results stored in pcis->vec1_B */
421: /* Action of B_delta */
422: MatMult(fetidpmat_ctx->B_delta,test_vec,fetidpmat_ctx->lambda_local);
423: VecSet(fetidp_global,0.0);
424: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
425: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
426: /* Action of B_Ddelta^T */
427: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
428: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
429: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
431: /* Average operator E_D : results stored in pcis->vec2_B */
432: PCBDDCScalingExtension(fetidpmat_ctx->pc,test_vec,pcis->vec1_global);
433: VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
434: VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
436: /* test E_D=I-P_D */
437: VecAXPY(pcis->vec1_B,1.0,pcis->vec2_B);
438: VecAXPY(pcis->vec1_B,-1.0,test_vec);
439: VecNorm(pcis->vec1_B,NORM_INFINITY,&val);
440: VecDestroy(&test_vec);
441: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
442: PetscViewerASCIIPrintf(viewer,"D: CHECK infty norm of E_D + P_D - I: % 1.14e\n",PetscGlobalRank,val);
444: /******************************************************************/
445: /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W} */
446: /* eq.48 Mandel Tezaur and Dohrmann 2005 */
447: /******************************************************************/
449: VecSetRandom(pcis->vec1_N,NULL);
450: /* set zero at vertices and essential dofs */
451: VecGetArray(pcis->vec1_N,&array);
452: for (i=0;i<n_vertices;i++) array[vertex_indices[i]] = 0.0;
453: if (dirdofs) {
454: const PetscInt *idxs;
455: PetscInt ndir;
457: ISGetLocalSize(dirdofs,&ndir);
458: ISGetIndices(dirdofs,&idxs);
459: for (i=0;i<ndir;i++) array[idxs[i]] = 0.0;
460: ISRestoreIndices(dirdofs,&idxs);
461: }
462: VecRestoreArray(pcis->vec1_N,&array);
464: /* Jump operator P_D : results stored in pcis->vec1_B */
466: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
467: VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
468: /* Action of B_delta */
469: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
470: VecSet(fetidp_global,0.0);
471: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
472: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
473: /* Action of B_Ddelta^T */
474: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
475: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
476: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
477: /* scaling */
478: PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);
479: VecNorm(pcis->vec1_global,NORM_INFINITY,&val);
480: PetscViewerASCIIPrintf(viewer,"E: CHECK infty norm of R^T_D P_D: % 1.14e\n",val);
482: if (!fetidp->fully_redundant) {
483: /******************************************************************/
484: /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */
485: /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */
486: /******************************************************************/
487: VecDuplicate(fetidp_global,&test_vec);
488: VecSetRandom(fetidp_global,NULL);
489: if (fetidpmat_ctx->l2g_p) {
490: VecSet(fetidpmat_ctx->vP,0.);
491: VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
492: VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
493: }
494: /* Action of B_Ddelta^T */
495: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
496: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
497: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
498: /* Action of B_delta */
499: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
500: VecSet(test_vec,0.0);
501: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
502: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
503: VecAXPY(fetidp_global,-1.,test_vec);
504: VecNorm(fetidp_global,NORM_INFINITY,&val);
505: PetscViewerASCIIPrintf(viewer,"E: CHECK infty norm of P^T_D - I: % 1.14e\n",val);
506: VecDestroy(&test_vec);
507: }
508: PetscViewerASCIIPrintf(viewer,"-------------------------------------\n");
509: PetscViewerFlush(viewer);
510: VecDestroy(&test_vec_p);
511: ISDestroy(&dirdofs);
512: VecDestroy(&fetidp_global);
513: ISRestoreIndices(isvert,&vertex_indices);
514: PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);
515: return(0);
516: }
518: static PetscErrorCode KSPFETIDPSetUpOperators(KSP ksp)
519: {
520: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
521: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
522: Mat A,Ap;
523: PetscInt fid = -1;
524: PetscMPIInt size;
525: PetscBool ismatis,pisz,allp,schp;
526: PetscBool flip; /* Usually, Stokes is written (B = -\int_\Omega \nabla \cdot u q)
527: | A B'| | v | = | f |
528: | B 0 | | p | = | g |
529: If -ksp_fetidp_saddlepoint_flip is true, the code assumes it is written as
530: | A B'| | v | = | f |
531: |-B 0 | | p | = |-g |
532: */
533: PetscObjectState matstate, matnnzstate;
534: PetscErrorCode ierr;
537: pisz = PETSC_FALSE;
538: flip = PETSC_FALSE;
539: allp = PETSC_FALSE;
540: schp = PETSC_FALSE;
541: PetscOptionsBegin(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"FETI-DP options","PC");
542: PetscOptionsInt("-ksp_fetidp_pressure_field","Field id for pressures for saddle-point problems",NULL,fid,&fid,NULL);
543: PetscOptionsBool("-ksp_fetidp_pressure_all","Use the whole pressure set instead of just that at the interface",NULL,allp,&allp,NULL);
544: PetscOptionsBool("-ksp_fetidp_saddlepoint_flip","Flip the sign of the pressure-velocity (lower-left) block",NULL,flip,&flip,NULL);
545: PetscOptionsBool("-ksp_fetidp_pressure_schur","Use a BDDC solver for pressure",NULL,schp,&schp,NULL);
546: PetscOptionsEnd();
548: MPI_Comm_size(PetscObjectComm((PetscObject)ksp),&size);
549: fetidp->saddlepoint = (fid >= 0 ? PETSC_TRUE : fetidp->saddlepoint);
550: if (size == 1) fetidp->saddlepoint = PETSC_FALSE;
552: KSPGetOperators(ksp,&A,&Ap);
553: PetscObjectTypeCompare((PetscObject)A,MATIS,&ismatis);
554: if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Amat should be of type MATIS");
556: /* Quiet return if the matrix states are unchanged.
557: Needed only for the saddle point case since it uses MatZeroRows
558: on a matrix that may not have changed */
559: PetscObjectStateGet((PetscObject)A,&matstate);
560: MatGetNonzeroState(A,&matnnzstate);
561: if (matstate == fetidp->matstate && matnnzstate == fetidp->matnnzstate) return(0);
562: fetidp->matstate = matstate;
563: fetidp->matnnzstate = matnnzstate;
564: fetidp->statechanged = fetidp->saddlepoint;
566: /* see if we have some fields attached */
567: if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
568: DM dm;
569: PetscContainer c;
571: KSPGetDM(ksp,&dm);
572: PetscObjectQuery((PetscObject)A,"_convert_nest_lfields",(PetscObject*)&c);
573: if (dm) {
574: IS *fields;
575: PetscInt nf,i;
577: DMCreateFieldDecomposition(dm,&nf,NULL,&fields,NULL);
578: PCBDDCSetDofsSplitting(fetidp->innerbddc,nf,fields);
579: for (i=0;i<nf;i++) {
580: ISDestroy(&fields[i]);
581: }
582: PetscFree(fields);
583: } else if (c) {
584: MatISLocalFields lf;
586: PetscContainerGetPointer(c,(void**)&lf);
587: PCBDDCSetDofsSplittingLocal(fetidp->innerbddc,lf->nr,lf->rf);
588: }
589: }
591: if (!fetidp->saddlepoint) {
592: PCSetOperators(fetidp->innerbddc,A,A);
593: } else {
594: Mat nA,lA,PPmat;
595: MatNullSpace nnsp;
596: IS pP;
597: PetscInt totP;
599: MatISGetLocalMat(A,&lA);
600: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",(PetscObject)lA);
602: pP = fetidp->pP;
603: if (!pP) { /* first time, need to compute pressure dofs */
604: PC_IS *pcis = (PC_IS*)fetidp->innerbddc->data;
605: Mat_IS *matis = (Mat_IS*)(A->data);
606: ISLocalToGlobalMapping l2g;
607: IS lP = NULL,II,pII,lPall,Pall,is1,is2;
608: const PetscInt *idxs;
609: PetscInt nl,ni,*widxs;
610: PetscInt i,j,n_neigh,*neigh,*n_shared,**shared,*count;
611: PetscInt rst,ren,n;
612: PetscBool ploc;
614: MatGetLocalSize(A,&nl,NULL);
615: MatGetOwnershipRange(A,&rst,&ren);
616: MatGetLocalSize(lA,&n,NULL);
617: MatGetLocalToGlobalMapping(A,&l2g,NULL);
619: if (!pcis->is_I_local) { /* need to compute interior dofs */
620: PetscCalloc1(n,&count);
621: ISLocalToGlobalMappingGetInfo(l2g,&n_neigh,&neigh,&n_shared,&shared);
622: for (i=1;i<n_neigh;i++)
623: for (j=0;j<n_shared[i];j++)
624: count[shared[i][j]] += 1;
625: for (i=0,j=0;i<n;i++) if (!count[i]) count[j++] = i;
626: ISLocalToGlobalMappingRestoreInfo(l2g,&n_neigh,&neigh,&n_shared,&shared);
627: ISCreateGeneral(PETSC_COMM_SELF,j,count,PETSC_OWN_POINTER,&II);
628: } else {
629: PetscObjectReference((PetscObject)pcis->is_I_local);
630: II = pcis->is_I_local;
631: }
633: /* interior dofs in layout */
634: MatISSetUpSF(A);
635: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
636: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
637: ISGetLocalSize(II,&ni);
638: ISGetIndices(II,&idxs);
639: for (i=0;i<ni;i++) matis->sf_leafdata[idxs[i]] = 1;
640: ISRestoreIndices(II,&idxs);
641: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
642: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
643: PetscMalloc1(PetscMax(nl,n),&widxs);
644: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
645: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&pII);
647: /* pressure dofs */
648: Pall = NULL;
649: lPall = NULL;
650: ploc = PETSC_FALSE;
651: if (fid < 0) { /* zero pressure block */
652: PetscInt np;
654: MatFindZeroDiagonals(A,&Pall);
655: ISGetSize(Pall,&np);
656: if (!np) { /* zero-block not found, defaults to last field (if set) */
657: fid = pcbddc->n_ISForDofsLocal ? pcbddc->n_ISForDofsLocal - 1 : pcbddc->n_ISForDofs - 1;
658: ISDestroy(&Pall);
659: } else if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
660: PCBDDCSetDofsSplitting(fetidp->innerbddc,1,&Pall);
661: }
662: }
663: if (!Pall) { /* look for registered fields */
664: if (pcbddc->n_ISForDofsLocal) {
665: PetscInt np;
667: if (fid < 0 || fid >= pcbddc->n_ISForDofsLocal) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Invalid field id for pressure %D, max %D",fid,pcbddc->n_ISForDofsLocal);
668: /* need a sequential IS */
669: ISGetLocalSize(pcbddc->ISForDofsLocal[fid],&np);
670: ISGetIndices(pcbddc->ISForDofsLocal[fid],&idxs);
671: ISCreateGeneral(PETSC_COMM_SELF,np,idxs,PETSC_COPY_VALUES,&lPall);
672: ISRestoreIndices(pcbddc->ISForDofsLocal[fid],&idxs);
673: ploc = PETSC_TRUE;
674: } else if (pcbddc->n_ISForDofs) {
675: if (fid < 0 || fid >= pcbddc->n_ISForDofs) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Invalid field id for pressure %D, max %D",fid,pcbddc->n_ISForDofs);
676: PetscObjectReference((PetscObject)pcbddc->ISForDofs[fid]);
677: Pall = pcbddc->ISForDofs[fid];
678: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Cannot detect pressure field! Use KSPFETIDPGetInnerBDDC() + PCBDDCSetDofsSplitting or PCBDDCSetDofsSplittingLocal");
679: }
681: /* if the user requested the entire pressure,
682: remove the interior pressure dofs from II (or pII) */
683: if (allp) {
684: if (ploc) {
685: IS nII;
686: ISDifference(II,lPall,&nII);
687: ISDestroy(&II);
688: II = nII;
689: } else {
690: IS nII;
691: ISDifference(pII,Pall,&nII);
692: ISDestroy(&pII);
693: pII = nII;
694: }
695: }
696: if (ploc) {
697: ISDifference(lPall,II,&lP);
698: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lP",(PetscObject)lP);
699: } else {
700: ISDifference(Pall,pII,&pP);
701: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",(PetscObject)pP);
702: /* need all local pressure dofs */
703: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
704: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
705: ISGetLocalSize(Pall,&ni);
706: ISGetIndices(Pall,&idxs);
707: for (i=0;i<ni;i++) matis->sf_rootdata[idxs[i]-rst] = 1;
708: ISRestoreIndices(Pall,&idxs);
709: PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
710: PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
711: for (i=0,ni=0;i<n;i++) if (matis->sf_leafdata[i]) widxs[ni++] = i;
712: ISCreateGeneral(PETSC_COMM_SELF,ni,widxs,PETSC_COPY_VALUES,&lPall);
713: }
715: if (!Pall) {
716: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
717: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
718: ISGetLocalSize(lPall,&ni);
719: ISGetIndices(lPall,&idxs);
720: for (i=0;i<ni;i++) matis->sf_leafdata[idxs[i]] = 1;
721: ISRestoreIndices(lPall,&idxs);
722: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
723: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
724: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
725: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&Pall);
726: }
727: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject)Pall);
729: if (flip) {
730: PetscInt npl;
731: ISGetLocalSize(Pall,&npl);
732: ISGetIndices(Pall,&idxs);
733: MatCreateVecs(A,NULL,&fetidp->rhs_flip);
734: VecSet(fetidp->rhs_flip,1.);
735: VecSetOption(fetidp->rhs_flip,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
736: for (i=0;i<npl;i++) {
737: VecSetValue(fetidp->rhs_flip,idxs[i],-1.,INSERT_VALUES);
738: }
739: VecAssemblyBegin(fetidp->rhs_flip);
740: VecAssemblyEnd(fetidp->rhs_flip);
741: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_flip",(PetscObject)fetidp->rhs_flip);
742: ISRestoreIndices(Pall,&idxs);
743: }
744: ISDestroy(&Pall);
745: ISDestroy(&pII);
747: /* local selected pressures in subdomain-wise and global ordering */
748: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
749: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
750: if (!ploc) {
751: PetscInt *widxs2;
753: if (!pP) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Missing parallel pressure IS");
754: ISGetLocalSize(pP,&ni);
755: ISGetIndices(pP,&idxs);
756: for (i=0;i<ni;i++) matis->sf_rootdata[idxs[i]-rst] = 1;
757: ISRestoreIndices(pP,&idxs);
758: PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
759: PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
760: for (i=0,ni=0;i<n;i++) if (matis->sf_leafdata[i]) widxs[ni++] = i;
761: PetscMalloc1(ni,&widxs2);
762: ISLocalToGlobalMappingApply(l2g,ni,widxs,widxs2);
763: ISCreateGeneral(PETSC_COMM_SELF,ni,widxs,PETSC_COPY_VALUES,&lP);
764: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lP",(PetscObject)lP);
765: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs2,PETSC_OWN_POINTER,&is1);
766: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_gP",(PetscObject)is1);
767: ISDestroy(&is1);
768: } else {
769: if (!lP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing sequential pressure IS");
770: ISGetLocalSize(lP,&ni);
771: ISGetIndices(lP,&idxs);
772: for (i=0;i<ni;i++)
773: if (idxs[i] >=0 && idxs[i] < n)
774: matis->sf_leafdata[idxs[i]] = 1;
775: ISRestoreIndices(lP,&idxs);
776: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
777: ISLocalToGlobalMappingApply(l2g,ni,idxs,widxs);
778: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&is1);
779: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_gP",(PetscObject)is1);
780: ISDestroy(&is1);
781: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
782: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
783: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&pP);
784: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",(PetscObject)pP);
785: }
786: PetscFree(widxs);
788: /* If there's any "interior pressure",
789: we may want to use a discrete harmonic solver instead
790: of a Stokes harmonic for the Dirichlet preconditioner
791: Need to extract the interior velocity dofs in interior dofs ordering (iV)
792: and interior pressure dofs in local ordering (iP) */
793: if (!allp) {
794: ISLocalToGlobalMapping l2g_t;
796: ISDifference(lPall,lP,&is1);
797: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_iP",(PetscObject)is1);
798: ISDifference(II,is1,&is2);
799: ISDestroy(&is1);
800: ISLocalToGlobalMappingCreateIS(II,&l2g_t);
801: ISGlobalToLocalMappingApplyIS(l2g_t,IS_GTOLM_DROP,is2,&is1);
802: ISGetLocalSize(is1,&i);
803: ISGetLocalSize(is2,&j);
804: if (i != j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local sizes %D and %D for iV",i,j);
805: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_iV",(PetscObject)is1);
806: ISLocalToGlobalMappingDestroy(&l2g_t);
807: ISDestroy(&is1);
808: ISDestroy(&is2);
809: }
810: ISDestroy(&II);
812: /* exclude selected pressures from the inner BDDC */
813: if (pcbddc->DirichletBoundariesLocal) {
814: IS list[2],plP,isout;
815: PetscInt np;
817: /* need a parallel IS */
818: ISGetLocalSize(lP,&np);
819: ISGetIndices(lP,&idxs);
820: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),np,idxs,PETSC_USE_POINTER,&plP);
821: list[0] = plP;
822: list[1] = pcbddc->DirichletBoundariesLocal;
823: ISConcatenate(PetscObjectComm((PetscObject)ksp),2,list,&isout);
824: ISSortRemoveDups(isout);
825: ISDestroy(&plP);
826: ISRestoreIndices(lP,&idxs);
827: PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc,isout);
828: ISDestroy(&isout);
829: } else if (pcbddc->DirichletBoundaries) {
830: IS list[2],isout;
832: list[0] = pP;
833: list[1] = pcbddc->DirichletBoundaries;
834: ISConcatenate(PetscObjectComm((PetscObject)ksp),2,list,&isout);
835: ISSortRemoveDups(isout);
836: PCBDDCSetDirichletBoundaries(fetidp->innerbddc,isout);
837: ISDestroy(&isout);
838: } else {
839: IS plP;
840: PetscInt np;
842: /* need a parallel IS */
843: ISGetLocalSize(lP,&np);
844: ISGetIndices(lP,&idxs);
845: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),np,idxs,PETSC_COPY_VALUES,&plP);
846: PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc,plP);
847: ISDestroy(&plP);
848: ISRestoreIndices(lP,&idxs);
849: }
851: /* save CSR information for the pressure BDDC solver (if any) */
852: if (schp) {
853: PetscInt np,nt;
855: MatGetSize(matis->A,&nt,NULL);
856: ISGetLocalSize(lP,&np);
857: if (np) {
858: PetscInt *xadj = pcbddc->mat_graph->xadj;
859: PetscInt *adjn = pcbddc->mat_graph->adjncy;
860: PetscInt nv = pcbddc->mat_graph->nvtxs_csr;
862: if (nv && nv == nt) {
863: ISLocalToGlobalMapping pmap;
864: PetscInt *schp_csr,*schp_xadj,*schp_adjn,p;
865: PetscContainer c;
867: ISLocalToGlobalMappingCreateIS(lPall,&pmap);
868: ISGetIndices(lPall,&idxs);
869: for (p = 0, nv = 0; p < np; p++) {
870: PetscInt x,n = idxs[p];
872: ISGlobalToLocalMappingApply(pmap,IS_GTOLM_DROP,xadj[n+1]-xadj[n],adjn+xadj[n],&x,NULL);
873: nv += x;
874: }
875: PetscMalloc1(np + 1 + nv,&schp_csr);
876: schp_xadj = schp_csr;
877: schp_adjn = schp_csr + np + 1;
878: for (p = 0, schp_xadj[0] = 0; p < np; p++) {
879: PetscInt x,n = idxs[p];
881: ISGlobalToLocalMappingApply(pmap,IS_GTOLM_DROP,xadj[n+1]-xadj[n],adjn+xadj[n],&x,schp_adjn + schp_xadj[p]);
882: schp_xadj[p+1] = schp_xadj[p] + x;
883: }
884: ISRestoreIndices(lPall,&idxs);
885: ISLocalToGlobalMappingDestroy(&pmap);
886: PetscContainerCreate(PETSC_COMM_SELF,&c);
887: PetscContainerSetPointer(c,schp_csr);
888: PetscContainerSetUserDestroy(c,PetscContainerUserDestroyDefault);
889: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pCSR",(PetscObject)c);
890: PetscContainerDestroy(&c);
892: }
893: }
894: }
895: ISDestroy(&lPall);
896: ISDestroy(&lP);
897: fetidp->pP = pP;
898: }
900: /* total number of selected pressure dofs */
901: ISGetSize(fetidp->pP,&totP);
903: /* Set operator for inner BDDC */
904: if (totP || fetidp->rhs_flip) {
905: MatDuplicate(A,MAT_COPY_VALUES,&nA);
906: } else {
907: PetscObjectReference((PetscObject)A);
908: nA = A;
909: }
910: if (fetidp->rhs_flip) {
911: MatDiagonalScale(nA,fetidp->rhs_flip,NULL);
912: if (totP) {
913: Mat lA2;
915: MatISGetLocalMat(nA,&lA);
916: MatDuplicate(lA,MAT_COPY_VALUES,&lA2);
917: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",(PetscObject)lA2);
918: MatDestroy(&lA2);
919: }
920: }
922: if (totP) {
923: MatSetOption(nA,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
924: MatZeroRowsColumnsIS(nA,fetidp->pP,1.,NULL,NULL);
925: } else {
926: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",NULL);
927: }
928: MatGetNearNullSpace(Ap,&nnsp);
929: if (!nnsp) {
930: MatGetNullSpace(Ap,&nnsp);
931: }
932: if (!nnsp) {
933: MatGetNearNullSpace(A,&nnsp);
934: }
935: if (!nnsp) {
936: MatGetNullSpace(A,&nnsp);
937: }
938: MatSetNearNullSpace(nA,nnsp);
939: PCSetOperators(fetidp->innerbddc,nA,nA);
940: MatDestroy(&nA);
942: /* non-zero rhs on interior dofs when applying the preconditioner */
943: if (totP) pcbddc->switch_static = PETSC_TRUE;
945: /* if there are no interface pressures, set inner bddc flag for benign saddle point */
946: if (!totP) {
947: pcbddc->benign_saddle_point = PETSC_TRUE;
948: pcbddc->compute_nonetflux = PETSC_TRUE;
949: }
951: /* Operators for pressure preconditioner */
952: if (totP) {
953: /* Extract pressure block if needed */
954: if (!pisz) {
955: Mat C;
956: IS nzrows = NULL;
958: MatCreateSubMatrix(A,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
959: MatFindNonzeroRows(C,&nzrows);
960: if (nzrows) {
961: PetscInt i;
963: ISGetSize(nzrows,&i);
964: ISDestroy(&nzrows);
965: if (!i) pisz = PETSC_TRUE;
966: }
967: if (!pisz) {
968: MatScale(C,-1.); /* i.e. Almost Incompressible Elasticity, Stokes discretized with Q1xQ1_stabilized */
969: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_C",(PetscObject)C);
970: }
971: MatDestroy(&C);
972: }
973: /* Divergence mat */
974: if (!pcbddc->divudotp) {
975: Mat B;
976: IS P;
977: IS l2l = NULL;
978: PetscBool save;
980: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&P);
981: if (!pisz) {
982: IS F,V;
983: PetscInt m,M;
985: MatGetOwnershipRange(A,&m,&M);
986: ISCreateStride(PetscObjectComm((PetscObject)A),M-m,m,1,&F);
987: ISComplement(P,m,M,&V);
988: MatCreateSubMatrix(A,P,V,MAT_INITIAL_MATRIX,&B);
989: {
990: Mat_IS *Bmatis = (Mat_IS*)B->data;
991: PetscObjectReference((PetscObject)Bmatis->getsub_cis);
992: l2l = Bmatis->getsub_cis;
993: }
994: ISDestroy(&V);
995: ISDestroy(&F);
996: } else {
997: MatCreateSubMatrix(A,P,NULL,MAT_INITIAL_MATRIX,&B);
998: }
999: save = pcbddc->compute_nonetflux; /* SetDivergenceMat activates nonetflux computation */
1000: PCBDDCSetDivergenceMat(fetidp->innerbddc,B,PETSC_FALSE,l2l);
1001: pcbddc->compute_nonetflux = save;
1002: MatDestroy(&B);
1003: ISDestroy(&l2l);
1004: }
1005: if (A != Ap) { /* user has provided a different Pmat, this always superseeds the setter (TODO: is it OK?) */
1006: /* use monolithic operator, we restrict later */
1007: KSPFETIDPSetPressureOperator(ksp,Ap);
1008: }
1009: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);
1011: /* PPmat not present, use some default choice */
1012: if (!PPmat) {
1013: Mat C;
1015: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_C",(PetscObject*)&C);
1016: if (!schp && C) { /* non-zero pressure block, most likely Almost Incompressible Elasticity */
1017: KSPFETIDPSetPressureOperator(ksp,C);
1018: } else if (!pisz && schp) { /* we need the whole pressure mass matrix to define the interface BDDC */
1019: IS P;
1021: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&P);
1022: MatCreateSubMatrix(A,P,P,MAT_INITIAL_MATRIX,&C);
1023: MatScale(C,-1.);
1024: KSPFETIDPSetPressureOperator(ksp,C);
1025: MatDestroy(&C);
1026: } else { /* identity (need to be scaled properly by the user using e.g. a Richardson method */
1027: PetscInt nl;
1029: ISGetLocalSize(fetidp->pP,&nl);
1030: MatCreate(PetscObjectComm((PetscObject)ksp),&C);
1031: MatSetSizes(C,nl,nl,totP,totP);
1032: MatSetType(C,MATAIJ);
1033: MatMPIAIJSetPreallocation(C,1,NULL,0,NULL);
1034: MatSeqAIJSetPreallocation(C,1,NULL);
1035: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1036: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1037: MatShift(C,1.);
1038: KSPFETIDPSetPressureOperator(ksp,C);
1039: MatDestroy(&C);
1040: }
1041: }
1043: /* Preconditioned operator for the pressure block */
1044: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);
1045: if (PPmat) {
1046: Mat C;
1047: IS Pall;
1048: PetscInt AM,PAM,PAN,pam,pan,am,an,pl,pIl,pAg,pIg;
1050: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&Pall);
1051: MatGetSize(A,&AM,NULL);
1052: MatGetSize(PPmat,&PAM,&PAN);
1053: ISGetSize(Pall,&pAg);
1054: ISGetSize(fetidp->pP,&pIg);
1055: MatGetLocalSize(PPmat,&pam,&pan);
1056: MatGetLocalSize(A,&am,&an);
1057: ISGetLocalSize(Pall,&pIl);
1058: ISGetLocalSize(fetidp->pP,&pl);
1059: if (PAM != PAN) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Pressure matrix must be square, unsupported %D x %D",PAM,PAN);
1060: if (pam != pan) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Local sizes of pressure matrix must be equal, unsupported %D x %D",pam,pan);
1061: if (pam != am && pam != pl && pam != pIl) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Invalid number of local rows %D for pressure matrix! Supported are %D, %D or %D",pam,am,pl,pIl);
1062: if (pan != an && pan != pl && pan != pIl) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_USER,"Invalid number of local columns %D for pressure matrix! Supported are %D, %D or %D",pan,an,pl,pIl);
1063: if (PAM == AM) { /* monolithic ordering, restrict to pressure */
1064: if (schp) {
1065: MatCreateSubMatrix(PPmat,Pall,Pall,MAT_INITIAL_MATRIX,&C);
1066: } else {
1067: MatCreateSubMatrix(PPmat,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
1068: }
1069: } else if (pAg == PAM) { /* global ordering for pressure only */
1070: if (!allp && !schp) { /* solving for interface pressure only */
1071: IS restr;
1073: ISRenumber(fetidp->pP,NULL,NULL,&restr);
1074: MatCreateSubMatrix(PPmat,restr,restr,MAT_INITIAL_MATRIX,&C);
1075: ISDestroy(&restr);
1076: } else {
1077: PetscObjectReference((PetscObject)PPmat);
1078: C = PPmat;
1079: }
1080: } else if (pIg == PAM) { /* global ordering for selected pressure only */
1081: if (schp) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Need the entire matrix");
1082: PetscObjectReference((PetscObject)PPmat);
1083: C = PPmat;
1084: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Unable to use the pressure matrix");
1086: KSPFETIDPSetPressureOperator(ksp,C);
1087: MatDestroy(&C);
1088: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Missing Pmat for pressure block");
1089: } else { /* totP == 0 */
1090: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",NULL);
1091: }
1092: }
1093: return(0);
1094: }
1096: static PetscErrorCode KSPSetUp_FETIDP(KSP ksp)
1097: {
1098: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1099: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1100: PetscBool flg;
1104: KSPFETIDPSetUpOperators(ksp);
1105: /* set up BDDC */
1106: PCSetErrorIfFailure(fetidp->innerbddc,ksp->errorifnotconverged);
1107: PCSetUp(fetidp->innerbddc);
1108: /* FETI-DP as it is implemented needs an exact coarse solver */
1109: if (pcbddc->coarse_ksp) {
1110: KSPSetTolerances(pcbddc->coarse_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,1000);
1111: KSPSetNormType(pcbddc->coarse_ksp,KSP_NORM_DEFAULT);
1112: }
1113: /* FETI-DP as it is implemented needs exact local Neumann solvers */
1114: KSPSetTolerances(pcbddc->ksp_R,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,1000);
1115: KSPSetNormType(pcbddc->ksp_R,KSP_NORM_DEFAULT);
1117: /* setup FETI-DP operators
1118: If fetidp->statechanged is true, we need to update the operators
1119: needed in the saddle-point case. This should be replaced
1120: by a better logic when the FETI-DP matrix and preconditioner will
1121: have their own classes */
1122: if (pcbddc->new_primal_space || fetidp->statechanged) {
1123: Mat F; /* the FETI-DP matrix */
1124: PC D; /* the FETI-DP preconditioner */
1125: KSPReset(fetidp->innerksp);
1126: PCBDDCCreateFETIDPOperators(fetidp->innerbddc,fetidp->fully_redundant,((PetscObject)ksp)->prefix,&F,&D);
1127: KSPSetOperators(fetidp->innerksp,F,F);
1128: KSPSetTolerances(fetidp->innerksp,ksp->rtol,ksp->abstol,ksp->divtol,ksp->max_it);
1129: KSPSetPC(fetidp->innerksp,D);
1130: PetscObjectIncrementTabLevel((PetscObject)D,(PetscObject)fetidp->innerksp,0);
1131: KSPSetFromOptions(fetidp->innerksp);
1132: MatCreateVecs(F,&(fetidp->innerksp)->vec_rhs,&(fetidp->innerksp)->vec_sol);
1133: MatDestroy(&F);
1134: PCDestroy(&D);
1135: if (fetidp->check) {
1136: PetscViewer viewer;
1138: if (!pcbddc->dbg_viewer) {
1139: viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
1140: } else {
1141: viewer = pcbddc->dbg_viewer;
1142: }
1143: KSPFETIDPCheckOperators(ksp,viewer);
1144: }
1145: }
1146: fetidp->statechanged = PETSC_FALSE;
1147: pcbddc->new_primal_space = PETSC_FALSE;
1149: /* propagate settings to the inner solve */
1150: KSPGetComputeSingularValues(ksp,&flg);
1151: KSPSetComputeSingularValues(fetidp->innerksp,flg);
1152: if (ksp->res_hist) {
1153: KSPSetResidualHistory(fetidp->innerksp,ksp->res_hist,ksp->res_hist_max,ksp->res_hist_reset);
1154: }
1155: KSPSetErrorIfNotConverged(fetidp->innerksp,ksp->errorifnotconverged);
1156: KSPSetUp(fetidp->innerksp);
1157: return(0);
1158: }
1160: static PetscErrorCode KSPSolve_FETIDP(KSP ksp)
1161: {
1163: Mat F,A;
1164: MatNullSpace nsp;
1165: Vec X,B,Xl,Bl;
1166: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1167: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1170: PetscCitationsRegister(citation,&cited);
1171: if (fetidp->saddlepoint) {
1172: PetscCitationsRegister(citation2,&cited2);
1173: }
1174: KSPGetOperators(ksp,&A,NULL);
1175: KSPGetRhs(ksp,&B);
1176: KSPGetSolution(ksp,&X);
1177: KSPGetOperators(fetidp->innerksp,&F,NULL);
1178: KSPGetRhs(fetidp->innerksp,&Bl);
1179: KSPGetSolution(fetidp->innerksp,&Xl);
1180: PCBDDCMatFETIDPGetRHS(F,B,Bl);
1181: if (ksp->transpose_solve) {
1182: KSPSolveTranspose(fetidp->innerksp,Bl,Xl);
1183: } else {
1184: KSPSolve(fetidp->innerksp,Bl,Xl);
1185: }
1186: PCBDDCMatFETIDPGetSolution(F,Xl,X);
1187: MatGetNullSpace(A,&nsp);
1188: if (nsp) {
1189: MatNullSpaceRemove(nsp,X);
1190: }
1191: /* update ksp with stats from inner ksp */
1192: KSPGetConvergedReason(fetidp->innerksp,&ksp->reason);
1193: KSPGetIterationNumber(fetidp->innerksp,&ksp->its);
1194: ksp->totalits += ksp->its;
1195: KSPGetResidualHistory(fetidp->innerksp,NULL,&ksp->res_hist_len);
1196: /* restore defaults for inner BDDC (Pre/PostSolve flags) */
1197: pcbddc->temp_solution_used = PETSC_FALSE;
1198: pcbddc->rhs_change = PETSC_FALSE;
1199: pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1200: return(0);
1201: }
1203: static PetscErrorCode KSPReset_FETIDP(KSP ksp)
1204: {
1205: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1206: PC_BDDC *pcbddc;
1210: ISDestroy(&fetidp->pP);
1211: VecDestroy(&fetidp->rhs_flip);
1212: /* avoid PCReset that does not take into account ref counting */
1213: PCDestroy(&fetidp->innerbddc);
1214: PCCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerbddc);
1215: PCSetType(fetidp->innerbddc,PCBDDC);
1216: pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1217: pcbddc->symmetric_primal = PETSC_FALSE;
1218: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerbddc);
1219: KSPDestroy(&fetidp->innerksp);
1220: fetidp->saddlepoint = PETSC_FALSE;
1221: fetidp->matstate = -1;
1222: fetidp->matnnzstate = -1;
1223: fetidp->statechanged = PETSC_TRUE;
1224: return(0);
1225: }
1227: static PetscErrorCode KSPDestroy_FETIDP(KSP ksp)
1228: {
1229: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1233: KSPReset_FETIDP(ksp);
1234: PCDestroy(&fetidp->innerbddc);
1235: KSPDestroy(&fetidp->innerksp);
1236: PetscFree(fetidp->monctx);
1237: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetInnerBDDC_C",NULL);
1238: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerBDDC_C",NULL);
1239: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerKSP_C",NULL);
1240: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetPressureOperator_C",NULL);
1241: PetscFree(ksp->data);
1242: return(0);
1243: }
1245: static PetscErrorCode KSPView_FETIDP(KSP ksp,PetscViewer viewer)
1246: {
1247: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1249: PetscBool iascii;
1252: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1253: if (iascii) {
1254: PetscViewerASCIIPrintf(viewer," fully redundant: %d\n",fetidp->fully_redundant);
1255: PetscViewerASCIIPrintf(viewer," saddle point: %d\n",fetidp->saddlepoint);
1256: PetscViewerASCIIPrintf(viewer,"Inner KSP solver details\n");
1257: }
1258: PetscViewerASCIIPushTab(viewer);
1259: KSPView(fetidp->innerksp,viewer);
1260: PetscViewerASCIIPopTab(viewer);
1261: if (iascii) {
1262: PetscViewerASCIIPrintf(viewer,"Inner BDDC solver details\n");
1263: }
1264: PetscViewerASCIIPushTab(viewer);
1265: PCView(fetidp->innerbddc,viewer);
1266: PetscViewerASCIIPopTab(viewer);
1267: return(0);
1268: }
1270: static PetscErrorCode KSPSetFromOptions_FETIDP(PetscOptionItems *PetscOptionsObject,KSP ksp)
1271: {
1272: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1276: /* set options prefixes for the inner objects, since the parent prefix will be valid at this point */
1277: PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerksp,((PetscObject)ksp)->prefix);
1278: PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerksp,"fetidp_");
1279: if (!fetidp->userbddc) {
1280: PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerbddc,((PetscObject)ksp)->prefix);
1281: PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerbddc,"fetidp_bddc_");
1282: }
1283: PetscOptionsHead(PetscOptionsObject,"KSP FETIDP options");
1284: PetscOptionsBool("-ksp_fetidp_fullyredundant","Use fully redundant multipliers","none",fetidp->fully_redundant,&fetidp->fully_redundant,NULL);
1285: PetscOptionsBool("-ksp_fetidp_saddlepoint","Activates support for saddle-point problems",NULL,fetidp->saddlepoint,&fetidp->saddlepoint,NULL);
1286: PetscOptionsBool("-ksp_fetidp_check","Activates verbose debugging output FETI-DP operators",NULL,fetidp->check,&fetidp->check,NULL);
1287: PetscOptionsTail();
1288: PCSetFromOptions(fetidp->innerbddc);
1289: return(0);
1290: }
1292: /*MC
1293: KSPFETIDP - The FETI-DP method
1295: This class implements the FETI-DP method [1].
1296: The matrix for the KSP must be of type MATIS.
1297: The FETI-DP linear system (automatically generated constructing an internal PCBDDC object) is solved using an internal KSP object.
1299: Options Database Keys:
1300: + -ksp_fetidp_fullyredundant <false> : use a fully redundant set of Lagrange multipliers
1301: . -ksp_fetidp_saddlepoint <false> : activates support for saddle point problems, see [2]
1302: . -ksp_fetidp_saddlepoint_flip <false> : usually, an incompressible Stokes problem is written as
1303: | A B^T | | v | = | f |
1304: | B 0 | | p | = | g |
1305: with B representing -\int_\Omega \nabla \cdot u q.
1306: If -ksp_fetidp_saddlepoint_flip is true, the code assumes that the user provides it as
1307: | A B^T | | v | = | f |
1308: |-B 0 | | p | = |-g |
1309: . -ksp_fetidp_pressure_field <-1> : activates support for saddle point problems, and identifies the pressure field id.
1310: If this information is not provided, the pressure field is detected by using MatFindZeroDiagonals().
1311: - -ksp_fetidp_pressure_all <false> : if false, uses the interface pressures, as described in [2]. If true, uses the entire pressure field.
1313: Level: Advanced
1315: Notes:
1316: 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.,
1317: .vb
1318: -fetidp_ksp_type gmres -fetidp_bddc_pc_bddc_symmetric false
1319: .ve
1320: will use GMRES for the solution of the linear system on the Lagrange multipliers, generated using a non-symmetric PCBDDC.
1322: For saddle point problems with continuous pressures, the preconditioned operator for the pressure solver can be specified with KSPFETIDPSetPressureOperator().
1323: Alternatively, the pressure operator is extracted from the precondioned matrix (if it is different from the linear solver matrix).
1324: If none of the above, an identity matrix will be created; the user then needs to scale it through a Richardson solver.
1325: Options for the pressure solver can be prefixed with -fetidp_fielsplit_p_, E.g.
1326: .vb
1327: -fetidp_fielsplit_p_ksp_type preonly -fetidp_fielsplit_p_pc_type lu -fetidp_fielsplit_p_pc_factor_mat_solver_type mumps
1328: .ve
1329: 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
1330: non-redundant multipliers, i.e. -ksp_fetidp_fullyredundant false. Options for the scaling solver are prefixed by -fetidp_bddelta_, E.g.
1331: .vb
1332: -fetidp_bddelta_pc_factor_mat_solver_type mumps -fetidp_bddelta_pc_type lu
1333: .ve
1335: 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.
1337: The converged reason and number of iterations computed are passed from the inner KSP to this KSP at the end of the solution.
1339: Developer Notes:
1340: Even though this method does not directly use any norms, the user is allowed to set the KSPNormType to any value.
1341: This is so users do not have to change KSPNormType options when they switch from other KSP methods to this one.
1343: References:
1344: .vb
1345: . [1] - 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
1346: . [2] - 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
1347: .ve
1349: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC(), KSPFETIDPGetInnerBDDC(), KSPFETIDPGetInnerKSP()
1350: M*/
1351: PETSC_EXTERN PetscErrorCode KSPCreate_FETIDP(KSP ksp)
1352: {
1354: KSP_FETIDP *fetidp;
1355: KSP_FETIDPMon *monctx;
1356: PC_BDDC *pcbddc;
1357: PC pc;
1360: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
1361: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
1362: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
1363: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
1364: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
1365: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
1366: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
1368: PetscNewLog(ksp,&fetidp);
1369: fetidp->matstate = -1;
1370: fetidp->matnnzstate = -1;
1371: fetidp->statechanged = PETSC_TRUE;
1373: ksp->data = (void*)fetidp;
1374: ksp->ops->setup = KSPSetUp_FETIDP;
1375: ksp->ops->solve = KSPSolve_FETIDP;
1376: ksp->ops->destroy = KSPDestroy_FETIDP;
1377: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_FETIDP;
1378: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_FETIDP;
1379: ksp->ops->view = KSPView_FETIDP;
1380: ksp->ops->setfromoptions = KSPSetFromOptions_FETIDP;
1381: ksp->ops->buildsolution = KSPBuildSolution_FETIDP;
1382: ksp->ops->buildresidual = KSPBuildResidualDefault;
1383: KSPGetPC(ksp,&pc);
1384: PCSetType(pc,PCNONE);
1385: /* create the inner KSP for the Lagrange multipliers */
1386: KSPCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerksp);
1387: KSPGetPC(fetidp->innerksp,&pc);
1388: PCSetType(pc,PCNONE);
1389: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerksp);
1390: /* monitor */
1391: PetscNew(&monctx);
1392: monctx->parentksp = ksp;
1393: fetidp->monctx = monctx;
1394: KSPMonitorSet(fetidp->innerksp,KSPMonitor_FETIDP,fetidp->monctx,NULL);
1395: /* create the inner BDDC */
1396: PCCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerbddc);
1397: PCSetType(fetidp->innerbddc,PCBDDC);
1398: /* make sure we always obtain a consistent FETI-DP matrix
1399: for symmetric problems, the user can always customize it through the command line */
1400: pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1401: pcbddc->symmetric_primal = PETSC_FALSE;
1402: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerbddc);
1403: /* composed functions */
1404: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetInnerBDDC_C",KSPFETIDPSetInnerBDDC_FETIDP);
1405: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerBDDC_C",KSPFETIDPGetInnerBDDC_FETIDP);
1406: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerKSP_C",KSPFETIDPGetInnerKSP_FETIDP);
1407: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetPressureOperator_C",KSPFETIDPSetPressureOperator_FETIDP);
1408: /* need to call KSPSetUp_FETIDP even with KSP_SETUP_NEWMATRIX */
1409: ksp->setupnewmatrix = PETSC_TRUE;
1410: return(0);
1411: }