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