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