1 2 /* 3 Defines a block Jacobi preconditioner. 4 */ 5 6 #include <../src/ksp/pc/impls/bjacobi/bjacobi.h> /*I "petscpc.h" I*/ 7 8 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC, Mat, Mat); 9 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC, Mat, Mat); 10 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC); 11 12 static PetscErrorCode PCSetUp_BJacobi(PC pc) 13 { 14 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 15 Mat mat = pc->mat, pmat = pc->pmat; 16 PetscBool hasop; 17 PetscInt N, M, start, i, sum, end; 18 PetscInt bs, i_start = -1, i_end = -1; 19 PetscMPIInt rank, size; 20 21 PetscFunctionBegin; 22 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 23 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 24 PetscCall(MatGetLocalSize(pc->pmat, &M, &N)); 25 PetscCall(MatGetBlockSize(pc->pmat, &bs)); 26 27 if (jac->n > 0 && jac->n < size) { 28 PetscCall(PCSetUp_BJacobi_Multiproc(pc)); 29 PetscFunctionReturn(0); 30 } 31 32 /* Determines the number of blocks assigned to each processor */ 33 /* local block count given */ 34 if (jac->n_local > 0 && jac->n < 0) { 35 PetscCall(MPIU_Allreduce(&jac->n_local, &jac->n, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc))); 36 if (jac->l_lens) { /* check that user set these correctly */ 37 sum = 0; 38 for (i = 0; i < jac->n_local; i++) { 39 PetscCheck(jac->l_lens[i] / bs * bs == jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout"); 40 sum += jac->l_lens[i]; 41 } 42 PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Local lens set incorrectly"); 43 } else { 44 PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 45 for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)); 46 } 47 } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */ 48 /* global blocks given: determine which ones are local */ 49 if (jac->g_lens) { 50 /* check if the g_lens is has valid entries */ 51 for (i = 0; i < jac->n; i++) { 52 PetscCheck(jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Zero block not allowed"); 53 PetscCheck(jac->g_lens[i] / bs * bs == jac->g_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat blocksize doesn't match block Jacobi layout"); 54 } 55 if (size == 1) { 56 jac->n_local = jac->n; 57 PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 58 PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens, jac->n_local)); 59 /* check that user set these correctly */ 60 sum = 0; 61 for (i = 0; i < jac->n_local; i++) sum += jac->l_lens[i]; 62 PetscCheck(sum == M, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Global lens set incorrectly"); 63 } else { 64 PetscCall(MatGetOwnershipRange(pc->pmat, &start, &end)); 65 /* loop over blocks determing first one owned by me */ 66 sum = 0; 67 for (i = 0; i < jac->n + 1; i++) { 68 if (sum == start) { 69 i_start = i; 70 goto start_1; 71 } 72 if (i < jac->n) sum += jac->g_lens[i]; 73 } 74 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout"); 75 start_1: 76 for (i = i_start; i < jac->n + 1; i++) { 77 if (sum == end) { 78 i_end = i; 79 goto end_1; 80 } 81 if (i < jac->n) sum += jac->g_lens[i]; 82 } 83 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Block sizes used in PCBJacobiSetTotalBlocks()\nare not compatible with parallel matrix layout"); 84 end_1: 85 jac->n_local = i_end - i_start; 86 PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 87 PetscCall(PetscArraycpy(jac->l_lens, jac->g_lens + i_start, jac->n_local)); 88 } 89 } else { /* no global blocks given, determine then using default layout */ 90 jac->n_local = jac->n / size + ((jac->n % size) > rank); 91 PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 92 for (i = 0; i < jac->n_local; i++) { 93 jac->l_lens[i] = ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)) * bs; 94 PetscCheck(jac->l_lens[i], PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Too many blocks given"); 95 } 96 } 97 } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */ 98 jac->n = size; 99 jac->n_local = 1; 100 PetscCall(PetscMalloc1(1, &jac->l_lens)); 101 jac->l_lens[0] = M; 102 } else { /* jac->n > 0 && jac->n_local > 0 */ 103 if (!jac->l_lens) { 104 PetscCall(PetscMalloc1(jac->n_local, &jac->l_lens)); 105 for (i = 0; i < jac->n_local; i++) jac->l_lens[i] = bs * ((M / bs) / jac->n_local + (((M / bs) % jac->n_local) > i)); 106 } 107 } 108 PetscCheck(jac->n_local >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of blocks is less than number of processors"); 109 110 /* Determines mat and pmat */ 111 PetscCall(MatHasOperation(pc->mat, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 112 if (!hasop && size == 1) { 113 mat = pc->mat; 114 pmat = pc->pmat; 115 } else { 116 if (pc->useAmat) { 117 /* use block from Amat matrix, not Pmat for local MatMult() */ 118 PetscCall(MatGetDiagonalBlock(pc->mat, &mat)); 119 } 120 if (pc->pmat != pc->mat || !pc->useAmat) { 121 PetscCall(MatGetDiagonalBlock(pc->pmat, &pmat)); 122 } else pmat = mat; 123 } 124 125 /* 126 Setup code depends on the number of blocks 127 */ 128 if (jac->n_local == 1) { 129 PetscCall(PCSetUp_BJacobi_Singleblock(pc, mat, pmat)); 130 } else { 131 PetscCall(PCSetUp_BJacobi_Multiblock(pc, mat, pmat)); 132 } 133 PetscFunctionReturn(0); 134 } 135 136 /* Default destroy, if it has never been setup */ 137 static PetscErrorCode PCDestroy_BJacobi(PC pc) 138 { 139 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 140 141 PetscFunctionBegin; 142 PetscCall(PetscFree(jac->g_lens)); 143 PetscCall(PetscFree(jac->l_lens)); 144 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", NULL)); 145 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", NULL)); 146 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", NULL)); 147 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", NULL)); 148 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", NULL)); 149 PetscCall(PetscFree(pc->data)); 150 PetscFunctionReturn(0); 151 } 152 153 static PetscErrorCode PCSetFromOptions_BJacobi(PC pc, PetscOptionItems *PetscOptionsObject) 154 { 155 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 156 PetscInt blocks, i; 157 PetscBool flg; 158 159 PetscFunctionBegin; 160 PetscOptionsHeadBegin(PetscOptionsObject, "Block Jacobi options"); 161 PetscCall(PetscOptionsInt("-pc_bjacobi_blocks", "Total number of blocks", "PCBJacobiSetTotalBlocks", jac->n, &blocks, &flg)); 162 if (flg) PetscCall(PCBJacobiSetTotalBlocks(pc, blocks, NULL)); 163 PetscCall(PetscOptionsInt("-pc_bjacobi_local_blocks", "Local number of blocks", "PCBJacobiSetLocalBlocks", jac->n_local, &blocks, &flg)); 164 if (flg) PetscCall(PCBJacobiSetLocalBlocks(pc, blocks, NULL)); 165 if (jac->ksp) { 166 /* The sub-KSP has already been set up (e.g., PCSetUp_BJacobi_Singleblock), but KSPSetFromOptions was not called 167 * unless we had already been called. */ 168 for (i = 0; i < jac->n_local; i++) PetscCall(KSPSetFromOptions(jac->ksp[i])); 169 } 170 PetscOptionsHeadEnd(); 171 PetscFunctionReturn(0); 172 } 173 174 #include <petscdraw.h> 175 static PetscErrorCode PCView_BJacobi(PC pc, PetscViewer viewer) 176 { 177 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 178 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 179 PetscMPIInt rank; 180 PetscInt i; 181 PetscBool iascii, isstring, isdraw; 182 PetscViewer sviewer; 183 PetscViewerFormat format; 184 const char *prefix; 185 186 PetscFunctionBegin; 187 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 188 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 189 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 190 if (iascii) { 191 if (pc->useAmat) PetscCall(PetscViewerASCIIPrintf(viewer, " using Amat local matrix, number of blocks = %" PetscInt_FMT "\n", jac->n)); 192 PetscCall(PetscViewerASCIIPrintf(viewer, " number of blocks = %" PetscInt_FMT "\n", jac->n)); 193 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 194 PetscCall(PetscViewerGetFormat(viewer, &format)); 195 if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 196 PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block is in the following KSP and PC objects on rank 0:\n")); 197 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 198 PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 199 if (jac->ksp && !jac->psubcomm) { 200 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 201 if (rank == 0) { 202 PetscCall(PetscViewerASCIIPushTab(viewer)); 203 PetscCall(KSPView(jac->ksp[0], sviewer)); 204 PetscCall(PetscViewerASCIIPopTab(viewer)); 205 } 206 PetscCall(PetscViewerFlush(sviewer)); 207 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 208 PetscCall(PetscViewerFlush(viewer)); 209 /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 210 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 211 } else if (mpjac && jac->ksp && mpjac->psubcomm) { 212 PetscCall(PetscViewerGetSubViewer(viewer, mpjac->psubcomm->child, &sviewer)); 213 if (!mpjac->psubcomm->color) { 214 PetscCall(PetscViewerASCIIPushTab(viewer)); 215 PetscCall(KSPView(*(jac->ksp), sviewer)); 216 PetscCall(PetscViewerASCIIPopTab(viewer)); 217 } 218 PetscCall(PetscViewerFlush(sviewer)); 219 PetscCall(PetscViewerRestoreSubViewer(viewer, mpjac->psubcomm->child, &sviewer)); 220 PetscCall(PetscViewerFlush(viewer)); 221 /* extra call needed because of the two calls to PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 222 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 223 } else { 224 PetscCall(PetscViewerFlush(viewer)); 225 } 226 } else { 227 PetscInt n_global; 228 PetscCall(MPIU_Allreduce(&jac->n_local, &n_global, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc))); 229 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 230 PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following KSP and PC objects:\n")); 231 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] number of local blocks = %" PetscInt_FMT ", first local block number = %" PetscInt_FMT "\n", rank, jac->n_local, jac->first_local)); 232 PetscCall(PetscViewerASCIIPushTab(viewer)); 233 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 234 for (i = 0; i < jac->n_local; i++) { 235 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] local block number %" PetscInt_FMT "\n", rank, i)); 236 PetscCall(KSPView(jac->ksp[i], sviewer)); 237 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 238 } 239 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 240 PetscCall(PetscViewerASCIIPopTab(viewer)); 241 PetscCall(PetscViewerFlush(viewer)); 242 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 243 } 244 } else if (isstring) { 245 PetscCall(PetscViewerStringSPrintf(viewer, " blks=%" PetscInt_FMT, jac->n)); 246 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 247 if (jac->ksp) PetscCall(KSPView(jac->ksp[0], sviewer)); 248 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 249 } else if (isdraw) { 250 PetscDraw draw; 251 char str[25]; 252 PetscReal x, y, bottom, h; 253 254 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 255 PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 256 PetscCall(PetscSNPrintf(str, 25, "Number blocks %" PetscInt_FMT, jac->n)); 257 PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 258 bottom = y - h; 259 PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 260 /* warning the communicator on viewer is different then on ksp in parallel */ 261 if (jac->ksp) PetscCall(KSPView(jac->ksp[0], viewer)); 262 PetscCall(PetscDrawPopCurrentPoint(draw)); 263 } 264 PetscFunctionReturn(0); 265 } 266 267 static PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) 268 { 269 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 270 271 PetscFunctionBegin; 272 PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call KSPSetUp() or PCSetUp() first"); 273 274 if (n_local) *n_local = jac->n_local; 275 if (first_local) *first_local = jac->first_local; 276 if (ksp) *ksp = jac->ksp; 277 PetscFunctionReturn(0); 278 } 279 280 static PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc, PetscInt blocks, PetscInt *lens) 281 { 282 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 283 284 PetscFunctionBegin; 285 PetscCheck(pc->setupcalled <= 0 || jac->n == blocks, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Cannot alter number of blocks after PCSetUp()/KSPSetUp() has been called"); 286 jac->n = blocks; 287 if (!lens) jac->g_lens = NULL; 288 else { 289 PetscCall(PetscMalloc1(blocks, &jac->g_lens)); 290 PetscCall(PetscArraycpy(jac->g_lens, lens, blocks)); 291 } 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) 296 { 297 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 298 299 PetscFunctionBegin; 300 *blocks = jac->n; 301 if (lens) *lens = jac->g_lens; 302 PetscFunctionReturn(0); 303 } 304 305 static PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc, PetscInt blocks, const PetscInt lens[]) 306 { 307 PC_BJacobi *jac; 308 309 PetscFunctionBegin; 310 jac = (PC_BJacobi *)pc->data; 311 312 jac->n_local = blocks; 313 if (!lens) jac->l_lens = NULL; 314 else { 315 PetscCall(PetscMalloc1(blocks, &jac->l_lens)); 316 PetscCall(PetscArraycpy(jac->l_lens, lens, blocks)); 317 } 318 PetscFunctionReturn(0); 319 } 320 321 static PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) 322 { 323 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 324 325 PetscFunctionBegin; 326 *blocks = jac->n_local; 327 if (lens) *lens = jac->l_lens; 328 PetscFunctionReturn(0); 329 } 330 331 /*@C 332 PCBJacobiGetSubKSP - Gets the local `KSP` contexts for all blocks on 333 this processor. 334 335 Not Collective 336 337 Input Parameter: 338 . pc - the preconditioner context 339 340 Output Parameters: 341 + n_local - the number of blocks on this processor, or NULL 342 . first_local - the global number of the first block on this processor, or NULL 343 - ksp - the array of KSP contexts 344 345 Notes: 346 After `PCBJacobiGetSubKSP()` the array of `KSP` contexts is not to be freed. 347 348 Currently for some matrix implementations only 1 block per processor 349 is supported. 350 351 You must call `KSPSetUp()` or `PCSetUp()` before calling `PCBJacobiGetSubKSP()`. 352 353 Fortran Usage: 354 You must pass in a `KSP` array that is large enough to contain all the local `KSP`s. 355 356 You can call `PCBJacobiGetSubKSP`(pc,nlocal,firstlocal,`PETSC_NULL_KSP`,ierr) to determine how large the 357 `KSP` array must be. 358 359 Level: advanced 360 361 .seealso: `PCBJACOBI`, `PCASM`, `PCASMGetSubKSP()` 362 @*/ 363 PetscErrorCode PCBJacobiGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) 364 { 365 PetscFunctionBegin; 366 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 367 PetscUseMethod(pc, "PCBJacobiGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp)); 368 PetscFunctionReturn(0); 369 } 370 371 /*@ 372 PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block 373 Jacobi preconditioner. 374 375 Collective on pc 376 377 Input Parameters: 378 + pc - the preconditioner context 379 . blocks - the number of blocks 380 - lens - [optional] integer array containing the size of each block 381 382 Options Database Key: 383 . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks 384 385 Note: 386 Currently only a limited number of blocking configurations are supported. 387 All processors sharing the `PC` must call this routine with the same data. 388 389 Level: intermediate 390 391 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetLocalBlocks()` 392 @*/ 393 PetscErrorCode PCBJacobiSetTotalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) 394 { 395 PetscFunctionBegin; 396 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 397 PetscCheck(blocks > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Must have positive blocks"); 398 PetscTryMethod(pc, "PCBJacobiSetTotalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens)); 399 PetscFunctionReturn(0); 400 } 401 402 /*@C 403 PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block 404 Jacobi, `PCBJACOBI`, preconditioner. 405 406 Not Collective 407 408 Input Parameter: 409 . pc - the preconditioner context 410 411 Output parameters: 412 + blocks - the number of blocks 413 - lens - integer array containing the size of each block 414 415 Level: intermediate 416 417 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetLocalBlocks()` 418 @*/ 419 PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) 420 { 421 PetscFunctionBegin; 422 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 423 PetscValidIntPointer(blocks, 2); 424 PetscUseMethod(pc, "PCBJacobiGetTotalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens)); 425 PetscFunctionReturn(0); 426 } 427 428 /*@ 429 PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block 430 Jacobi, `PCBJACOBI`, preconditioner. 431 432 Not Collective 433 434 Input Parameters: 435 + pc - the preconditioner context 436 . blocks - the number of blocks 437 - lens - [optional] integer array containing size of each block 438 439 Options Database Key: 440 . -pc_bjacobi_local_blocks <blocks> - Sets the number of local blocks 441 442 Note: 443 Currently only a limited number of blocking configurations are supported. 444 445 Level: intermediate 446 447 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiSetTotalBlocks()` 448 @*/ 449 PetscErrorCode PCBJacobiSetLocalBlocks(PC pc, PetscInt blocks, const PetscInt lens[]) 450 { 451 PetscFunctionBegin; 452 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 453 PetscCheck(blocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have nonegative blocks"); 454 PetscTryMethod(pc, "PCBJacobiSetLocalBlocks_C", (PC, PetscInt, const PetscInt[]), (pc, blocks, lens)); 455 PetscFunctionReturn(0); 456 } 457 458 /*@C 459 PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block 460 Jacobi, `PCBJACOBI`, preconditioner. 461 462 Not Collective 463 464 Input Parameters: 465 + pc - the preconditioner context 466 . blocks - the number of blocks 467 - lens - [optional] integer array containing size of each block 468 469 Note: 470 Currently only a limited number of blocking configurations are supported. 471 472 Level: intermediate 473 474 .seealso: `PCBJACOBI`, `PCSetUseAmat()`, `PCBJacobiGetTotalBlocks()` 475 @*/ 476 PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) 477 { 478 PetscFunctionBegin; 479 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 480 PetscValidIntPointer(blocks, 2); 481 PetscUseMethod(pc, "PCBJacobiGetLocalBlocks_C", (PC, PetscInt *, const PetscInt *[]), (pc, blocks, lens)); 482 PetscFunctionReturn(0); 483 } 484 485 /*MC 486 PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with 487 its own `KSP` object. 488 489 Options Database Keys: 490 + -pc_use_amat - use Amat to apply block of operator in inner Krylov method 491 - -pc_bjacobi_blocks <n> - use n total blocks 492 493 Notes: 494 See `PCJACOBI` for diagonal Jacobi, `PCVPBJACOBI` for variable point block, and `PCPBJACOBI` for fixed size point block 495 496 Each processor can have one or more blocks, or a single block can be shared by several processes. Defaults to one block per processor. 497 498 To set options on the solvers for each block append -sub_ to all the `KSP` and `PC` 499 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 500 501 To set the options on the solvers separate for each block call `PCBJacobiGetSubKSP()` 502 and set the options directly on the resulting `KSP` object (you can access its `PC` 503 `KSPGetPC())` 504 505 For GPU-based vectors (`VECCUDA`, `VECViennaCL`) it is recommended to use exactly one block per MPI process for best 506 performance. Different block partitioning may lead to additional data transfers 507 between host and GPU that lead to degraded performance. 508 509 When multiple processes share a single block, each block encompasses exactly all the unknowns owned its set of processes. 510 511 Level: beginner 512 513 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCType`, 514 `PCASM`, `PCSetUseAmat()`, `PCGetUseAmat()`, `PCBJacobiGetSubKSP()`, `PCBJacobiSetTotalBlocks()`, 515 `PCBJacobiSetLocalBlocks()`, `PCSetModifySubMatrices()`, `PCJACOBI`, `PCVPBJACOBI`, `PCPBJACOBI` 516 M*/ 517 518 PETSC_EXTERN PetscErrorCode PCCreate_BJacobi(PC pc) 519 { 520 PetscMPIInt rank; 521 PC_BJacobi *jac; 522 523 PetscFunctionBegin; 524 PetscCall(PetscNew(&jac)); 525 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 526 527 pc->ops->apply = NULL; 528 pc->ops->matapply = NULL; 529 pc->ops->applytranspose = NULL; 530 pc->ops->setup = PCSetUp_BJacobi; 531 pc->ops->destroy = PCDestroy_BJacobi; 532 pc->ops->setfromoptions = PCSetFromOptions_BJacobi; 533 pc->ops->view = PCView_BJacobi; 534 pc->ops->applyrichardson = NULL; 535 536 pc->data = (void *)jac; 537 jac->n = -1; 538 jac->n_local = -1; 539 jac->first_local = rank; 540 jac->ksp = NULL; 541 jac->g_lens = NULL; 542 jac->l_lens = NULL; 543 jac->psubcomm = NULL; 544 545 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetSubKSP_C", PCBJacobiGetSubKSP_BJacobi)); 546 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetTotalBlocks_C", PCBJacobiSetTotalBlocks_BJacobi)); 547 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetTotalBlocks_C", PCBJacobiGetTotalBlocks_BJacobi)); 548 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiSetLocalBlocks_C", PCBJacobiSetLocalBlocks_BJacobi)); 549 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJacobiGetLocalBlocks_C", PCBJacobiGetLocalBlocks_BJacobi)); 550 PetscFunctionReturn(0); 551 } 552 553 /* 554 These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI 555 */ 556 static PetscErrorCode PCReset_BJacobi_Singleblock(PC pc) 557 { 558 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 559 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 560 561 PetscFunctionBegin; 562 PetscCall(KSPReset(jac->ksp[0])); 563 PetscCall(VecDestroy(&bjac->x)); 564 PetscCall(VecDestroy(&bjac->y)); 565 PetscFunctionReturn(0); 566 } 567 568 static PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc) 569 { 570 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 571 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 572 573 PetscFunctionBegin; 574 PetscCall(PCReset_BJacobi_Singleblock(pc)); 575 PetscCall(KSPDestroy(&jac->ksp[0])); 576 PetscCall(PetscFree(jac->ksp)); 577 PetscCall(PetscFree(bjac)); 578 PetscCall(PCDestroy_BJacobi(pc)); 579 PetscFunctionReturn(0); 580 } 581 582 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc) 583 { 584 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 585 KSP subksp = jac->ksp[0]; 586 KSPConvergedReason reason; 587 588 PetscFunctionBegin; 589 PetscCall(KSPSetUp(subksp)); 590 PetscCall(KSPGetConvergedReason(subksp, &reason)); 591 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 592 PetscFunctionReturn(0); 593 } 594 595 static PetscErrorCode PCApply_BJacobi_Singleblock(PC pc, Vec x, Vec y) 596 { 597 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 598 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 599 600 PetscFunctionBegin; 601 PetscCall(VecGetLocalVectorRead(x, bjac->x)); 602 PetscCall(VecGetLocalVector(y, bjac->y)); 603 /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner 604 matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild 605 of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/ 606 PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner)); 607 PetscCall(KSPSolve(jac->ksp[0], bjac->x, bjac->y)); 608 PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y)); 609 PetscCall(VecRestoreLocalVectorRead(x, bjac->x)); 610 PetscCall(VecRestoreLocalVector(y, bjac->y)); 611 PetscFunctionReturn(0); 612 } 613 614 static PetscErrorCode PCMatApply_BJacobi_Singleblock(PC pc, Mat X, Mat Y) 615 { 616 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 617 Mat sX, sY; 618 619 PetscFunctionBegin; 620 /* Since the inner KSP matrix may point directly to the diagonal block of an MPI matrix the inner 621 matrix may change even if the outer KSP/PC has not updated the preconditioner, this will trigger a rebuild 622 of the inner preconditioner automatically unless we pass down the outer preconditioners reuse flag.*/ 623 PetscCall(KSPSetReusePreconditioner(jac->ksp[0], pc->reusepreconditioner)); 624 PetscCall(MatDenseGetLocalMatrix(X, &sX)); 625 PetscCall(MatDenseGetLocalMatrix(Y, &sY)); 626 PetscCall(KSPMatSolve(jac->ksp[0], sX, sY)); 627 PetscFunctionReturn(0); 628 } 629 630 static PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc, Vec x, Vec y) 631 { 632 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 633 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 634 PetscScalar *y_array; 635 const PetscScalar *x_array; 636 PC subpc; 637 638 PetscFunctionBegin; 639 /* 640 The VecPlaceArray() is to avoid having to copy the 641 y vector into the bjac->x vector. The reason for 642 the bjac->x vector is that we need a sequential vector 643 for the sequential solve. 644 */ 645 PetscCall(VecGetArrayRead(x, &x_array)); 646 PetscCall(VecGetArray(y, &y_array)); 647 PetscCall(VecPlaceArray(bjac->x, x_array)); 648 PetscCall(VecPlaceArray(bjac->y, y_array)); 649 /* apply the symmetric left portion of the inner PC operator */ 650 /* note this by-passes the inner KSP and its options completely */ 651 PetscCall(KSPGetPC(jac->ksp[0], &subpc)); 652 PetscCall(PCApplySymmetricLeft(subpc, bjac->x, bjac->y)); 653 PetscCall(VecResetArray(bjac->x)); 654 PetscCall(VecResetArray(bjac->y)); 655 PetscCall(VecRestoreArrayRead(x, &x_array)); 656 PetscCall(VecRestoreArray(y, &y_array)); 657 PetscFunctionReturn(0); 658 } 659 660 static PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc, Vec x, Vec y) 661 { 662 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 663 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 664 PetscScalar *y_array; 665 const PetscScalar *x_array; 666 PC subpc; 667 668 PetscFunctionBegin; 669 /* 670 The VecPlaceArray() is to avoid having to copy the 671 y vector into the bjac->x vector. The reason for 672 the bjac->x vector is that we need a sequential vector 673 for the sequential solve. 674 */ 675 PetscCall(VecGetArrayRead(x, &x_array)); 676 PetscCall(VecGetArray(y, &y_array)); 677 PetscCall(VecPlaceArray(bjac->x, x_array)); 678 PetscCall(VecPlaceArray(bjac->y, y_array)); 679 680 /* apply the symmetric right portion of the inner PC operator */ 681 /* note this by-passes the inner KSP and its options completely */ 682 683 PetscCall(KSPGetPC(jac->ksp[0], &subpc)); 684 PetscCall(PCApplySymmetricRight(subpc, bjac->x, bjac->y)); 685 686 PetscCall(VecRestoreArrayRead(x, &x_array)); 687 PetscCall(VecRestoreArray(y, &y_array)); 688 PetscFunctionReturn(0); 689 } 690 691 static PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc, Vec x, Vec y) 692 { 693 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 694 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock *)jac->data; 695 PetscScalar *y_array; 696 const PetscScalar *x_array; 697 698 PetscFunctionBegin; 699 /* 700 The VecPlaceArray() is to avoid having to copy the 701 y vector into the bjac->x vector. The reason for 702 the bjac->x vector is that we need a sequential vector 703 for the sequential solve. 704 */ 705 PetscCall(VecGetArrayRead(x, &x_array)); 706 PetscCall(VecGetArray(y, &y_array)); 707 PetscCall(VecPlaceArray(bjac->x, x_array)); 708 PetscCall(VecPlaceArray(bjac->y, y_array)); 709 PetscCall(KSPSolveTranspose(jac->ksp[0], bjac->x, bjac->y)); 710 PetscCall(KSPCheckSolve(jac->ksp[0], pc, bjac->y)); 711 PetscCall(VecResetArray(bjac->x)); 712 PetscCall(VecResetArray(bjac->y)); 713 PetscCall(VecRestoreArrayRead(x, &x_array)); 714 PetscCall(VecRestoreArray(y, &y_array)); 715 PetscFunctionReturn(0); 716 } 717 718 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc, Mat mat, Mat pmat) 719 { 720 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 721 PetscInt m; 722 KSP ksp; 723 PC_BJacobi_Singleblock *bjac; 724 PetscBool wasSetup = PETSC_TRUE; 725 VecType vectype; 726 const char *prefix; 727 728 PetscFunctionBegin; 729 if (!pc->setupcalled) { 730 if (!jac->ksp) { 731 wasSetup = PETSC_FALSE; 732 733 PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 734 PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure)); 735 PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1)); 736 PetscCall(KSPSetType(ksp, KSPPREONLY)); 737 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 738 PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 739 PetscCall(KSPAppendOptionsPrefix(ksp, "sub_")); 740 741 pc->ops->reset = PCReset_BJacobi_Singleblock; 742 pc->ops->destroy = PCDestroy_BJacobi_Singleblock; 743 pc->ops->apply = PCApply_BJacobi_Singleblock; 744 pc->ops->matapply = PCMatApply_BJacobi_Singleblock; 745 pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock; 746 pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock; 747 pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock; 748 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock; 749 750 PetscCall(PetscMalloc1(1, &jac->ksp)); 751 jac->ksp[0] = ksp; 752 753 PetscCall(PetscNew(&bjac)); 754 jac->data = (void *)bjac; 755 } else { 756 ksp = jac->ksp[0]; 757 bjac = (PC_BJacobi_Singleblock *)jac->data; 758 } 759 760 /* 761 The reason we need to generate these vectors is to serve 762 as the right-hand side and solution vector for the solve on the 763 block. We do not need to allocate space for the vectors since 764 that is provided via VecPlaceArray() just before the call to 765 KSPSolve() on the block. 766 */ 767 PetscCall(MatGetSize(pmat, &m, &m)); 768 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->x)); 769 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &bjac->y)); 770 PetscCall(MatGetVecType(pmat, &vectype)); 771 PetscCall(VecSetType(bjac->x, vectype)); 772 PetscCall(VecSetType(bjac->y, vectype)); 773 } else { 774 ksp = jac->ksp[0]; 775 bjac = (PC_BJacobi_Singleblock *)jac->data; 776 } 777 PetscCall(KSPGetOptionsPrefix(ksp, &prefix)); 778 if (pc->useAmat) { 779 PetscCall(KSPSetOperators(ksp, mat, pmat)); 780 PetscCall(MatSetOptionsPrefix(mat, prefix)); 781 } else { 782 PetscCall(KSPSetOperators(ksp, pmat, pmat)); 783 } 784 PetscCall(MatSetOptionsPrefix(pmat, prefix)); 785 if (!wasSetup && pc->setfromoptionscalled) { 786 /* If PCSetFromOptions_BJacobi is called later, KSPSetFromOptions will be called at that time. */ 787 PetscCall(KSPSetFromOptions(ksp)); 788 } 789 PetscFunctionReturn(0); 790 } 791 792 static PetscErrorCode PCReset_BJacobi_Multiblock(PC pc) 793 { 794 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 795 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 796 PetscInt i; 797 798 PetscFunctionBegin; 799 if (bjac && bjac->pmat) { 800 PetscCall(MatDestroyMatrices(jac->n_local, &bjac->pmat)); 801 if (pc->useAmat) PetscCall(MatDestroyMatrices(jac->n_local, &bjac->mat)); 802 } 803 804 for (i = 0; i < jac->n_local; i++) { 805 PetscCall(KSPReset(jac->ksp[i])); 806 if (bjac && bjac->x) { 807 PetscCall(VecDestroy(&bjac->x[i])); 808 PetscCall(VecDestroy(&bjac->y[i])); 809 PetscCall(ISDestroy(&bjac->is[i])); 810 } 811 } 812 PetscCall(PetscFree(jac->l_lens)); 813 PetscCall(PetscFree(jac->g_lens)); 814 PetscFunctionReturn(0); 815 } 816 817 static PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc) 818 { 819 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 820 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 821 PetscInt i; 822 823 PetscFunctionBegin; 824 PetscCall(PCReset_BJacobi_Multiblock(pc)); 825 if (bjac) { 826 PetscCall(PetscFree2(bjac->x, bjac->y)); 827 PetscCall(PetscFree(bjac->starts)); 828 PetscCall(PetscFree(bjac->is)); 829 } 830 PetscCall(PetscFree(jac->data)); 831 for (i = 0; i < jac->n_local; i++) PetscCall(KSPDestroy(&jac->ksp[i])); 832 PetscCall(PetscFree(jac->ksp)); 833 PetscCall(PCDestroy_BJacobi(pc)); 834 PetscFunctionReturn(0); 835 } 836 837 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc) 838 { 839 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 840 PetscInt i, n_local = jac->n_local; 841 KSPConvergedReason reason; 842 843 PetscFunctionBegin; 844 for (i = 0; i < n_local; i++) { 845 PetscCall(KSPSetUp(jac->ksp[i])); 846 PetscCall(KSPGetConvergedReason(jac->ksp[i], &reason)); 847 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 848 } 849 PetscFunctionReturn(0); 850 } 851 852 /* 853 Preconditioner for block Jacobi 854 */ 855 static PetscErrorCode PCApply_BJacobi_Multiblock(PC pc, Vec x, Vec y) 856 { 857 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 858 PetscInt i, n_local = jac->n_local; 859 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 860 PetscScalar *yin; 861 const PetscScalar *xin; 862 863 PetscFunctionBegin; 864 PetscCall(VecGetArrayRead(x, &xin)); 865 PetscCall(VecGetArray(y, &yin)); 866 for (i = 0; i < n_local; i++) { 867 /* 868 To avoid copying the subvector from x into a workspace we instead 869 make the workspace vector array point to the subpart of the array of 870 the global vector. 871 */ 872 PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i])); 873 PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i])); 874 875 PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 876 PetscCall(KSPSolve(jac->ksp[i], bjac->x[i], bjac->y[i])); 877 PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i])); 878 PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 879 880 PetscCall(VecResetArray(bjac->x[i])); 881 PetscCall(VecResetArray(bjac->y[i])); 882 } 883 PetscCall(VecRestoreArrayRead(x, &xin)); 884 PetscCall(VecRestoreArray(y, &yin)); 885 PetscFunctionReturn(0); 886 } 887 888 /* 889 Preconditioner for block Jacobi 890 */ 891 static PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc, Vec x, Vec y) 892 { 893 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 894 PetscInt i, n_local = jac->n_local; 895 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 896 PetscScalar *yin; 897 const PetscScalar *xin; 898 899 PetscFunctionBegin; 900 PetscCall(VecGetArrayRead(x, &xin)); 901 PetscCall(VecGetArray(y, &yin)); 902 for (i = 0; i < n_local; i++) { 903 /* 904 To avoid copying the subvector from x into a workspace we instead 905 make the workspace vector array point to the subpart of the array of 906 the global vector. 907 */ 908 PetscCall(VecPlaceArray(bjac->x[i], xin + bjac->starts[i])); 909 PetscCall(VecPlaceArray(bjac->y[i], yin + bjac->starts[i])); 910 911 PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 912 PetscCall(KSPSolveTranspose(jac->ksp[i], bjac->x[i], bjac->y[i])); 913 PetscCall(KSPCheckSolve(jac->ksp[i], pc, bjac->y[i])); 914 PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, jac->ksp[i], bjac->x[i], bjac->y[i], 0)); 915 916 PetscCall(VecResetArray(bjac->x[i])); 917 PetscCall(VecResetArray(bjac->y[i])); 918 } 919 PetscCall(VecRestoreArrayRead(x, &xin)); 920 PetscCall(VecRestoreArray(y, &yin)); 921 PetscFunctionReturn(0); 922 } 923 924 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc, Mat mat, Mat pmat) 925 { 926 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 927 PetscInt m, n_local, N, M, start, i; 928 const char *prefix; 929 KSP ksp; 930 Vec x, y; 931 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock *)jac->data; 932 PC subpc; 933 IS is; 934 MatReuse scall; 935 VecType vectype; 936 937 PetscFunctionBegin; 938 PetscCall(MatGetLocalSize(pc->pmat, &M, &N)); 939 940 n_local = jac->n_local; 941 942 if (pc->useAmat) { 943 PetscBool same; 944 PetscCall(PetscObjectTypeCompare((PetscObject)mat, ((PetscObject)pmat)->type_name, &same)); 945 PetscCheck(same, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Matrices not of same type"); 946 } 947 948 if (!pc->setupcalled) { 949 scall = MAT_INITIAL_MATRIX; 950 951 if (!jac->ksp) { 952 pc->ops->reset = PCReset_BJacobi_Multiblock; 953 pc->ops->destroy = PCDestroy_BJacobi_Multiblock; 954 pc->ops->apply = PCApply_BJacobi_Multiblock; 955 pc->ops->matapply = NULL; 956 pc->ops->applytranspose = PCApplyTranspose_BJacobi_Multiblock; 957 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock; 958 959 PetscCall(PetscNew(&bjac)); 960 PetscCall(PetscMalloc1(n_local, &jac->ksp)); 961 PetscCall(PetscMalloc2(n_local, &bjac->x, n_local, &bjac->y)); 962 PetscCall(PetscMalloc1(n_local, &bjac->starts)); 963 964 jac->data = (void *)bjac; 965 PetscCall(PetscMalloc1(n_local, &bjac->is)); 966 967 for (i = 0; i < n_local; i++) { 968 PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 969 PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure)); 970 PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1)); 971 PetscCall(KSPSetType(ksp, KSPPREONLY)); 972 PetscCall(KSPGetPC(ksp, &subpc)); 973 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 974 PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 975 PetscCall(KSPAppendOptionsPrefix(ksp, "sub_")); 976 977 jac->ksp[i] = ksp; 978 } 979 } else { 980 bjac = (PC_BJacobi_Multiblock *)jac->data; 981 } 982 983 start = 0; 984 PetscCall(MatGetVecType(pmat, &vectype)); 985 for (i = 0; i < n_local; i++) { 986 m = jac->l_lens[i]; 987 /* 988 The reason we need to generate these vectors is to serve 989 as the right-hand side and solution vector for the solve on the 990 block. We do not need to allocate space for the vectors since 991 that is provided via VecPlaceArray() just before the call to 992 KSPSolve() on the block. 993 994 */ 995 PetscCall(VecCreateSeq(PETSC_COMM_SELF, m, &x)); 996 PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, m, NULL, &y)); 997 PetscCall(VecSetType(x, vectype)); 998 PetscCall(VecSetType(y, vectype)); 999 1000 bjac->x[i] = x; 1001 bjac->y[i] = y; 1002 bjac->starts[i] = start; 1003 1004 PetscCall(ISCreateStride(PETSC_COMM_SELF, m, start, 1, &is)); 1005 bjac->is[i] = is; 1006 1007 start += m; 1008 } 1009 } else { 1010 bjac = (PC_BJacobi_Multiblock *)jac->data; 1011 /* 1012 Destroy the blocks from the previous iteration 1013 */ 1014 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1015 PetscCall(MatDestroyMatrices(n_local, &bjac->pmat)); 1016 if (pc->useAmat) PetscCall(MatDestroyMatrices(n_local, &bjac->mat)); 1017 scall = MAT_INITIAL_MATRIX; 1018 } else scall = MAT_REUSE_MATRIX; 1019 } 1020 1021 PetscCall(MatCreateSubMatrices(pmat, n_local, bjac->is, bjac->is, scall, &bjac->pmat)); 1022 if (pc->useAmat) PetscCall(MatCreateSubMatrices(mat, n_local, bjac->is, bjac->is, scall, &bjac->mat)); 1023 /* Return control to the user so that the submatrices can be modified (e.g., to apply 1024 different boundary conditions for the submatrices than for the global problem) */ 1025 PetscCall(PCModifySubMatrices(pc, n_local, bjac->is, bjac->is, bjac->pmat, pc->modifysubmatricesP)); 1026 1027 for (i = 0; i < n_local; i++) { 1028 PetscCall(KSPGetOptionsPrefix(jac->ksp[i], &prefix)); 1029 if (pc->useAmat) { 1030 PetscCall(KSPSetOperators(jac->ksp[i], bjac->mat[i], bjac->pmat[i])); 1031 PetscCall(MatSetOptionsPrefix(bjac->mat[i], prefix)); 1032 } else { 1033 PetscCall(KSPSetOperators(jac->ksp[i], bjac->pmat[i], bjac->pmat[i])); 1034 } 1035 PetscCall(MatSetOptionsPrefix(bjac->pmat[i], prefix)); 1036 if (pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[i])); 1037 } 1038 PetscFunctionReturn(0); 1039 } 1040 1041 /* 1042 These are for a single block with multiple processes 1043 */ 1044 static PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiproc(PC pc) 1045 { 1046 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1047 KSP subksp = jac->ksp[0]; 1048 KSPConvergedReason reason; 1049 1050 PetscFunctionBegin; 1051 PetscCall(KSPSetUp(subksp)); 1052 PetscCall(KSPGetConvergedReason(subksp, &reason)); 1053 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1054 PetscFunctionReturn(0); 1055 } 1056 1057 static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc) 1058 { 1059 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1060 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 1061 1062 PetscFunctionBegin; 1063 PetscCall(VecDestroy(&mpjac->ysub)); 1064 PetscCall(VecDestroy(&mpjac->xsub)); 1065 PetscCall(MatDestroy(&mpjac->submats)); 1066 if (jac->ksp) PetscCall(KSPReset(jac->ksp[0])); 1067 PetscFunctionReturn(0); 1068 } 1069 1070 static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc) 1071 { 1072 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1073 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 1074 1075 PetscFunctionBegin; 1076 PetscCall(PCReset_BJacobi_Multiproc(pc)); 1077 PetscCall(KSPDestroy(&jac->ksp[0])); 1078 PetscCall(PetscFree(jac->ksp)); 1079 PetscCall(PetscSubcommDestroy(&mpjac->psubcomm)); 1080 1081 PetscCall(PetscFree(mpjac)); 1082 PetscCall(PCDestroy_BJacobi(pc)); 1083 PetscFunctionReturn(0); 1084 } 1085 1086 static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc, Vec x, Vec y) 1087 { 1088 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1089 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 1090 PetscScalar *yarray; 1091 const PetscScalar *xarray; 1092 KSPConvergedReason reason; 1093 1094 PetscFunctionBegin; 1095 /* place x's and y's local arrays into xsub and ysub */ 1096 PetscCall(VecGetArrayRead(x, &xarray)); 1097 PetscCall(VecGetArray(y, &yarray)); 1098 PetscCall(VecPlaceArray(mpjac->xsub, xarray)); 1099 PetscCall(VecPlaceArray(mpjac->ysub, yarray)); 1100 1101 /* apply preconditioner on each matrix block */ 1102 PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0)); 1103 PetscCall(KSPSolve(jac->ksp[0], mpjac->xsub, mpjac->ysub)); 1104 PetscCall(KSPCheckSolve(jac->ksp[0], pc, mpjac->ysub)); 1105 PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], mpjac->xsub, mpjac->ysub, 0)); 1106 PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason)); 1107 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1108 1109 PetscCall(VecResetArray(mpjac->xsub)); 1110 PetscCall(VecResetArray(mpjac->ysub)); 1111 PetscCall(VecRestoreArrayRead(x, &xarray)); 1112 PetscCall(VecRestoreArray(y, &yarray)); 1113 PetscFunctionReturn(0); 1114 } 1115 1116 static PetscErrorCode PCMatApply_BJacobi_Multiproc(PC pc, Mat X, Mat Y) 1117 { 1118 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1119 KSPConvergedReason reason; 1120 Mat sX, sY; 1121 const PetscScalar *x; 1122 PetscScalar *y; 1123 PetscInt m, N, lda, ldb; 1124 1125 PetscFunctionBegin; 1126 /* apply preconditioner on each matrix block */ 1127 PetscCall(MatGetLocalSize(X, &m, NULL)); 1128 PetscCall(MatGetSize(X, NULL, &N)); 1129 PetscCall(MatDenseGetLDA(X, &lda)); 1130 PetscCall(MatDenseGetLDA(Y, &ldb)); 1131 PetscCall(MatDenseGetArrayRead(X, &x)); 1132 PetscCall(MatDenseGetArrayWrite(Y, &y)); 1133 PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, (PetscScalar *)x, &sX)); 1134 PetscCall(MatCreateDense(PetscObjectComm((PetscObject)jac->ksp[0]), m, PETSC_DECIDE, PETSC_DECIDE, N, y, &sY)); 1135 PetscCall(MatDenseSetLDA(sX, lda)); 1136 PetscCall(MatDenseSetLDA(sY, ldb)); 1137 PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0)); 1138 PetscCall(KSPMatSolve(jac->ksp[0], sX, sY)); 1139 PetscCall(KSPCheckSolve(jac->ksp[0], pc, NULL)); 1140 PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, jac->ksp[0], X, Y, 0)); 1141 PetscCall(MatDestroy(&sY)); 1142 PetscCall(MatDestroy(&sX)); 1143 PetscCall(MatDenseRestoreArrayWrite(Y, &y)); 1144 PetscCall(MatDenseRestoreArrayRead(X, &x)); 1145 PetscCall(KSPGetConvergedReason(jac->ksp[0], &reason)); 1146 if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 1147 PetscFunctionReturn(0); 1148 } 1149 1150 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc) 1151 { 1152 PC_BJacobi *jac = (PC_BJacobi *)pc->data; 1153 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc *)jac->data; 1154 PetscInt m, n; 1155 MPI_Comm comm, subcomm = 0; 1156 const char *prefix; 1157 PetscBool wasSetup = PETSC_TRUE; 1158 VecType vectype; 1159 1160 PetscFunctionBegin; 1161 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 1162 PetscCheck(jac->n_local <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Only a single block in a subcommunicator is supported"); 1163 jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */ 1164 if (!pc->setupcalled) { 1165 wasSetup = PETSC_FALSE; 1166 PetscCall(PetscNew(&mpjac)); 1167 jac->data = (void *)mpjac; 1168 1169 /* initialize datastructure mpjac */ 1170 if (!jac->psubcomm) { 1171 /* Create default contiguous subcommunicatiors if user does not provide them */ 1172 PetscCall(PetscSubcommCreate(comm, &jac->psubcomm)); 1173 PetscCall(PetscSubcommSetNumber(jac->psubcomm, jac->n)); 1174 PetscCall(PetscSubcommSetType(jac->psubcomm, PETSC_SUBCOMM_CONTIGUOUS)); 1175 } 1176 mpjac->psubcomm = jac->psubcomm; 1177 subcomm = PetscSubcommChild(mpjac->psubcomm); 1178 1179 /* Get matrix blocks of pmat */ 1180 PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats)); 1181 1182 /* create a new PC that processors in each subcomm have copy of */ 1183 PetscCall(PetscMalloc1(1, &jac->ksp)); 1184 PetscCall(KSPCreate(subcomm, &jac->ksp[0])); 1185 PetscCall(KSPSetErrorIfNotConverged(jac->ksp[0], pc->erroriffailure)); 1186 PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp[0], (PetscObject)pc, 1)); 1187 PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats)); 1188 PetscCall(KSPGetPC(jac->ksp[0], &mpjac->pc)); 1189 1190 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1191 PetscCall(KSPSetOptionsPrefix(jac->ksp[0], prefix)); 1192 PetscCall(KSPAppendOptionsPrefix(jac->ksp[0], "sub_")); 1193 PetscCall(KSPGetOptionsPrefix(jac->ksp[0], &prefix)); 1194 PetscCall(MatSetOptionsPrefix(mpjac->submats, prefix)); 1195 1196 /* create dummy vectors xsub and ysub */ 1197 PetscCall(MatGetLocalSize(mpjac->submats, &m, &n)); 1198 PetscCall(VecCreateMPIWithArray(subcomm, 1, n, PETSC_DECIDE, NULL, &mpjac->xsub)); 1199 PetscCall(VecCreateMPIWithArray(subcomm, 1, m, PETSC_DECIDE, NULL, &mpjac->ysub)); 1200 PetscCall(MatGetVecType(mpjac->submats, &vectype)); 1201 PetscCall(VecSetType(mpjac->xsub, vectype)); 1202 PetscCall(VecSetType(mpjac->ysub, vectype)); 1203 1204 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiproc; 1205 pc->ops->reset = PCReset_BJacobi_Multiproc; 1206 pc->ops->destroy = PCDestroy_BJacobi_Multiproc; 1207 pc->ops->apply = PCApply_BJacobi_Multiproc; 1208 pc->ops->matapply = PCMatApply_BJacobi_Multiproc; 1209 } else { /* pc->setupcalled */ 1210 subcomm = PetscSubcommChild(mpjac->psubcomm); 1211 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1212 /* destroy old matrix blocks, then get new matrix blocks */ 1213 if (mpjac->submats) PetscCall(MatDestroy(&mpjac->submats)); 1214 PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_INITIAL_MATRIX, &mpjac->submats)); 1215 } else { 1216 PetscCall(MatGetMultiProcBlock(pc->pmat, subcomm, MAT_REUSE_MATRIX, &mpjac->submats)); 1217 } 1218 PetscCall(KSPSetOperators(jac->ksp[0], mpjac->submats, mpjac->submats)); 1219 } 1220 1221 if (!wasSetup && pc->setfromoptionscalled) PetscCall(KSPSetFromOptions(jac->ksp[0])); 1222 PetscFunctionReturn(0); 1223 } 1224