Actual source code: fetidp.c
petsc-3.9.3 2018-07-02
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: The operator can be either passed in a) monolithic global ordering, b) pressure-only global ordering
90: or c) interface pressure ordering (if -ksp_fetidp_pressure_all false).
91: In cases b) and c), the pressure ordering of dofs needs to satisfy
92: pid_1 < pid_2 iff gid_1 < gid_2
93: where pid_1 and pid_2 are two different pressure dof numbers and gid_1 and gid_2 the corresponding
94: id in the monolithic global ordering.
96: .seealso: MATIS, PCBDDC, KSPFETIDPGetInnerBDDC, KSPFETIDPGetInnerKSP, KSPSetOperators
97: @*/
98: PetscErrorCode KSPFETIDPSetPressureOperator(KSP ksp, Mat P)
99: {
105: PetscTryMethod(ksp,"KSPFETIDPSetPressureOperator_C",(KSP,Mat),(ksp,P));
106: return(0);
107: }
109: static PetscErrorCode KSPFETIDPGetInnerKSP_FETIDP(KSP ksp, KSP* innerksp)
110: {
111: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
114: *innerksp = fetidp->innerksp;
115: return(0);
116: }
118: /*@
119: KSPFETIDPGetInnerKSP - Gets the KSP object for the Lagrange multipliers
121: Input Parameters:
122: + ksp - the FETI-DP KSP
123: - innerksp - the KSP for the multipliers
125: Level: advanced
127: Notes:
129: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC, KSPFETIDPGetInnerBDDC
130: @*/
131: PetscErrorCode KSPFETIDPGetInnerKSP(KSP ksp, KSP* innerksp)
132: {
138: PetscUseMethod(ksp,"KSPFETIDPGetInnerKSP_C",(KSP,KSP*),(ksp,innerksp));
139: return(0);
140: }
142: static PetscErrorCode KSPFETIDPGetInnerBDDC_FETIDP(KSP ksp, PC* pc)
143: {
144: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
147: *pc = fetidp->innerbddc;
148: return(0);
149: }
151: /*@
152: KSPFETIDPGetInnerBDDC - Gets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers
154: Input Parameters:
155: + ksp - the FETI-DP Krylov solver
156: - pc - the BDDC preconditioner
158: Level: advanced
160: Notes:
162: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC, KSPFETIDPGetInnerKSP
163: @*/
164: PetscErrorCode KSPFETIDPGetInnerBDDC(KSP ksp, PC* pc)
165: {
171: PetscUseMethod(ksp,"KSPFETIDPGetInnerBDDC_C",(KSP,PC*),(ksp,pc));
172: return(0);
173: }
175: static PetscErrorCode KSPFETIDPSetInnerBDDC_FETIDP(KSP ksp, PC pc)
176: {
177: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
181: PetscObjectReference((PetscObject)pc);
182: PCDestroy(&fetidp->innerbddc);
183: fetidp->innerbddc = pc;
184: fetidp->userbddc = PETSC_TRUE;
185: return(0);
186: }
188: /*@
189: KSPFETIDPSetInnerBDDC - Sets the BDDC preconditioner used to setup the FETI-DP matrix for the Lagrange multipliers
191: Collective on KSP
193: Input Parameters:
194: + ksp - the FETI-DP Krylov solver
195: - pc - the BDDC preconditioner
197: Level: advanced
199: Notes:
201: .seealso: MATIS, PCBDDC, KSPFETIDPGetInnerBDDC, KSPFETIDPGetInnerKSP
202: @*/
203: PetscErrorCode KSPFETIDPSetInnerBDDC(KSP ksp, PC pc)
204: {
205: PetscBool isbddc;
211: PetscObjectTypeCompare((PetscObject)pc,PCBDDC,&isbddc);
212: if (!isbddc) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"KSPFETIDPSetInnerBDDC need a PCBDDC preconditioner");
213: PetscTryMethod(ksp,"KSPFETIDPSetInnerBDDC_C",(KSP,PC),(ksp,pc));
214: return(0);
215: }
217: static PetscErrorCode KSPBuildSolution_FETIDP(KSP ksp,Vec v,Vec *V)
218: {
219: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
220: Mat F;
221: Vec Xl;
225: KSPGetOperators(fetidp->innerksp,&F,NULL);
226: KSPBuildSolution(fetidp->innerksp,NULL,&Xl);
227: if (v) {
228: PCBDDCMatFETIDPGetSolution(F,Xl,v);
229: *V = v;
230: } else {
231: PCBDDCMatFETIDPGetSolution(F,Xl,*V);
232: }
233: return(0);
234: }
236: static PetscErrorCode KSPMonitor_FETIDP(KSP ksp,PetscInt it,PetscReal rnorm,void* ctx)
237: {
238: KSP_FETIDPMon *monctx = (KSP_FETIDPMon*)ctx;
242: KSPMonitor(monctx->parentksp,it,rnorm);
243: return(0);
244: }
246: static PetscErrorCode KSPComputeEigenvalues_FETIDP(KSP ksp,PetscInt nmax,PetscReal *r,PetscReal *c,PetscInt *neig)
247: {
248: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
252: KSPComputeEigenvalues(fetidp->innerksp,nmax,r,c,neig);
253: return(0);
254: }
256: static PetscErrorCode KSPComputeExtremeSingularValues_FETIDP(KSP ksp,PetscReal *emax,PetscReal *emin)
257: {
258: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
262: KSPComputeExtremeSingularValues(fetidp->innerksp,emax,emin);
263: return(0);
264: }
266: static PetscErrorCode KSPFETIDPCheckOperators(KSP ksp, PetscViewer viewer)
267: {
268: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
269: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
270: PC_IS *pcis = (PC_IS*)fetidp->innerbddc->data;
271: Mat_IS *matis = (Mat_IS*)fetidp->innerbddc->pmat->data;
272: Mat F;
273: FETIDPMat_ctx fetidpmat_ctx;
274: Vec test_vec,test_vec_p = NULL,fetidp_global;
275: IS dirdofs,isvert;
276: MPI_Comm comm = PetscObjectComm((PetscObject)ksp);
277: PetscScalar sval,*array;
278: PetscReal val,rval;
279: const PetscInt *vertex_indices;
280: PetscInt i,n_vertices;
281: PetscBool isascii;
286: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
287: if (!isascii) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported viewer");
288: PetscViewerASCIIPrintf(viewer,"----------FETI-DP MAT --------------\n");
289: PetscViewerASCIIAddTab(viewer,2);
290: KSPGetOperators(fetidp->innerksp,&F,NULL);
291: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
292: MatView(F,viewer);
293: PetscViewerPopFormat(viewer);
294: PetscViewerASCIISubtractTab(viewer,2);
295: MatShellGetContext(F,(void**)&fetidpmat_ctx);
296: PetscViewerASCIIPrintf(viewer,"----------FETI-DP TESTS--------------\n");
297: PetscViewerASCIIPrintf(viewer,"All tests should return zero!\n");
298: PetscViewerASCIIPrintf(viewer,"FETIDP MAT context in the ");
299: if (fetidp->fully_redundant) {
300: PetscViewerASCIIPrintf(viewer,"fully redundant case for lagrange multipliers.\n");
301: } else {
302: PetscViewerASCIIPrintf(viewer,"Non-fully redundant case for lagrange multiplier.\n");
303: }
304: PetscViewerFlush(viewer);
306: /* Get Vertices used to define the BDDC */
307: PCBDDCGraphGetCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);
308: ISGetLocalSize(isvert,&n_vertices);
309: ISGetIndices(isvert,&vertex_indices);
311: /******************************************************************/
312: /* TEST A/B: Test numbering of global fetidp dofs */
313: /******************************************************************/
314: MatCreateVecs(F,&fetidp_global,NULL);
315: VecDuplicate(fetidpmat_ctx->lambda_local,&test_vec);
316: VecSet(fetidp_global,1.0);
317: VecSet(test_vec,1.);
318: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
319: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
320: if (fetidpmat_ctx->l2g_p) {
321: VecDuplicate(fetidpmat_ctx->vP,&test_vec_p);
322: VecSet(test_vec_p,1.);
323: VecScatterBegin(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
324: VecScatterEnd(fetidpmat_ctx->l2g_p,fetidp_global,fetidpmat_ctx->vP,INSERT_VALUES,SCATTER_REVERSE);
325: }
326: VecAXPY(test_vec,-1.0,fetidpmat_ctx->lambda_local);
327: VecNorm(test_vec,NORM_INFINITY,&val);
328: VecDestroy(&test_vec);
329: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
330: PetscViewerASCIIPrintf(viewer,"A: CHECK glob to loc: % 1.14e\n",rval);
332: if (fetidpmat_ctx->l2g_p) {
333: VecAXPY(test_vec_p,-1.0,fetidpmat_ctx->vP);
334: VecNorm(test_vec_p,NORM_INFINITY,&val);
335: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
336: PetscViewerASCIIPrintf(viewer,"A: CHECK glob to loc (p): % 1.14e\n",rval);
337: }
339: if (fetidp->fully_redundant) {
340: VecSet(fetidp_global,0.0);
341: VecSet(fetidpmat_ctx->lambda_local,0.5);
342: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
343: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
344: VecSum(fetidp_global,&sval);
345: val = PetscRealPart(sval)-fetidpmat_ctx->n_lambda;
346: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
347: PetscViewerASCIIPrintf(viewer,"B: CHECK loc to glob: % 1.14e\n",rval);
348: }
350: if (fetidpmat_ctx->l2g_p) {
351: VecSet(pcis->vec1_N,1.0);
352: VecSet(pcis->vec1_global,0.0);
353: VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
354: VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
356: VecSet(fetidp_global,0.0);
357: VecSet(fetidpmat_ctx->vP,-1.0);
358: VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
359: VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
360: VecScatterBegin(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
361: VecScatterEnd(fetidpmat_ctx->g2g_p,fetidp_global,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
362: VecScatterBegin(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
363: VecScatterEnd(fetidpmat_ctx->g2g_p,pcis->vec1_global,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
364: VecSum(fetidp_global,&sval);
365: val = PetscRealPart(sval);
366: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
367: PetscViewerASCIIPrintf(viewer,"B: CHECK loc to glob (p): % 1.14e\n",rval);
368: }
370: /******************************************************************/
371: /* TEST C: It should hold B_delta*w=0, w\in\widehat{W} */
372: /* This is the meaning of the B matrix */
373: /******************************************************************/
375: VecSetRandom(pcis->vec1_N,NULL);
376: VecSet(pcis->vec1_global,0.0);
377: VecScatterBegin(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
378: VecScatterEnd(matis->rctx,pcis->vec1_N,pcis->vec1_global,ADD_VALUES,SCATTER_REVERSE);
379: VecScatterBegin(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
380: VecScatterEnd(matis->rctx,pcis->vec1_global,pcis->vec1_N,INSERT_VALUES,SCATTER_FORWARD);
381: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
382: VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
383: /* Action of B_delta */
384: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
385: VecSet(fetidp_global,0.0);
386: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
387: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
388: VecNorm(fetidp_global,NORM_INFINITY,&val);
389: PetscViewerASCIIPrintf(viewer,"C: CHECK infty norm of B_delta*w (w continuous): % 1.14e\n",val);
391: /******************************************************************/
392: /* TEST D: It should hold E_Dw = w - P_Dw w\in\widetilde{W} */
393: /* E_D = R_D^TR */
394: /* P_D = B_{D,delta}^T B_{delta} */
395: /* eq.44 Mandel Tezaur and Dohrmann 2005 */
396: /******************************************************************/
398: /* compute a random vector in \widetilde{W} */
399: VecSetRandom(pcis->vec1_N,NULL);
400: /* set zero at vertices and essential dofs */
401: VecGetArray(pcis->vec1_N,&array);
402: for (i=0;i<n_vertices;i++) array[vertex_indices[i]] = 0.0;
403: PCBDDCGraphGetDirichletDofs(pcbddc->mat_graph,&dirdofs);
404: if (dirdofs) {
405: const PetscInt *idxs;
406: PetscInt ndir;
408: ISGetLocalSize(dirdofs,&ndir);
409: ISGetIndices(dirdofs,&idxs);
410: for (i=0;i<ndir;i++) array[idxs[i]] = 0.0;
411: ISRestoreIndices(dirdofs,&idxs);
412: }
413: VecRestoreArray(pcis->vec1_N,&array);
414: /* store w for final comparison */
415: VecDuplicate(pcis->vec1_B,&test_vec);
416: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
417: VecScatterEnd(pcis->N_to_B,pcis->vec1_N,test_vec,INSERT_VALUES,SCATTER_FORWARD);
419: /* Jump operator P_D : results stored in pcis->vec1_B */
420: /* Action of B_delta */
421: MatMult(fetidpmat_ctx->B_delta,test_vec,fetidpmat_ctx->lambda_local);
422: VecSet(fetidp_global,0.0);
423: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
424: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
425: /* Action of B_Ddelta^T */
426: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
427: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
428: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
430: /* Average operator E_D : results stored in pcis->vec2_B */
431: PCBDDCScalingExtension(fetidpmat_ctx->pc,test_vec,pcis->vec1_global);
432: VecScatterBegin(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
433: VecScatterEnd(pcis->global_to_B,pcis->vec1_global,pcis->vec2_B,INSERT_VALUES,SCATTER_FORWARD);
435: /* test E_D=I-P_D */
436: VecAXPY(pcis->vec1_B,1.0,pcis->vec2_B);
437: VecAXPY(pcis->vec1_B,-1.0,test_vec);
438: VecNorm(pcis->vec1_B,NORM_INFINITY,&val);
439: VecDestroy(&test_vec);
440: MPI_Reduce(&val,&rval,1,MPIU_REAL,MPI_MAX,0,comm);
441: PetscViewerASCIIPrintf(viewer,"D: CHECK infty norm of E_D + P_D - I: % 1.14e\n",PetscGlobalRank,val);
443: /******************************************************************/
444: /* TEST E: It should hold R_D^TP_Dw=0 w\in\widetilde{W} */
445: /* eq.48 Mandel Tezaur and Dohrmann 2005 */
446: /******************************************************************/
448: VecSetRandom(pcis->vec1_N,NULL);
449: /* set zero at vertices and essential dofs */
450: VecGetArray(pcis->vec1_N,&array);
451: for (i=0;i<n_vertices;i++) array[vertex_indices[i]] = 0.0;
452: if (dirdofs) {
453: const PetscInt *idxs;
454: PetscInt ndir;
456: ISGetLocalSize(dirdofs,&ndir);
457: ISGetIndices(dirdofs,&idxs);
458: for (i=0;i<ndir;i++) array[idxs[i]] = 0.0;
459: ISRestoreIndices(dirdofs,&idxs);
460: }
461: VecRestoreArray(pcis->vec1_N,&array);
463: /* Jump operator P_D : results stored in pcis->vec1_B */
465: VecScatterBegin(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
466: VecScatterEnd(pcis->N_to_B,pcis->vec1_N,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
467: /* Action of B_delta */
468: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
469: VecSet(fetidp_global,0.0);
470: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
471: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,fetidp_global,ADD_VALUES,SCATTER_FORWARD);
472: /* Action of B_Ddelta^T */
473: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
474: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
475: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
476: /* scaling */
477: PCBDDCScalingExtension(fetidpmat_ctx->pc,pcis->vec1_B,pcis->vec1_global);
478: VecNorm(pcis->vec1_global,NORM_INFINITY,&val);
479: PetscViewerASCIIPrintf(viewer,"E: CHECK infty norm of R^T_D P_D: % 1.14e\n",val);
481: if (!fetidp->fully_redundant) {
482: /******************************************************************/
483: /* TEST F: It should holds B_{delta}B^T_{D,delta}=I */
484: /* Corollary thm 14 Mandel Tezaur and Dohrmann 2005 */
485: /******************************************************************/
486: VecDuplicate(fetidp_global,&test_vec);
487: VecSetRandom(fetidp_global,NULL);
488: if (fetidpmat_ctx->l2g_p) {
489: VecSet(fetidpmat_ctx->vP,0.);
490: VecScatterBegin(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
491: VecScatterEnd(fetidpmat_ctx->l2g_p,fetidpmat_ctx->vP,fetidp_global,INSERT_VALUES,SCATTER_FORWARD);
492: }
493: /* Action of B_Ddelta^T */
494: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
495: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidp_global,fetidpmat_ctx->lambda_local,INSERT_VALUES,SCATTER_REVERSE);
496: MatMultTranspose(fetidpmat_ctx->B_Ddelta,fetidpmat_ctx->lambda_local,pcis->vec1_B);
497: /* Action of B_delta */
498: MatMult(fetidpmat_ctx->B_delta,pcis->vec1_B,fetidpmat_ctx->lambda_local);
499: VecSet(test_vec,0.0);
500: VecScatterBegin(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
501: VecScatterEnd(fetidpmat_ctx->l2g_lambda,fetidpmat_ctx->lambda_local,test_vec,ADD_VALUES,SCATTER_FORWARD);
502: VecAXPY(fetidp_global,-1.,test_vec);
503: VecNorm(fetidp_global,NORM_INFINITY,&val);
504: PetscViewerASCIIPrintf(viewer,"E: CHECK infty norm of P^T_D - I: % 1.14e\n",val);
505: VecDestroy(&test_vec);
506: }
507: PetscViewerASCIIPrintf(viewer,"-------------------------------------\n");
508: PetscViewerFlush(viewer);
509: VecDestroy(&test_vec_p);
510: ISDestroy(&dirdofs);
511: VecDestroy(&fetidp_global);
512: ISRestoreIndices(isvert,&vertex_indices);
513: PCBDDCGraphRestoreCandidatesIS(pcbddc->mat_graph,NULL,NULL,NULL,NULL,&isvert);
514: return(0);
515: }
517: static PetscErrorCode KSPFETIDPSetUpOperators(KSP ksp)
518: {
519: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
520: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
521: Mat A,Ap;
522: PetscInt fid = -1;
523: PetscMPIInt size;
524: PetscBool ismatis,pisz,allp;
525: PetscBool flip; /* Usually, Stokes is written (B = -\int_\Omega \nabla \cdot u q)
526: | A B'| | v | = | f |
527: | B 0 | | p | = | g |
528: If -ksp_fetidp_saddlepoint_flip is true, the code assumes it is written as
529: | A B'| | v | = | f |
530: |-B 0 | | p | = |-g |
531: */
532: PetscObjectState matstate, matnnzstate;
533: PetscErrorCode ierr;
536: pisz = PETSC_FALSE;
537: flip = PETSC_FALSE;
538: allp = PETSC_FALSE;
539: PetscOptionsBegin(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"FETI-DP options","PC");
540: PetscOptionsInt("-ksp_fetidp_pressure_field","Field id for pressures for saddle-point problems",NULL,fid,&fid,NULL);
541: PetscOptionsBool("-ksp_fetidp_pressure_all","Use the whole pressure set instead of just that at the interface",NULL,allp,&allp,NULL);
542: PetscOptionsBool("-ksp_fetidp_saddlepoint_flip","Flip the sign of the pressure-velocity (lower-left) block",NULL,flip,&flip,NULL);
543: PetscOptionsEnd();
545: MPI_Comm_size(PetscObjectComm((PetscObject)ksp),&size);
546: fetidp->saddlepoint = (fid >= 0 ? PETSC_TRUE : fetidp->saddlepoint);
547: if (size == 1) fetidp->saddlepoint = PETSC_FALSE;
549: KSPGetOperators(ksp,&A,&Ap);
550: PetscObjectTypeCompare((PetscObject)A,MATIS,&ismatis);
551: if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Amat should be of type MATIS");
553: /* Quiet return if the matrix states are unchanged.
554: Needed only for the saddle point case since it uses MatZeroRows
555: on a matrix that may not have changed */
556: PetscObjectStateGet((PetscObject)A,&matstate);
557: MatGetNonzeroState(A,&matnnzstate);
558: if (matstate == fetidp->matstate && matnnzstate == fetidp->matnnzstate) return(0);
559: fetidp->matstate = matstate;
560: fetidp->matnnzstate = matnnzstate;
561: fetidp->statechanged = fetidp->saddlepoint;
563: /* see if we have some fields attached */
564: if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
565: DM dm;
566: PetscContainer c;
568: KSPGetDM(ksp,&dm);
569: PetscObjectQuery((PetscObject)A,"_convert_nest_lfields",(PetscObject*)&c);
570: if (dm) {
571: IS *fields;
572: PetscInt nf,i;
574: DMCreateFieldDecomposition(dm,&nf,NULL,&fields,NULL);
575: PCBDDCSetDofsSplitting(fetidp->innerbddc,nf,fields);
576: for (i=0;i<nf;i++) {
577: ISDestroy(&fields[i]);
578: }
579: PetscFree(fields);
580: } else if (c) {
581: MatISLocalFields lf;
582: PetscContainerGetPointer(c,(void**)&lf);
583: PCBDDCSetDofsSplittingLocal(fetidp->innerbddc,lf->nr,lf->rf);
584: }
585: }
587: if (!fetidp->saddlepoint) {
588: PCSetOperators(fetidp->innerbddc,A,A);
589: } else {
590: Mat nA,lA;
591: Mat PPmat;
592: IS pP;
593: PetscInt totP;
595: MatISGetLocalMat(A,&lA);
596: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",(PetscObject)lA);
598: pP = fetidp->pP;
599: if (!pP) { /* first time, need to compute pressure dofs */
600: PC_IS *pcis = (PC_IS*)fetidp->innerbddc->data;
601: Mat_IS *matis = (Mat_IS*)(A->data);
602: ISLocalToGlobalMapping l2g;
603: IS lP = NULL,II,pII,lPall,Pall,is1,is2;
604: const PetscInt *idxs;
605: PetscInt nl,ni,*widxs;
606: PetscInt i,j,n_neigh,*neigh,*n_shared,**shared,*count;
607: PetscInt rst,ren,n;
608: PetscBool ploc;
610: MatGetLocalSize(A,&nl,NULL);
611: MatGetOwnershipRange(A,&rst,&ren);
612: MatGetLocalSize(lA,&n,NULL);
613: MatGetLocalToGlobalMapping(A,&l2g,NULL);
615: if (!pcis->is_I_local) { /* need to compute interior dofs */
616: PetscCalloc1(n,&count);
617: ISLocalToGlobalMappingGetInfo(l2g,&n_neigh,&neigh,&n_shared,&shared);
618: for (i=1;i<n_neigh;i++)
619: for (j=0;j<n_shared[i];j++)
620: count[shared[i][j]] += 1;
621: for (i=0,j=0;i<n;i++) if (!count[i]) count[j++] = i;
622: ISLocalToGlobalMappingRestoreInfo(l2g,&n_neigh,&neigh,&n_shared,&shared);
623: ISCreateGeneral(PETSC_COMM_SELF,j,count,PETSC_OWN_POINTER,&II);
624: } else {
625: PetscObjectReference((PetscObject)pcis->is_I_local);
626: II = pcis->is_I_local;
627: }
629: /* interior dofs in layout */
630: MatISSetUpSF(A);
631: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
632: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
633: ISGetLocalSize(II,&ni);
634: ISGetIndices(II,&idxs);
635: for (i=0;i<ni;i++) matis->sf_leafdata[idxs[i]] = 1;
636: ISRestoreIndices(II,&idxs);
637: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
638: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
639: PetscMalloc1(PetscMax(nl,n),&widxs);
640: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
641: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&pII);
643: /* pressure dofs */
644: Pall = NULL;
645: lPall = NULL;
646: ploc = PETSC_FALSE;
647: if (fid >= 0) {
648: if (pcbddc->n_ISForDofsLocal) {
649: PetscInt np;
651: if (fid >= pcbddc->n_ISForDofsLocal) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Invalid field id for pressure %D, max %D",fid,pcbddc->n_ISForDofsLocal);
652: /* need a sequential IS */
653: ISGetLocalSize(pcbddc->ISForDofsLocal[fid],&np);
654: ISGetIndices(pcbddc->ISForDofsLocal[fid],&idxs);
655: ISCreateGeneral(PETSC_COMM_SELF,np,idxs,PETSC_COPY_VALUES,&lPall);
656: ISRestoreIndices(pcbddc->ISForDofsLocal[fid],&idxs);
657: ploc = PETSC_TRUE;
658: } else if (pcbddc->n_ISForDofs) {
659: if (fid >= pcbddc->n_ISForDofs) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Invalid field id for pressure %D, max %D",fid,pcbddc->n_ISForDofs);
660: PetscObjectReference((PetscObject)pcbddc->ISForDofs[fid]);
661: Pall = pcbddc->ISForDofs[fid];
662: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Missing fields! Use PCBDDCSetDofsSplitting/Local");
663: } else { /* fallback to zero pressure block */
664: IS list[2];
666: MatFindZeroDiagonals(A,&list[1]);
667: ISComplement(list[1],rst,ren,&list[0]);
668: PCBDDCSetDofsSplitting(fetidp->innerbddc,2,list);
669: ISDestroy(&list[0]);
670: Pall = list[1];
671: }
672: /* if the user requested the entire pressure,
673: remove the interior pressure dofs from II (or pII) */
674: if (allp) {
675: if (ploc) {
676: IS nII;
677: ISDifference(II,lPall,&nII);
678: ISDestroy(&II);
679: II = nII;
680: } else {
681: IS nII;
682: ISDifference(pII,Pall,&nII);
683: ISDestroy(&pII);
684: pII = nII;
685: }
686: }
687: if (ploc) {
688: ISDifference(lPall,II,&lP);
689: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lP",(PetscObject)lP);
690: } else {
691: ISDifference(Pall,pII,&pP);
692: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",(PetscObject)pP);
693: /* need all local pressure dofs */
694: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
695: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
696: ISGetLocalSize(Pall,&ni);
697: ISGetIndices(Pall,&idxs);
698: for (i=0;i<ni;i++) matis->sf_rootdata[idxs[i]-rst] = 1;
699: ISRestoreIndices(Pall,&idxs);
700: PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
701: PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
702: for (i=0,ni=0;i<n;i++) if (matis->sf_leafdata[i]) widxs[ni++] = i;
703: ISCreateGeneral(PETSC_COMM_SELF,ni,widxs,PETSC_COPY_VALUES,&lPall);
704: }
706: if (!Pall) {
707: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
708: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
709: ISGetLocalSize(lPall,&ni);
710: ISGetIndices(lPall,&idxs);
711: for (i=0;i<ni;i++) matis->sf_leafdata[idxs[i]] = 1;
712: ISRestoreIndices(lPall,&idxs);
713: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
714: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
715: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
716: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&Pall);
717: }
718: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject)Pall);
720: if (flip) {
721: PetscInt npl;
722: ISGetLocalSize(Pall,&npl);
723: ISGetIndices(Pall,&idxs);
724: MatCreateVecs(A,NULL,&fetidp->rhs_flip);
725: VecSet(fetidp->rhs_flip,1.);
726: VecSetOption(fetidp->rhs_flip,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
727: for (i=0;i<npl;i++) {
728: VecSetValue(fetidp->rhs_flip,idxs[i],-1.,INSERT_VALUES);
729: }
730: VecAssemblyBegin(fetidp->rhs_flip);
731: VecAssemblyEnd(fetidp->rhs_flip);
732: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_flip",(PetscObject)fetidp->rhs_flip);
733: ISRestoreIndices(Pall,&idxs);
734: }
735: ISDestroy(&Pall);
736: ISDestroy(&pII);
738: /* local selected pressures in subdomain-wise and global ordering */
739: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
740: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
741: if (!ploc) {
742: PetscInt *widxs2;
744: if (!pP) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Missing parallel pressure IS");
745: ISGetLocalSize(pP,&ni);
746: ISGetIndices(pP,&idxs);
747: for (i=0;i<ni;i++) matis->sf_rootdata[idxs[i]-rst] = 1;
748: ISRestoreIndices(pP,&idxs);
749: PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
750: PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
751: for (i=0,ni=0;i<n;i++) if (matis->sf_leafdata[i]) widxs[ni++] = i;
752: PetscMalloc1(ni,&widxs2);
753: ISLocalToGlobalMappingApply(l2g,ni,widxs,widxs2);
754: ISCreateGeneral(PETSC_COMM_SELF,ni,widxs,PETSC_COPY_VALUES,&lP);
755: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lP",(PetscObject)lP);
756: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs2,PETSC_OWN_POINTER,&is1);
757: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_gP",(PetscObject)is1);
758: ISDestroy(&is1);
759: } else {
760: if (!lP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing sequential pressure IS");
761: ISGetLocalSize(lP,&ni);
762: ISGetIndices(lP,&idxs);
763: for (i=0;i<ni;i++)
764: if (idxs[i] >=0 && idxs[i] < n)
765: matis->sf_leafdata[idxs[i]] = 1;
766: ISRestoreIndices(lP,&idxs);
767: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
768: ISLocalToGlobalMappingApply(l2g,ni,idxs,widxs);
769: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&is1);
770: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_gP",(PetscObject)is1);
771: ISDestroy(&is1);
772: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
773: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
774: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&pP);
775: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",(PetscObject)pP);
776: }
777: PetscFree(widxs);
779: /* If there's any "interior pressure",
780: we may want to use a discrete harmonic solver instead
781: of a Stokes harmonic for the Dirichlet preconditioner
782: Need to extract the interior velocity dofs in interior dofs ordering (iV)
783: and interior pressure dofs in local ordering (iP) */
784: if (!allp) {
785: ISLocalToGlobalMapping l2g_t;
787: ISDifference(lPall,lP,&is1);
788: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_iP",(PetscObject)is1);
789: ISDifference(II,is1,&is2);
790: ISDestroy(&is1);
791: ISLocalToGlobalMappingCreateIS(II,&l2g_t);
792: ISGlobalToLocalMappingApplyIS(l2g_t,IS_GTOLM_DROP,is2,&is1);
793: ISGetLocalSize(is1,&i);
794: ISGetLocalSize(is2,&j);
795: if (i != j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local sizes %D and %D for iV",i,j);
796: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_iV",(PetscObject)is1);
797: ISLocalToGlobalMappingDestroy(&l2g_t);
798: ISDestroy(&is1);
799: ISDestroy(&is2);
800: }
801: ISDestroy(&lPall);
802: ISDestroy(&II);
804: /* exclude selected pressures from the inner BDDC */
805: if (pcbddc->DirichletBoundariesLocal) {
806: IS list[2],plP,isout;
807: PetscInt np;
809: /* need a parallel IS */
810: ISGetLocalSize(lP,&np);
811: ISGetIndices(lP,&idxs);
812: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),np,idxs,PETSC_USE_POINTER,&plP);
813: list[0] = plP;
814: list[1] = pcbddc->DirichletBoundariesLocal;
815: ISConcatenate(PetscObjectComm((PetscObject)ksp),2,list,&isout);
816: ISSortRemoveDups(isout);
817: ISDestroy(&plP);
818: ISRestoreIndices(lP,&idxs);
819: PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc,isout);
820: ISDestroy(&isout);
821: } else if (pcbddc->DirichletBoundaries) {
822: IS list[2],isout;
824: list[0] = pP;
825: list[1] = pcbddc->DirichletBoundaries;
826: ISConcatenate(PetscObjectComm((PetscObject)ksp),2,list,&isout);
827: ISSortRemoveDups(isout);
828: PCBDDCSetDirichletBoundaries(fetidp->innerbddc,isout);
829: ISDestroy(&isout);
830: } else {
831: IS plP;
832: PetscInt np;
834: /* need a parallel IS */
835: ISGetLocalSize(lP,&np);
836: ISGetIndices(lP,&idxs);
837: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),np,idxs,PETSC_COPY_VALUES,&plP);
838: PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc,plP);
839: ISDestroy(&plP);
840: ISRestoreIndices(lP,&idxs);
841: }
842: ISDestroy(&lP);
843: fetidp->pP = pP;
844: }
846: /* total number of selected pressure dofs */
847: ISGetSize(fetidp->pP,&totP);
849: /* Set operator for inner BDDC */
850: if (totP || fetidp->rhs_flip) {
851: MatDuplicate(A,MAT_COPY_VALUES,&nA);
852: } else {
853: PetscObjectReference((PetscObject)A);
854: nA = A;
855: }
856: if (fetidp->rhs_flip) {
857: MatDiagonalScale(nA,fetidp->rhs_flip,NULL);
858: if (totP) {
859: Mat lA2;
861: MatISGetLocalMat(nA,&lA);
862: MatDuplicate(lA,MAT_COPY_VALUES,&lA2);
863: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",(PetscObject)lA2);
864: MatDestroy(&lA2);
865: }
866: }
868: if (totP) {
869: MatSetOption(nA,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
870: MatZeroRowsColumnsIS(nA,fetidp->pP,1.,NULL,NULL);
871: } else {
872: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",NULL);
873: }
874: PCSetOperators(fetidp->innerbddc,nA,nA);
875: MatDestroy(&nA);
877: /* non-zero rhs on interior dofs when applying the preconditioner */
878: if (totP) pcbddc->switch_static = PETSC_TRUE;
880: /* if there are no interface pressures, set inner bddc flag for benign saddle point */
881: if (!totP) {
882: pcbddc->benign_saddle_point = PETSC_TRUE;
883: pcbddc->compute_nonetflux = PETSC_TRUE;
884: }
886: /* Divergence mat */
887: if (totP) {
888: Mat B;
889: IS P;
890: PetscBool save;
892: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&P);
893: MatCreateSubMatrix(A,P,NULL,MAT_INITIAL_MATRIX,&B);
894: save = pcbddc->compute_nonetflux; /* SetDivergenceMat activates nonetflux computation */
895: PCBDDCSetDivergenceMat(fetidp->innerbddc,B,PETSC_FALSE,NULL);
896: pcbddc->compute_nonetflux = save;
897: MatDestroy(&B);
898: }
900: /* Operators for pressure preconditioner */
901: if (totP) {
902: /* Extract pressure block if needed */
903: if (!pisz) {
904: Mat C;
905: IS nzrows = NULL;
907: MatCreateSubMatrix(A,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
908: MatFindNonzeroRows(C,&nzrows);
909: if (nzrows) {
910: PetscInt i;
912: ISGetSize(nzrows,&i);
913: ISDestroy(&nzrows);
914: if (!i) pisz = PETSC_TRUE;
915: }
916: if (!pisz) {
917: MatScale(C,-1.); /* i.e. Almost Incompressible Elasticity, Stokes discretized with Q1xQ1_stabilized */
918: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_C",(PetscObject)C);
919: }
920: MatDestroy(&C);
921: }
922: if (A != Ap) { /* user has provided a different Pmat, use it to extract the pressure preconditioner */
923: Mat C;
925: MatCreateSubMatrix(Ap,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
926: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
927: MatDestroy(&C);
928: }
929: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);
931: /* Preconditioned operator for the pressure block */
932: if (PPmat) {
933: Mat C;
934: IS Pall;
935: PetscInt AM,PAM,PAN,pam,pan,am,an,pl,pIl,pAg,pIg;
936: PetscBool ismatis;
938: PetscObjectTypeCompare((PetscObject)PPmat,MATIS,&ismatis);
939: if (ismatis) {
940: MatISGetMPIXAIJ(PPmat,MAT_INITIAL_MATRIX,&C);
941: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
942: MatDestroy(&C);
943: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);
944: }
945: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&Pall);
946: MatGetSize(A,&AM,NULL);
947: MatGetSize(PPmat,&PAM,&PAN);
948: ISGetSize(Pall,&pAg);
949: ISGetSize(fetidp->pP,&pIg);
950: MatGetLocalSize(PPmat,&pam,&pan);
951: MatGetLocalSize(A,&am,&an);
952: ISGetLocalSize(Pall,&pIl);
953: ISGetLocalSize(fetidp->pP,&pl);
954: if (PAM != PAN) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Pressure matrix must be square, unsupported %D x %D",PAM,PAN);
955: if (pam != pan) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Local sizes of pressure matrix must be equal, unsupported %D x %D",pam,pan);
956: 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);
957: 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);
958: if (PAM == AM) { /* monolithic ordering, restrict to pressure */
959: MatCreateSubMatrix(PPmat,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
960: } else if (pAg == PAM) { /* global ordering for pressure only */
961: if (!allp) { /* solving for interface pressure only */
962: IS restr;
964: ISRenumber(fetidp->pP,NULL,NULL,&restr);
965: MatCreateSubMatrix(PPmat,restr,restr,MAT_INITIAL_MATRIX,&C);
966: ISDestroy(&restr);
967: } else {
968: PetscObjectReference((PetscObject)PPmat);
969: C = PPmat;
970: }
971: } else if (pIg == PAM) { /* global ordering for selected pressure only */
972: PetscObjectReference((PetscObject)PPmat);
973: C = PPmat;
974: } else {
975: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Unable to use the pressure matrix");
976: }
977: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
978: MatDestroy(&C);
979: } else {
980: Mat C;
982: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_C",(PetscObject*)&C);
983: if (C) { /* non-zero pressure block, most likely Almost Incompressible Elasticity */
984: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
985: } else { /* identity (need to be scaled properly by the user using e.g. a Richardson method */
986: PetscInt nl;
988: ISGetLocalSize(fetidp->pP,&nl);
989: MatCreate(PetscObjectComm((PetscObject)ksp),&PPmat);
990: MatSetSizes(PPmat,nl,nl,totP,totP);
991: MatSetType(PPmat,MATAIJ);
992: MatMPIAIJSetPreallocation(PPmat,1,NULL,0,NULL);
993: MatSeqAIJSetPreallocation(PPmat,1,NULL);
994: MatAssemblyBegin(PPmat,MAT_FINAL_ASSEMBLY);
995: MatAssemblyEnd(PPmat,MAT_FINAL_ASSEMBLY);
996: MatShift(PPmat,1.);
997: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)PPmat);
998: MatDestroy(&PPmat);
999: }
1000: }
1001: } else { /* totP == 0 */
1002: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",NULL);
1003: }
1004: }
1005: return(0);
1006: }
1008: static PetscErrorCode KSPSetUp_FETIDP(KSP ksp)
1009: {
1010: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1011: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1012: PetscBool flg;
1016: KSPFETIDPSetUpOperators(ksp);
1017: /* set up BDDC */
1018: PCSetErrorIfFailure(fetidp->innerbddc,ksp->errorifnotconverged);
1019: PCSetUp(fetidp->innerbddc);
1020: /* FETI-DP as it is implemented needs an exact coarse solver */
1021: if (pcbddc->coarse_ksp) {
1022: KSPSetTolerances(pcbddc->coarse_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,1000);
1023: KSPSetNormType(pcbddc->coarse_ksp,KSP_NORM_DEFAULT);
1024: }
1025: /* FETI-DP as it is implemented needs exact local Neumann solvers */
1026: KSPSetTolerances(pcbddc->ksp_R,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,1000);
1027: KSPSetNormType(pcbddc->ksp_R,KSP_NORM_DEFAULT);
1029: /* setup FETI-DP operators
1030: If fetidp->statechanged is true, we need update the operators
1031: that are needed in the saddle-point case. This should be replaced
1032: by a better logic when the FETI-DP matrix and preconditioner will
1033: have their own classes */
1034: if (pcbddc->new_primal_space || fetidp->statechanged) {
1035: Mat F; /* the FETI-DP matrix */
1036: PC D; /* the FETI-DP preconditioner */
1037: KSPReset(fetidp->innerksp);
1038: PCBDDCCreateFETIDPOperators(fetidp->innerbddc,fetidp->fully_redundant,((PetscObject)ksp)->prefix,&F,&D);
1039: KSPSetOperators(fetidp->innerksp,F,F);
1040: KSPSetTolerances(fetidp->innerksp,ksp->rtol,ksp->abstol,ksp->divtol,ksp->max_it);
1041: KSPSetPC(fetidp->innerksp,D);
1042: KSPSetFromOptions(fetidp->innerksp);
1043: MatCreateVecs(F,&(fetidp->innerksp)->vec_rhs,&(fetidp->innerksp)->vec_sol);
1044: MatDestroy(&F);
1045: PCDestroy(&D);
1046: if (fetidp->check) {
1047: PetscViewer viewer;
1049: if (!pcbddc->dbg_viewer) {
1050: viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
1051: } else {
1052: viewer = pcbddc->dbg_viewer;
1053: }
1054: KSPFETIDPCheckOperators(ksp,viewer);
1055: }
1056: }
1057: fetidp->statechanged = PETSC_FALSE;
1058: pcbddc->new_primal_space = PETSC_FALSE;
1060: /* propagate settings to the inner solve */
1061: KSPGetComputeSingularValues(ksp,&flg);
1062: KSPSetComputeSingularValues(fetidp->innerksp,flg);
1063: if (ksp->res_hist) {
1064: KSPSetResidualHistory(fetidp->innerksp,ksp->res_hist,ksp->res_hist_max,ksp->res_hist_reset);
1065: }
1066: KSPSetErrorIfNotConverged(fetidp->innerksp,ksp->errorifnotconverged);
1067: KSPSetUp(fetidp->innerksp);
1068: return(0);
1069: }
1071: static PetscErrorCode KSPSolve_FETIDP(KSP ksp)
1072: {
1074: Mat F,A;
1075: MatNullSpace nsp;
1076: Vec X,B,Xl,Bl;
1077: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1078: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1081: PetscCitationsRegister(citation,&cited);
1082: if (fetidp->saddlepoint) {
1083: PetscCitationsRegister(citation2,&cited2);
1084: }
1085: KSPGetOperators(ksp,&A,NULL);
1086: KSPGetRhs(ksp,&B);
1087: KSPGetSolution(ksp,&X);
1088: KSPGetOperators(fetidp->innerksp,&F,NULL);
1089: KSPGetRhs(fetidp->innerksp,&Bl);
1090: KSPGetSolution(fetidp->innerksp,&Xl);
1091: PCBDDCMatFETIDPGetRHS(F,B,Bl);
1092: if (ksp->transpose_solve) {
1093: KSPSolveTranspose(fetidp->innerksp,Bl,Xl);
1094: } else {
1095: KSPSolve(fetidp->innerksp,Bl,Xl);
1096: }
1097: PCBDDCMatFETIDPGetSolution(F,Xl,X);
1098: MatGetNullSpace(A,&nsp);
1099: if (nsp) {
1100: MatNullSpaceRemove(nsp,X);
1101: }
1102: /* update ksp with stats from inner ksp */
1103: KSPGetConvergedReason(fetidp->innerksp,&ksp->reason);
1104: KSPGetIterationNumber(fetidp->innerksp,&ksp->its);
1105: ksp->totalits += ksp->its;
1106: KSPGetResidualHistory(fetidp->innerksp,NULL,&ksp->res_hist_len);
1107: /* restore defaults for inner BDDC (Pre/PostSolve flags) */
1108: pcbddc->temp_solution_used = PETSC_FALSE;
1109: pcbddc->rhs_change = PETSC_FALSE;
1110: pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1111: return(0);
1112: }
1114: static PetscErrorCode KSPReset_FETIDP(KSP ksp)
1115: {
1116: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1117: PC_BDDC *pcbddc;
1121: ISDestroy(&fetidp->pP);
1122: VecDestroy(&fetidp->rhs_flip);
1123: /* avoid PCReset that does not take into account ref counting */
1124: PCDestroy(&fetidp->innerbddc);
1125: PCCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerbddc);
1126: PCSetType(fetidp->innerbddc,PCBDDC);
1127: pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1128: pcbddc->symmetric_primal = PETSC_FALSE;
1129: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerbddc);
1130: KSPDestroy(&fetidp->innerksp);
1131: fetidp->saddlepoint = PETSC_FALSE;
1132: fetidp->matstate = -1;
1133: fetidp->matnnzstate = -1;
1134: fetidp->statechanged = PETSC_TRUE;
1135: return(0);
1136: }
1138: static PetscErrorCode KSPDestroy_FETIDP(KSP ksp)
1139: {
1140: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1144: KSPReset_FETIDP(ksp);
1145: PCDestroy(&fetidp->innerbddc);
1146: KSPDestroy(&fetidp->innerksp);
1147: PetscFree(fetidp->monctx);
1148: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetInnerBDDC_C",NULL);
1149: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerBDDC_C",NULL);
1150: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerKSP_C",NULL);
1151: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetPressureOperator_C",NULL);
1152: PetscFree(ksp->data);
1153: return(0);
1154: }
1156: static PetscErrorCode KSPView_FETIDP(KSP ksp,PetscViewer viewer)
1157: {
1158: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1160: PetscBool iascii;
1163: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1164: if (iascii) {
1165: PetscViewerASCIIPrintf(viewer," fully redundant: %d\n",fetidp->fully_redundant);
1166: PetscViewerASCIIPrintf(viewer," saddle point: %d\n",fetidp->saddlepoint);
1167: PetscViewerASCIIPrintf(viewer," inner solver details\n");
1168: PetscViewerASCIIAddTab(viewer,2);
1169: }
1170: KSPView(fetidp->innerksp,viewer);
1171: if (iascii) {
1172: PetscViewerASCIISubtractTab(viewer,2);
1173: PetscViewerASCIIPrintf(viewer," BDDC solver details\n");
1174: PetscViewerASCIIAddTab(viewer,2);
1175: }
1176: PCView(fetidp->innerbddc,viewer);
1177: if (iascii) {
1178: PetscViewerASCIISubtractTab(viewer,2);
1179: }
1180: return(0);
1181: }
1183: static PetscErrorCode KSPSetFromOptions_FETIDP(PetscOptionItems *PetscOptionsObject,KSP ksp)
1184: {
1185: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1189: /* set options prefixes for the inner objects, since the parent prefix will be valid at this point */
1190: PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerksp,((PetscObject)ksp)->prefix);
1191: PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerksp,"fetidp_");
1192: if (!fetidp->userbddc) {
1193: PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerbddc,((PetscObject)ksp)->prefix);
1194: PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerbddc,"fetidp_bddc_");
1195: }
1196: PetscOptionsHead(PetscOptionsObject,"KSP FETIDP options");
1197: PetscOptionsBool("-ksp_fetidp_fullyredundant","Use fully redundant multipliers","none",fetidp->fully_redundant,&fetidp->fully_redundant,NULL);
1198: PetscOptionsBool("-ksp_fetidp_saddlepoint","Activates support for saddle-point problems",NULL,fetidp->saddlepoint,&fetidp->saddlepoint,NULL);
1199: PetscOptionsBool("-ksp_fetidp_check","Activates verbose debugging output FETI-DP operators",NULL,fetidp->check,&fetidp->check,NULL);
1200: PetscOptionsTail();
1201: PCSetFromOptions(fetidp->innerbddc);
1202: return(0);
1203: }
1205: /*MC
1206: KSPFETIDP - The FETI-DP method
1208: This class implements the FETI-DP method [1].
1209: The matrix for the KSP must be of type MATIS.
1210: The FETI-DP linear system (automatically generated constructing an internal PCBDDC object) is solved using an internal KSP object.
1212: Options Database Keys:
1213: + -ksp_fetidp_fullyredundant <false> : use a fully redundant set of Lagrange multipliers
1214: . -ksp_fetidp_saddlepoint <false> : activates support for saddle point problems, see [2]
1215: . -ksp_fetidp_saddlepoint_flip <false> : usually, an incompressible Stokes problem is written as
1216: | A B^T | | v | = | f |
1217: | B 0 | | p | = | g |
1218: with B representing -\int_\Omega \nabla \cdot u q.
1219: If -ksp_fetidp_saddlepoint_flip is true, the code assumes that the user provides it as
1220: | A B^T | | v | = | f |
1221: |-B 0 | | p | = |-g |
1222: . -ksp_fetidp_pressure_field <-1> : activates support for saddle point problems, and identifies the pressure field id.
1223: If this information is not provided, the pressure field is detected by using MatFindZeroDiagonals().
1224: - -ksp_fetidp_pressure_all <false> : if false, uses the interface pressures, as described in [2]. If true, uses the entire pressure field.
1226: Level: Advanced
1228: Notes: 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.,
1229: .vb
1230: -fetidp_ksp_type gmres -fetidp_bddc_pc_bddc_symmetric false
1231: .ve
1232: will use GMRES for the solution of the linear system on the Lagrange multipliers, generated using a non-symmetric PCBDDC.
1234: For saddle point problems with continuous pressures, the preconditioned operator for the pressure solver can be specified with KSPFETIDPSetPressureOperator().
1235: Alternatively, the pressure operator is extracted from the precondioned matrix (if it is different from the linear solver matrix).
1236: If none of the above, an identity matrix will be created; the user then needs to scale it through a Richardson solver.
1237: Options for the pressure solver can be prefixed with -fetidp_fielsplit_p_, E.g.
1238: .vb
1239: -fetidp_fielsplit_p_ksp_type preonly -fetidp_fielsplit_p_pc_type lu -fetidp_fielsplit_p_pc_factor_mat_solver_type mumps
1240: .ve
1241: 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
1242: non-redundant multipliers, i.e. -ksp_fetidp_fullyredundant false. Options for the scaling solver are prefixed by -fetidp_bddelta_, E.g.
1243: .vb
1244: -fetidp_bddelta_pc_factor_mat_solver_type mumps -fetidp_bddelta_pc_type lu
1245: .ve
1247: 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.
1249: The converged reason and number of iterations computed are passed from the inner KSP to this KSP at the end of the solution.
1251: Developer Notes: Even though this method does not directly use any norms, the user is allowed to set the KSPNormType to any value.
1252: This is so users do not have to change KSPNormType options when they switch from other KSP methods to this one.
1254: References:
1255: .vb
1256: . [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
1257: . [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
1258: .ve
1260: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC(), KSPFETIDPGetInnerBDDC(), KSPFETIDPGetInnerKSP()
1261: M*/
1262: PETSC_EXTERN PetscErrorCode KSPCreate_FETIDP(KSP ksp)
1263: {
1265: KSP_FETIDP *fetidp;
1266: KSP_FETIDPMon *monctx;
1267: PC_BDDC *pcbddc;
1268: PC pc;
1271: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
1272: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
1273: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
1274: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
1275: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
1276: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
1277: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
1279: PetscNewLog(ksp,&fetidp);
1280: fetidp->matstate = -1;
1281: fetidp->matnnzstate = -1;
1282: fetidp->statechanged = PETSC_TRUE;
1284: ksp->data = (void*)fetidp;
1285: ksp->ops->setup = KSPSetUp_FETIDP;
1286: ksp->ops->solve = KSPSolve_FETIDP;
1287: ksp->ops->destroy = KSPDestroy_FETIDP;
1288: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_FETIDP;
1289: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_FETIDP;
1290: ksp->ops->view = KSPView_FETIDP;
1291: ksp->ops->setfromoptions = KSPSetFromOptions_FETIDP;
1292: ksp->ops->buildsolution = KSPBuildSolution_FETIDP;
1293: ksp->ops->buildresidual = KSPBuildResidualDefault;
1294: KSPGetPC(ksp,&pc);
1295: PCSetType(pc,PCNONE);
1296: /* create the inner KSP for the Lagrange multipliers */
1297: KSPCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerksp);
1298: KSPGetPC(fetidp->innerksp,&pc);
1299: PCSetType(pc,PCNONE);
1300: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerksp);
1301: /* monitor */
1302: PetscNew(&monctx);
1303: monctx->parentksp = ksp;
1304: fetidp->monctx = monctx;
1305: KSPMonitorSet(fetidp->innerksp,KSPMonitor_FETIDP,fetidp->monctx,NULL);
1306: /* create the inner BDDC */
1307: PCCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerbddc);
1308: PCSetType(fetidp->innerbddc,PCBDDC);
1309: /* make sure we always obtain a consistent FETI-DP matrix
1310: for symmetric problems, the user can always customize it through the command line */
1311: pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1312: pcbddc->symmetric_primal = PETSC_FALSE;
1313: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerbddc);
1314: /* composed functions */
1315: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetInnerBDDC_C",KSPFETIDPSetInnerBDDC_FETIDP);
1316: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerBDDC_C",KSPFETIDPGetInnerBDDC_FETIDP);
1317: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerKSP_C",KSPFETIDPGetInnerKSP_FETIDP);
1318: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetPressureOperator_C",KSPFETIDPSetPressureOperator_FETIDP);
1319: /* need to call KSPSetUp_FETIDP even with KSP_SETUP_NEWMATRIX */
1320: ksp->setupnewmatrix = PETSC_TRUE;
1321: return(0);
1322: }