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