Actual source code: fetidp.c
petsc-3.9.2 2018-05-20
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: PetscBool ismatis,pisz,allp;
524: PetscBool flip; /* Usually, Stokes is written (B = -\int_\Omega \nabla \cdot u q)
525: | A B'| | v | = | f |
526: | B 0 | | p | = | g |
527: If -ksp_fetidp_saddlepoint_flip is true, the code assumes it is written as
528: | A B'| | v | = | f |
529: |-B 0 | | p | = |-g |
530: */
531: PetscObjectState matstate, matnnzstate;
532: PetscErrorCode ierr;
535: pisz = PETSC_FALSE;
536: flip = PETSC_FALSE;
537: allp = PETSC_FALSE;
538: PetscOptionsBegin(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->prefix,"FETI-DP options","PC");
539: PetscOptionsInt("-ksp_fetidp_pressure_field","Field id for pressures for saddle-point problems",NULL,fid,&fid,NULL);
540: PetscOptionsBool("-ksp_fetidp_pressure_all","Use the whole pressure set instead of just that at the interface",NULL,allp,&allp,NULL);
541: PetscOptionsBool("-ksp_fetidp_saddlepoint_flip","Flip the sign of the pressure-velocity (lower-left) block",NULL,flip,&flip,NULL);
542: PetscOptionsEnd();
544: fetidp->saddlepoint = (fid >= 0 ? PETSC_TRUE : fetidp->saddlepoint);
546: KSPGetOperators(ksp,&A,&Ap);
547: PetscObjectTypeCompare((PetscObject)A,MATIS,&ismatis);
548: if (!ismatis) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Amat should be of type MATIS");
550: /* Quiet return if the matrix states are unchanged.
551: Needed only for the saddle point case since it uses MatZeroRows
552: on a matrix that may not have changed */
553: PetscObjectStateGet((PetscObject)A,&matstate);
554: MatGetNonzeroState(A,&matnnzstate);
555: if (matstate == fetidp->matstate && matnnzstate == fetidp->matnnzstate) return(0);
556: fetidp->matstate = matstate;
557: fetidp->matnnzstate = matnnzstate;
558: fetidp->statechanged = fetidp->saddlepoint;
560: /* see if we have some fields attached */
561: if (!pcbddc->n_ISForDofsLocal && !pcbddc->n_ISForDofs) {
562: DM dm;
563: PetscContainer c;
565: KSPGetDM(ksp,&dm);
566: PetscObjectQuery((PetscObject)A,"_convert_nest_lfields",(PetscObject*)&c);
567: if (dm) {
568: IS *fields;
569: PetscInt nf,i;
571: DMCreateFieldDecomposition(dm,&nf,NULL,&fields,NULL);
572: PCBDDCSetDofsSplitting(fetidp->innerbddc,nf,fields);
573: for (i=0;i<nf;i++) {
574: ISDestroy(&fields[i]);
575: }
576: PetscFree(fields);
577: } else if (c) {
578: MatISLocalFields lf;
579: PetscContainerGetPointer(c,(void**)&lf);
580: PCBDDCSetDofsSplittingLocal(fetidp->innerbddc,lf->nr,lf->rf);
581: }
582: }
584: if (!fetidp->saddlepoint) {
585: PCSetOperators(fetidp->innerbddc,A,A);
586: } else {
587: Mat nA,lA;
588: Mat PPmat;
589: IS pP;
590: PetscInt totP;
592: MatISGetLocalMat(A,&lA);
593: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",(PetscObject)lA);
595: pP = fetidp->pP;
596: if (!pP) { /* first time, need to compute pressure dofs */
597: PC_IS *pcis = (PC_IS*)fetidp->innerbddc->data;
598: Mat_IS *matis = (Mat_IS*)(A->data);
599: ISLocalToGlobalMapping l2g;
600: IS lP = NULL,II,pII,lPall,Pall,is1,is2;
601: const PetscInt *idxs;
602: PetscInt nl,ni,*widxs;
603: PetscInt i,j,n_neigh,*neigh,*n_shared,**shared,*count;
604: PetscInt rst,ren,n;
605: PetscBool ploc;
607: MatGetLocalSize(A,&nl,NULL);
608: MatGetOwnershipRange(A,&rst,&ren);
609: MatGetLocalSize(lA,&n,NULL);
610: MatGetLocalToGlobalMapping(A,&l2g,NULL);
612: if (!pcis->is_I_local) { /* need to compute interior dofs */
613: PetscCalloc1(n,&count);
614: ISLocalToGlobalMappingGetInfo(l2g,&n_neigh,&neigh,&n_shared,&shared);
615: for (i=1;i<n_neigh;i++)
616: for (j=0;j<n_shared[i];j++)
617: count[shared[i][j]] += 1;
618: for (i=0,j=0;i<n;i++) if (!count[i]) count[j++] = i;
619: ISLocalToGlobalMappingRestoreInfo(l2g,&n_neigh,&neigh,&n_shared,&shared);
620: ISCreateGeneral(PETSC_COMM_SELF,j,count,PETSC_OWN_POINTER,&II);
621: } else {
622: PetscObjectReference((PetscObject)pcis->is_I_local);
623: II = pcis->is_I_local;
624: }
626: /* interior dofs in layout */
627: MatISSetUpSF(A);
628: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
629: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
630: ISGetLocalSize(II,&ni);
631: ISGetIndices(II,&idxs);
632: for (i=0;i<ni;i++) matis->sf_leafdata[idxs[i]] = 1;
633: ISRestoreIndices(II,&idxs);
634: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
635: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
636: PetscMalloc1(PetscMax(nl,n),&widxs);
637: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
638: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&pII);
640: /* pressure dofs */
641: Pall = NULL;
642: lPall = NULL;
643: ploc = PETSC_FALSE;
644: if (fid >= 0) {
645: if (pcbddc->n_ISForDofsLocal) {
646: PetscInt np;
648: if (fid >= pcbddc->n_ISForDofsLocal) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Invalid field id for pressure %D, max %D",fid,pcbddc->n_ISForDofsLocal);
649: /* need a sequential IS */
650: ISGetLocalSize(pcbddc->ISForDofsLocal[fid],&np);
651: ISGetIndices(pcbddc->ISForDofsLocal[fid],&idxs);
652: ISCreateGeneral(PETSC_COMM_SELF,np,idxs,PETSC_COPY_VALUES,&lPall);
653: ISRestoreIndices(pcbddc->ISForDofsLocal[fid],&idxs);
654: ploc = PETSC_TRUE;
655: } else if (pcbddc->n_ISForDofs) {
656: if (fid >= pcbddc->n_ISForDofs) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Invalid field id for pressure %D, max %D",fid,pcbddc->n_ISForDofs);
657: PetscObjectReference((PetscObject)pcbddc->ISForDofs[fid]);
658: Pall = pcbddc->ISForDofs[fid];
659: } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Missing fields! Use PCBDDCSetDofsSplitting/Local");
660: } else { /* fallback to zero pressure block */
661: IS list[2];
663: MatFindZeroDiagonals(A,&list[1]);
664: ISComplement(list[1],rst,ren,&list[0]);
665: PCBDDCSetDofsSplitting(fetidp->innerbddc,2,list);
666: ISDestroy(&list[0]);
667: Pall = list[1];
668: }
669: /* if the user requested the entire pressure,
670: remove the interior pressure dofs from II (or pII) */
671: if (allp) {
672: if (ploc) {
673: IS nII;
674: ISDifference(II,lPall,&nII);
675: ISDestroy(&II);
676: II = nII;
677: } else {
678: IS nII;
679: ISDifference(pII,Pall,&nII);
680: ISDestroy(&pII);
681: pII = nII;
682: }
683: }
684: if (ploc) {
685: ISDifference(lPall,II,&lP);
686: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lP",(PetscObject)lP);
687: } else {
688: ISDifference(Pall,pII,&pP);
689: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",(PetscObject)pP);
690: /* need all local pressure dofs */
691: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
692: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
693: ISGetLocalSize(Pall,&ni);
694: ISGetIndices(Pall,&idxs);
695: for (i=0;i<ni;i++) matis->sf_rootdata[idxs[i]-rst] = 1;
696: ISRestoreIndices(Pall,&idxs);
697: PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
698: PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
699: for (i=0,ni=0;i<n;i++) if (matis->sf_leafdata[i]) widxs[ni++] = i;
700: ISCreateGeneral(PETSC_COMM_SELF,ni,widxs,PETSC_COPY_VALUES,&lPall);
701: }
703: if (!Pall) {
704: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
705: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
706: ISGetLocalSize(lPall,&ni);
707: ISGetIndices(lPall,&idxs);
708: for (i=0;i<ni;i++) matis->sf_leafdata[idxs[i]] = 1;
709: ISRestoreIndices(lPall,&idxs);
710: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
711: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
712: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
713: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&Pall);
714: }
715: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject)Pall);
717: if (flip) {
718: PetscInt npl;
719: ISGetLocalSize(Pall,&npl);
720: ISGetIndices(Pall,&idxs);
721: MatCreateVecs(A,NULL,&fetidp->rhs_flip);
722: VecSet(fetidp->rhs_flip,1.);
723: VecSetOption(fetidp->rhs_flip,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);
724: for (i=0;i<npl;i++) {
725: VecSetValue(fetidp->rhs_flip,idxs[i],-1.,INSERT_VALUES);
726: }
727: VecAssemblyBegin(fetidp->rhs_flip);
728: VecAssemblyEnd(fetidp->rhs_flip);
729: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_flip",(PetscObject)fetidp->rhs_flip);
730: ISRestoreIndices(Pall,&idxs);
731: }
732: ISDestroy(&Pall);
733: ISDestroy(&pII);
735: /* local selected pressures in subdomain-wise and global ordering */
736: PetscMemzero(matis->sf_leafdata,n*sizeof(PetscInt));
737: PetscMemzero(matis->sf_rootdata,nl*sizeof(PetscInt));
738: if (!ploc) {
739: PetscInt *widxs2;
741: if (!pP) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Missing parallel pressure IS");
742: ISGetLocalSize(pP,&ni);
743: ISGetIndices(pP,&idxs);
744: for (i=0;i<ni;i++) matis->sf_rootdata[idxs[i]-rst] = 1;
745: ISRestoreIndices(pP,&idxs);
746: PetscSFBcastBegin(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
747: PetscSFBcastEnd(matis->sf,MPIU_INT,matis->sf_rootdata,matis->sf_leafdata);
748: for (i=0,ni=0;i<n;i++) if (matis->sf_leafdata[i]) widxs[ni++] = i;
749: PetscMalloc1(ni,&widxs2);
750: ISLocalToGlobalMappingApply(l2g,ni,widxs,widxs2);
751: ISCreateGeneral(PETSC_COMM_SELF,ni,widxs,PETSC_COPY_VALUES,&lP);
752: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lP",(PetscObject)lP);
753: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs2,PETSC_OWN_POINTER,&is1);
754: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_gP",(PetscObject)is1);
755: ISDestroy(&is1);
756: } else {
757: if (!lP) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing sequential pressure IS");
758: ISGetLocalSize(lP,&ni);
759: ISGetIndices(lP,&idxs);
760: for (i=0;i<ni;i++)
761: if (idxs[i] >=0 && idxs[i] < n)
762: matis->sf_leafdata[idxs[i]] = 1;
763: ISRestoreIndices(lP,&idxs);
764: PetscSFReduceBegin(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
765: ISLocalToGlobalMappingApply(l2g,ni,idxs,widxs);
766: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&is1);
767: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_gP",(PetscObject)is1);
768: ISDestroy(&is1);
769: PetscSFReduceEnd(matis->sf,MPIU_INT,matis->sf_leafdata,matis->sf_rootdata,MPIU_REPLACE);
770: for (i=0,ni=0;i<nl;i++) if (matis->sf_rootdata[i]) widxs[ni++] = i+rst;
771: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),ni,widxs,PETSC_COPY_VALUES,&pP);
772: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",(PetscObject)pP);
773: }
774: PetscFree(widxs);
776: /* If there's any "interior pressure",
777: we may want to use a discrete harmonic solver instead
778: of a Stokes harmonic for the Dirichlet preconditioner
779: Need to extract the interior velocity dofs in interior dofs ordering (iV)
780: and interior pressure dofs in local ordering (iP) */
781: if (!allp) {
782: ISLocalToGlobalMapping l2g_t;
784: ISDifference(lPall,lP,&is1);
785: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_iP",(PetscObject)is1);
786: ISDifference(II,is1,&is2);
787: ISDestroy(&is1);
788: ISLocalToGlobalMappingCreateIS(II,&l2g_t);
789: ISGlobalToLocalMappingApplyIS(l2g_t,IS_GTOLM_DROP,is2,&is1);
790: ISGetLocalSize(is1,&i);
791: ISGetLocalSize(is2,&j);
792: if (i != j) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local sizes %D and %D for iV",i,j);
793: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_iV",(PetscObject)is1);
794: ISLocalToGlobalMappingDestroy(&l2g_t);
795: ISDestroy(&is1);
796: ISDestroy(&is2);
797: }
798: ISDestroy(&lPall);
799: ISDestroy(&II);
801: /* exclude selected pressures from the inner BDDC */
802: if (pcbddc->DirichletBoundariesLocal) {
803: IS list[2],plP,isout;
804: PetscInt np;
806: /* need a parallel IS */
807: ISGetLocalSize(lP,&np);
808: ISGetIndices(lP,&idxs);
809: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),np,idxs,PETSC_USE_POINTER,&plP);
810: list[0] = plP;
811: list[1] = pcbddc->DirichletBoundariesLocal;
812: ISConcatenate(PetscObjectComm((PetscObject)ksp),2,list,&isout);
813: ISSortRemoveDups(isout);
814: ISDestroy(&plP);
815: ISRestoreIndices(lP,&idxs);
816: PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc,isout);
817: ISDestroy(&isout);
818: } else if (pcbddc->DirichletBoundaries) {
819: IS list[2],isout;
821: list[0] = pP;
822: list[1] = pcbddc->DirichletBoundaries;
823: ISConcatenate(PetscObjectComm((PetscObject)ksp),2,list,&isout);
824: ISSortRemoveDups(isout);
825: PCBDDCSetDirichletBoundaries(fetidp->innerbddc,isout);
826: ISDestroy(&isout);
827: } else {
828: IS plP;
829: PetscInt np;
831: /* need a parallel IS */
832: ISGetLocalSize(lP,&np);
833: ISGetIndices(lP,&idxs);
834: ISCreateGeneral(PetscObjectComm((PetscObject)ksp),np,idxs,PETSC_COPY_VALUES,&plP);
835: PCBDDCSetDirichletBoundariesLocal(fetidp->innerbddc,plP);
836: ISDestroy(&plP);
837: ISRestoreIndices(lP,&idxs);
838: }
839: ISDestroy(&lP);
840: fetidp->pP = pP;
841: }
843: /* total number of selected pressure dofs */
844: ISGetSize(fetidp->pP,&totP);
846: /* Set operator for inner BDDC */
847: if (totP || fetidp->rhs_flip) {
848: MatDuplicate(A,MAT_COPY_VALUES,&nA);
849: } else {
850: PetscObjectReference((PetscObject)A);
851: nA = A;
852: }
853: if (fetidp->rhs_flip) {
854: MatDiagonalScale(nA,fetidp->rhs_flip,NULL);
855: if (totP) {
856: Mat lA2;
858: MatISGetLocalMat(nA,&lA);
859: MatDuplicate(lA,MAT_COPY_VALUES,&lA2);
860: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",(PetscObject)lA2);
861: MatDestroy(&lA2);
862: }
863: }
865: if (totP) {
866: MatSetOption(nA,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE);
867: MatZeroRowsColumnsIS(nA,fetidp->pP,1.,NULL,NULL);
868: } else {
869: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_lA",NULL);
870: }
871: PCSetOperators(fetidp->innerbddc,nA,nA);
872: MatDestroy(&nA);
874: /* non-zero rhs on interior dofs when applying the preconditioner */
875: if (totP) pcbddc->switch_static = PETSC_TRUE;
877: /* if there are no interface pressures, set inner bddc flag for benign saddle point */
878: if (!totP) {
879: pcbddc->benign_saddle_point = PETSC_TRUE;
880: pcbddc->compute_nonetflux = PETSC_TRUE;
881: }
883: /* Divergence mat */
884: if (totP) {
885: Mat B;
886: IS P;
887: PetscBool save;
889: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&P);
890: MatCreateSubMatrix(A,P,NULL,MAT_INITIAL_MATRIX,&B);
891: save = pcbddc->compute_nonetflux; /* SetDivergenceMat activates nonetflux computation */
892: PCBDDCSetDivergenceMat(fetidp->innerbddc,B,PETSC_FALSE,NULL);
893: pcbddc->compute_nonetflux = save;
894: MatDestroy(&B);
895: }
897: /* Operators for pressure preconditioner */
898: if (totP) {
899: /* Extract pressure block if needed */
900: if (!pisz) {
901: Mat C;
902: IS nzrows = NULL;
904: MatCreateSubMatrix(A,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
905: MatFindNonzeroRows(C,&nzrows);
906: if (nzrows) {
907: PetscInt i;
909: ISGetSize(nzrows,&i);
910: ISDestroy(&nzrows);
911: if (!i) pisz = PETSC_TRUE;
912: }
913: if (!pisz) {
914: MatScale(C,-1.); /* i.e. Almost Incompressible Elasticity, Stokes discretized with Q1xQ1_stabilized */
915: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_C",(PetscObject)C);
916: }
917: MatDestroy(&C);
918: }
919: if (A != Ap) { /* user has provided a different Pmat, use it to extract the pressure preconditioner */
920: Mat C;
922: MatCreateSubMatrix(Ap,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
923: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
924: MatDestroy(&C);
925: }
926: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);
928: /* Preconditioned operator for the pressure block */
929: if (PPmat) {
930: Mat C;
931: IS Pall;
932: PetscInt AM,PAM,PAN,pam,pan,am,an,pl,pIl,pAg,pIg;
933: PetscBool ismatis;
935: PetscObjectTypeCompare((PetscObject)PPmat,MATIS,&ismatis);
936: if (ismatis) {
937: MatISGetMPIXAIJ(PPmat,MAT_INITIAL_MATRIX,&C);
938: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
939: MatDestroy(&C);
940: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject*)&PPmat);
941: }
942: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_aP",(PetscObject*)&Pall);
943: MatGetSize(A,&AM,NULL);
944: MatGetSize(PPmat,&PAM,&PAN);
945: ISGetSize(Pall,&pAg);
946: ISGetSize(fetidp->pP,&pIg);
947: MatGetLocalSize(PPmat,&pam,&pan);
948: MatGetLocalSize(A,&am,&an);
949: ISGetLocalSize(Pall,&pIl);
950: ISGetLocalSize(fetidp->pP,&pl);
951: if (PAM != PAN) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Pressure matrix must be square, unsupported %D x %D",PAM,PAN);
952: if (pam != pan) SETERRQ2(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Local sizes of pressure matrix must be equal, unsupported %D x %D",pam,pan);
953: 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);
954: 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);
955: if (PAM == AM) { /* monolithic ordering, restrict to pressure */
956: MatCreateSubMatrix(PPmat,fetidp->pP,fetidp->pP,MAT_INITIAL_MATRIX,&C);
957: } else if (pAg == PAM) { /* global ordering for pressure only */
958: if (!allp) { /* solving for interface pressure only */
959: IS restr;
961: ISRenumber(fetidp->pP,NULL,NULL,&restr);
962: MatCreateSubMatrix(PPmat,restr,restr,MAT_INITIAL_MATRIX,&C);
963: ISDestroy(&restr);
964: } else {
965: PetscObjectReference((PetscObject)PPmat);
966: C = PPmat;
967: }
968: } else if (pIg == PAM) { /* global ordering for selected pressure only */
969: PetscObjectReference((PetscObject)PPmat);
970: C = PPmat;
971: } else {
972: SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_USER,"Unable to use the pressure matrix");
973: }
974: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
975: MatDestroy(&C);
976: } else {
977: Mat C;
979: PetscObjectQuery((PetscObject)fetidp->innerbddc,"__KSPFETIDP_C",(PetscObject*)&C);
980: if (C) { /* non-zero pressure block, most likely Almost Incompressible Elasticity */
981: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)C);
982: } else { /* identity (need to be scaled properly by the user using e.g. a Richardson method */
983: PetscInt nl;
985: ISGetLocalSize(fetidp->pP,&nl);
986: MatCreate(PetscObjectComm((PetscObject)ksp),&PPmat);
987: MatSetSizes(PPmat,nl,nl,totP,totP);
988: MatSetType(PPmat,MATAIJ);
989: MatMPIAIJSetPreallocation(PPmat,1,NULL,0,NULL);
990: MatSeqAIJSetPreallocation(PPmat,1,NULL);
991: MatAssemblyBegin(PPmat,MAT_FINAL_ASSEMBLY);
992: MatAssemblyEnd(PPmat,MAT_FINAL_ASSEMBLY);
993: MatShift(PPmat,1.);
994: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_PPmat",(PetscObject)PPmat);
995: MatDestroy(&PPmat);
996: }
997: }
998: } else { /* totP == 0 */
999: PetscObjectCompose((PetscObject)fetidp->innerbddc,"__KSPFETIDP_pP",NULL);
1000: }
1001: }
1002: return(0);
1003: }
1005: static PetscErrorCode KSPSetUp_FETIDP(KSP ksp)
1006: {
1007: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1008: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1009: PetscBool flg;
1013: KSPFETIDPSetUpOperators(ksp);
1014: /* set up BDDC */
1015: PCSetErrorIfFailure(fetidp->innerbddc,ksp->errorifnotconverged);
1016: PCSetUp(fetidp->innerbddc);
1017: /* FETI-DP as it is implemented needs an exact coarse solver */
1018: if (pcbddc->coarse_ksp) {
1019: KSPSetTolerances(pcbddc->coarse_ksp,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,1000);
1020: KSPSetNormType(pcbddc->coarse_ksp,KSP_NORM_DEFAULT);
1021: }
1022: /* FETI-DP as it is implemented needs exact local Neumann solvers */
1023: KSPSetTolerances(pcbddc->ksp_R,PETSC_SMALL,PETSC_SMALL,PETSC_DEFAULT,1000);
1024: KSPSetNormType(pcbddc->ksp_R,KSP_NORM_DEFAULT);
1026: /* setup FETI-DP operators
1027: If fetidp->statechanged is true, we need update the operators
1028: that are needed in the saddle-point case. This should be replaced
1029: by a better logic when the FETI-DP matrix and preconditioner will
1030: have their own classes */
1031: if (pcbddc->new_primal_space || fetidp->statechanged) {
1032: Mat F; /* the FETI-DP matrix */
1033: PC D; /* the FETI-DP preconditioner */
1034: KSPReset(fetidp->innerksp);
1035: PCBDDCCreateFETIDPOperators(fetidp->innerbddc,fetidp->fully_redundant,((PetscObject)ksp)->prefix,&F,&D);
1036: KSPSetOperators(fetidp->innerksp,F,F);
1037: KSPSetTolerances(fetidp->innerksp,ksp->rtol,ksp->abstol,ksp->divtol,ksp->max_it);
1038: KSPSetPC(fetidp->innerksp,D);
1039: KSPSetFromOptions(fetidp->innerksp);
1040: MatCreateVecs(F,&(fetidp->innerksp)->vec_rhs,&(fetidp->innerksp)->vec_sol);
1041: MatDestroy(&F);
1042: PCDestroy(&D);
1043: if (fetidp->check) {
1044: PetscViewer viewer;
1046: if (!pcbddc->dbg_viewer) {
1047: viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ksp));
1048: } else {
1049: viewer = pcbddc->dbg_viewer;
1050: }
1051: KSPFETIDPCheckOperators(ksp,viewer);
1052: }
1053: }
1054: fetidp->statechanged = PETSC_FALSE;
1055: pcbddc->new_primal_space = PETSC_FALSE;
1057: /* propagate settings to the inner solve */
1058: KSPGetComputeSingularValues(ksp,&flg);
1059: KSPSetComputeSingularValues(fetidp->innerksp,flg);
1060: if (ksp->res_hist) {
1061: KSPSetResidualHistory(fetidp->innerksp,ksp->res_hist,ksp->res_hist_max,ksp->res_hist_reset);
1062: }
1063: KSPSetErrorIfNotConverged(fetidp->innerksp,ksp->errorifnotconverged);
1064: KSPSetUp(fetidp->innerksp);
1065: return(0);
1066: }
1068: static PetscErrorCode KSPSolve_FETIDP(KSP ksp)
1069: {
1071: Mat F,A;
1072: MatNullSpace nsp;
1073: Vec X,B,Xl,Bl;
1074: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1075: PC_BDDC *pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1078: PetscCitationsRegister(citation,&cited);
1079: if (fetidp->saddlepoint) {
1080: PetscCitationsRegister(citation2,&cited2);
1081: }
1082: KSPGetOperators(ksp,&A,NULL);
1083: KSPGetRhs(ksp,&B);
1084: KSPGetSolution(ksp,&X);
1085: KSPGetOperators(fetidp->innerksp,&F,NULL);
1086: KSPGetRhs(fetidp->innerksp,&Bl);
1087: KSPGetSolution(fetidp->innerksp,&Xl);
1088: PCBDDCMatFETIDPGetRHS(F,B,Bl);
1089: if (ksp->transpose_solve) {
1090: KSPSolveTranspose(fetidp->innerksp,Bl,Xl);
1091: } else {
1092: KSPSolve(fetidp->innerksp,Bl,Xl);
1093: }
1094: PCBDDCMatFETIDPGetSolution(F,Xl,X);
1095: MatGetNullSpace(A,&nsp);
1096: if (nsp) {
1097: MatNullSpaceRemove(nsp,X);
1098: }
1099: /* update ksp with stats from inner ksp */
1100: KSPGetConvergedReason(fetidp->innerksp,&ksp->reason);
1101: KSPGetIterationNumber(fetidp->innerksp,&ksp->its);
1102: ksp->totalits += ksp->its;
1103: KSPGetResidualHistory(fetidp->innerksp,NULL,&ksp->res_hist_len);
1104: /* restore defaults for inner BDDC (Pre/PostSolve flags) */
1105: pcbddc->temp_solution_used = PETSC_FALSE;
1106: pcbddc->rhs_change = PETSC_FALSE;
1107: pcbddc->exact_dirichlet_trick_app = PETSC_FALSE;
1108: return(0);
1109: }
1111: static PetscErrorCode KSPReset_FETIDP(KSP ksp)
1112: {
1113: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1114: PC_BDDC *pcbddc;
1118: ISDestroy(&fetidp->pP);
1119: VecDestroy(&fetidp->rhs_flip);
1120: /* avoid PCReset that does not take into account ref counting */
1121: PCDestroy(&fetidp->innerbddc);
1122: PCCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerbddc);
1123: PCSetType(fetidp->innerbddc,PCBDDC);
1124: pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1125: pcbddc->symmetric_primal = PETSC_FALSE;
1126: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerbddc);
1127: KSPDestroy(&fetidp->innerksp);
1128: fetidp->saddlepoint = PETSC_FALSE;
1129: fetidp->matstate = -1;
1130: fetidp->matnnzstate = -1;
1131: fetidp->statechanged = PETSC_TRUE;
1132: return(0);
1133: }
1135: static PetscErrorCode KSPDestroy_FETIDP(KSP ksp)
1136: {
1137: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1141: KSPReset_FETIDP(ksp);
1142: PCDestroy(&fetidp->innerbddc);
1143: KSPDestroy(&fetidp->innerksp);
1144: PetscFree(fetidp->monctx);
1145: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetInnerBDDC_C",NULL);
1146: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerBDDC_C",NULL);
1147: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerKSP_C",NULL);
1148: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetPressureOperator_C",NULL);
1149: PetscFree(ksp->data);
1150: return(0);
1151: }
1153: static PetscErrorCode KSPView_FETIDP(KSP ksp,PetscViewer viewer)
1154: {
1155: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1157: PetscBool iascii;
1160: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
1161: if (iascii) {
1162: PetscViewerASCIIPrintf(viewer," fully redundant: %d\n",fetidp->fully_redundant);
1163: PetscViewerASCIIPrintf(viewer," saddle point: %d\n",fetidp->saddlepoint);
1164: PetscViewerASCIIPrintf(viewer," inner solver details\n");
1165: PetscViewerASCIIAddTab(viewer,2);
1166: }
1167: KSPView(fetidp->innerksp,viewer);
1168: if (iascii) {
1169: PetscViewerASCIISubtractTab(viewer,2);
1170: PetscViewerASCIIPrintf(viewer," BDDC solver details\n");
1171: PetscViewerASCIIAddTab(viewer,2);
1172: }
1173: PCView(fetidp->innerbddc,viewer);
1174: if (iascii) {
1175: PetscViewerASCIISubtractTab(viewer,2);
1176: }
1177: return(0);
1178: }
1180: static PetscErrorCode KSPSetFromOptions_FETIDP(PetscOptionItems *PetscOptionsObject,KSP ksp)
1181: {
1182: KSP_FETIDP *fetidp = (KSP_FETIDP*)ksp->data;
1186: /* set options prefixes for the inner objects, since the parent prefix will be valid at this point */
1187: PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerksp,((PetscObject)ksp)->prefix);
1188: PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerksp,"fetidp_");
1189: if (!fetidp->userbddc) {
1190: PetscObjectSetOptionsPrefix((PetscObject)fetidp->innerbddc,((PetscObject)ksp)->prefix);
1191: PetscObjectAppendOptionsPrefix((PetscObject)fetidp->innerbddc,"fetidp_bddc_");
1192: }
1193: PetscOptionsHead(PetscOptionsObject,"KSP FETIDP options");
1194: PetscOptionsBool("-ksp_fetidp_fullyredundant","Use fully redundant multipliers","none",fetidp->fully_redundant,&fetidp->fully_redundant,NULL);
1195: PetscOptionsBool("-ksp_fetidp_saddlepoint","Activates support for saddle-point problems",NULL,fetidp->saddlepoint,&fetidp->saddlepoint,NULL);
1196: PetscOptionsBool("-ksp_fetidp_check","Activates verbose debugging output FETI-DP operators",NULL,fetidp->check,&fetidp->check,NULL);
1197: PetscOptionsTail();
1198: PCSetFromOptions(fetidp->innerbddc);
1199: return(0);
1200: }
1202: /*MC
1203: KSPFETIDP - The FETI-DP method
1205: This class implements the FETI-DP method [1].
1206: The matrix for the KSP must be of type MATIS.
1207: The FETI-DP linear system (automatically generated constructing an internal PCBDDC object) is solved using an internal KSP object.
1209: Options Database Keys:
1210: + -ksp_fetidp_fullyredundant <false> : use a fully redundant set of Lagrange multipliers
1211: . -ksp_fetidp_saddlepoint <false> : activates support for saddle point problems, see [2]
1212: . -ksp_fetidp_saddlepoint_flip <false> : usually, an incompressible Stokes problem is written as
1213: | A B^T | | v | = | f |
1214: | B 0 | | p | = | g |
1215: with B representing -\int_\Omega \nabla \cdot u q.
1216: If -ksp_fetidp_saddlepoint_flip is true, the code assumes that the user provides it as
1217: | A B^T | | v | = | f |
1218: |-B 0 | | p | = |-g |
1219: . -ksp_fetidp_pressure_field <-1> : activates support for saddle point problems, and identifies the pressure field id.
1220: If this information is not provided, the pressure field is detected by using MatFindZeroDiagonals().
1221: - -ksp_fetidp_pressure_all <false> : if false, uses the interface pressures, as described in [2]. If true, uses the entire pressure field.
1223: Level: Advanced
1225: 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.,
1226: .vb
1227: -fetidp_ksp_type gmres -fetidp_bddc_pc_bddc_symmetric false
1228: .ve
1229: will use GMRES for the solution of the linear system on the Lagrange multipliers, generated using a non-symmetric PCBDDC.
1231: For saddle point problems with continuous pressures, the preconditioned operator for the pressure solver can be specified with KSPFETIDPSetPressureOperator().
1232: Alternatively, the pressure operator is extracted from the precondioned matrix (if it is different from the linear solver matrix).
1233: If none of the above, an identity matrix will be created; the user then needs to scale it through a Richardson solver.
1234: Options for the pressure solver can be prefixed with -fetidp_fielsplit_p_, E.g.
1235: .vb
1236: -fetidp_fielsplit_p_ksp_type preonly -fetidp_fielsplit_p_pc_type lu -fetidp_fielsplit_p_pc_factor_mat_solver_type mumps
1237: .ve
1238: 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
1239: non-redundant multipliers, i.e. -ksp_fetidp_fullyredundant false. Options for the scaling solver are prefixed by -fetidp_bddelta_, E.g.
1240: .vb
1241: -fetidp_bddelta_pc_factor_mat_solver_type mumps -fetidp_bddelta_pc_type lu
1242: .ve
1244: 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.
1246: The converged reason and number of iterations computed are passed from the inner KSP to this KSP at the end of the solution.
1248: Developer Notes: Even though this method does not directly use any norms, the user is allowed to set the KSPNormType to any value.
1249: This is so users do not have to change KSPNormType options when they switch from other KSP methods to this one.
1251: References:
1252: .vb
1253: . [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
1254: . [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
1255: .ve
1257: .seealso: MATIS, PCBDDC, KSPFETIDPSetInnerBDDC(), KSPFETIDPGetInnerBDDC(), KSPFETIDPGetInnerKSP()
1258: M*/
1259: PETSC_EXTERN PetscErrorCode KSPCreate_FETIDP(KSP ksp)
1260: {
1262: KSP_FETIDP *fetidp;
1263: KSP_FETIDPMon *monctx;
1264: PC_BDDC *pcbddc;
1265: PC pc;
1268: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_LEFT,3);
1269: KSPSetSupportedNorm(ksp,KSP_NORM_NONE,PC_RIGHT,2);
1270: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_LEFT,2);
1271: KSPSetSupportedNorm(ksp,KSP_NORM_PRECONDITIONED,PC_RIGHT,2);
1272: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_LEFT,2);
1273: KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,2);
1274: KSPSetSupportedNorm(ksp,KSP_NORM_NATURAL,PC_LEFT,2);
1276: PetscNewLog(ksp,&fetidp);
1277: fetidp->matstate = -1;
1278: fetidp->matnnzstate = -1;
1279: fetidp->statechanged = PETSC_TRUE;
1281: ksp->data = (void*)fetidp;
1282: ksp->ops->setup = KSPSetUp_FETIDP;
1283: ksp->ops->solve = KSPSolve_FETIDP;
1284: ksp->ops->destroy = KSPDestroy_FETIDP;
1285: ksp->ops->computeeigenvalues = KSPComputeEigenvalues_FETIDP;
1286: ksp->ops->computeextremesingularvalues = KSPComputeExtremeSingularValues_FETIDP;
1287: ksp->ops->view = KSPView_FETIDP;
1288: ksp->ops->setfromoptions = KSPSetFromOptions_FETIDP;
1289: ksp->ops->buildsolution = KSPBuildSolution_FETIDP;
1290: ksp->ops->buildresidual = KSPBuildResidualDefault;
1291: KSPGetPC(ksp,&pc);
1292: PCSetType(pc,PCNONE);
1293: /* create the inner KSP for the Lagrange multipliers */
1294: KSPCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerksp);
1295: KSPGetPC(fetidp->innerksp,&pc);
1296: PCSetType(pc,PCNONE);
1297: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerksp);
1298: /* monitor */
1299: PetscNew(&monctx);
1300: monctx->parentksp = ksp;
1301: fetidp->monctx = monctx;
1302: KSPMonitorSet(fetidp->innerksp,KSPMonitor_FETIDP,fetidp->monctx,NULL);
1303: /* create the inner BDDC */
1304: PCCreate(PetscObjectComm((PetscObject)ksp),&fetidp->innerbddc);
1305: PCSetType(fetidp->innerbddc,PCBDDC);
1306: /* make sure we always obtain a consistent FETI-DP matrix
1307: for symmetric problems, the user can always customize it through the command line */
1308: pcbddc = (PC_BDDC*)fetidp->innerbddc->data;
1309: pcbddc->symmetric_primal = PETSC_FALSE;
1310: PetscLogObjectParent((PetscObject)ksp,(PetscObject)fetidp->innerbddc);
1311: /* composed functions */
1312: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetInnerBDDC_C",KSPFETIDPSetInnerBDDC_FETIDP);
1313: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerBDDC_C",KSPFETIDPGetInnerBDDC_FETIDP);
1314: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPGetInnerKSP_C",KSPFETIDPGetInnerKSP_FETIDP);
1315: PetscObjectComposeFunction((PetscObject)ksp,"KSPFETIDPSetPressureOperator_C",KSPFETIDPSetPressureOperator_FETIDP);
1316: /* need to call KSPSetUp_FETIDP even with KSP_SETUP_NEWMATRIX */
1317: ksp->setupnewmatrix = PETSC_TRUE;
1318: return(0);
1319: }