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