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