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