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