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