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