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