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