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