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 Developer Notes: This preconditioner does not currently work with CUDA/CUSP for a couple of reasons. 608 (1) It creates seq vectors as work vectors that should be cusp 609 (2) The use of VecPlaceArray() is not handled properly by CUSP (that is it will not know where 610 the ownership of the vector is so may use wrong values) even if it did know the ownership 611 it may induce extra copy ups and downs. Satish suggests a VecTransplantArray() to handle two 612 vectors sharing the same pointer and handling the CUSP side as well instead of VecGetArray()/VecPlaceArray(). 613 614 615 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 616 PCASM, PCBJacobiSetUseTrueLocal(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(), 617 PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices() 618 M*/ 619 620 EXTERN_C_BEGIN 621 #undef __FUNCT__ 622 #define __FUNCT__ "PCCreate_BJacobi" 623 PetscErrorCode PCCreate_BJacobi(PC pc) 624 { 625 PetscErrorCode ierr; 626 PetscMPIInt rank; 627 PC_BJacobi *jac; 628 629 PetscFunctionBegin; 630 ierr = PetscNewLog(pc,PC_BJacobi,&jac);CHKERRQ(ierr); 631 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 632 pc->ops->apply = 0; 633 pc->ops->applytranspose = 0; 634 pc->ops->setup = PCSetUp_BJacobi; 635 pc->ops->destroy = PCDestroy_BJacobi; 636 pc->ops->setfromoptions = PCSetFromOptions_BJacobi; 637 pc->ops->view = PCView_BJacobi; 638 pc->ops->applyrichardson = 0; 639 640 pc->data = (void*)jac; 641 jac->n = -1; 642 jac->n_local = -1; 643 jac->first_local = rank; 644 jac->ksp = 0; 645 jac->use_true_local = PETSC_FALSE; 646 jac->same_local_solves = PETSC_TRUE; 647 jac->g_lens = 0; 648 jac->l_lens = 0; 649 jac->psubcomm = 0; 650 651 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C", 652 "PCBJacobiSetUseTrueLocal_BJacobi", 653 PCBJacobiSetUseTrueLocal_BJacobi);CHKERRQ(ierr); 654 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetSubKSP_C","PCBJacobiGetSubKSP_BJacobi", 655 PCBJacobiGetSubKSP_BJacobi);CHKERRQ(ierr); 656 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetTotalBlocks_C","PCBJacobiSetTotalBlocks_BJacobi", 657 PCBJacobiSetTotalBlocks_BJacobi);CHKERRQ(ierr); 658 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetTotalBlocks_C","PCBJacobiGetTotalBlocks_BJacobi", 659 PCBJacobiGetTotalBlocks_BJacobi);CHKERRQ(ierr); 660 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetLocalBlocks_C","PCBJacobiSetLocalBlocks_BJacobi", 661 PCBJacobiSetLocalBlocks_BJacobi);CHKERRQ(ierr); 662 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetLocalBlocks_C","PCBJacobiGetLocalBlocks_BJacobi", 663 PCBJacobiGetLocalBlocks_BJacobi);CHKERRQ(ierr); 664 665 PetscFunctionReturn(0); 666 } 667 EXTERN_C_END 668 669 /* --------------------------------------------------------------------------------------------*/ 670 /* 671 These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI 672 */ 673 #undef __FUNCT__ 674 #define __FUNCT__ "PCReset_BJacobi_Singleblock" 675 PetscErrorCode PCReset_BJacobi_Singleblock(PC pc) 676 { 677 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 678 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 679 PetscErrorCode ierr; 680 681 PetscFunctionBegin; 682 ierr = KSPReset(jac->ksp[0]);CHKERRQ(ierr); 683 ierr = VecDestroy(&bjac->x);CHKERRQ(ierr); 684 ierr = VecDestroy(&bjac->y);CHKERRQ(ierr); 685 PetscFunctionReturn(0); 686 } 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "PCDestroy_BJacobi_Singleblock" 690 PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc) 691 { 692 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 693 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 694 PetscErrorCode ierr; 695 696 PetscFunctionBegin; 697 ierr = PCReset_BJacobi_Singleblock(pc);CHKERRQ(ierr); 698 ierr = KSPDestroy(&jac->ksp[0]);CHKERRQ(ierr); 699 ierr = PetscFree(jac->ksp);CHKERRQ(ierr); 700 ierr = PetscFree(jac->l_lens);CHKERRQ(ierr); 701 ierr = PetscFree(jac->g_lens);CHKERRQ(ierr); 702 ierr = PetscFree(bjac);CHKERRQ(ierr); 703 ierr = PetscFree(pc->data);CHKERRQ(ierr); 704 PetscFunctionReturn(0); 705 } 706 707 #undef __FUNCT__ 708 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Singleblock" 709 PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc) 710 { 711 PetscErrorCode ierr; 712 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 713 714 PetscFunctionBegin; 715 ierr = KSPSetUp(jac->ksp[0]);CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 #undef __FUNCT__ 720 #define __FUNCT__ "PCApply_BJacobi_Singleblock" 721 PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y) 722 { 723 PetscErrorCode ierr; 724 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 725 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 726 PetscScalar *x_array,*y_array; 727 728 PetscFunctionBegin; 729 /* 730 The VecPlaceArray() is to avoid having to copy the 731 y vector into the bjac->x vector. The reason for 732 the bjac->x vector is that we need a sequential vector 733 for the sequential solve. 734 */ 735 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 736 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 737 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 738 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 739 ierr = KSPSolve(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr); 740 ierr = VecResetArray(bjac->x);CHKERRQ(ierr); 741 ierr = VecResetArray(bjac->y);CHKERRQ(ierr); 742 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 743 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "PCApplySymmetricLeft_BJacobi_Singleblock" 749 PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y) 750 { 751 PetscErrorCode ierr; 752 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 753 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 754 PetscScalar *x_array,*y_array; 755 PC subpc; 756 757 PetscFunctionBegin; 758 /* 759 The VecPlaceArray() is to avoid having to copy the 760 y vector into the bjac->x vector. The reason for 761 the bjac->x vector is that we need a sequential vector 762 for the sequential solve. 763 */ 764 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 765 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 766 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 767 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 768 769 /* apply the symmetric left portion of the inner PC operator */ 770 /* note this by-passes the inner KSP and its options completely */ 771 772 ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr); 773 ierr = PCApplySymmetricLeft(subpc,bjac->x,bjac->y);CHKERRQ(ierr); 774 ierr = VecResetArray(bjac->x);CHKERRQ(ierr); 775 ierr = VecResetArray(bjac->y);CHKERRQ(ierr); 776 777 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 778 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 779 PetscFunctionReturn(0); 780 } 781 782 #undef __FUNCT__ 783 #define __FUNCT__ "PCApplySymmetricRight_BJacobi_Singleblock" 784 PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y) 785 { 786 PetscErrorCode ierr; 787 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 788 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 789 PetscScalar *x_array,*y_array; 790 PC subpc; 791 792 PetscFunctionBegin; 793 /* 794 The VecPlaceArray() is to avoid having to copy the 795 y vector into the bjac->x vector. The reason for 796 the bjac->x vector is that we need a sequential vector 797 for the sequential solve. 798 */ 799 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 800 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 801 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 802 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 803 804 /* apply the symmetric right portion of the inner PC operator */ 805 /* note this by-passes the inner KSP and its options completely */ 806 807 ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr); 808 ierr = PCApplySymmetricRight(subpc,bjac->x,bjac->y);CHKERRQ(ierr); 809 810 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 811 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 812 PetscFunctionReturn(0); 813 } 814 815 #undef __FUNCT__ 816 #define __FUNCT__ "PCApplyTranspose_BJacobi_Singleblock" 817 PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y) 818 { 819 PetscErrorCode ierr; 820 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 821 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 822 PetscScalar *x_array,*y_array; 823 824 PetscFunctionBegin; 825 /* 826 The VecPlaceArray() is to avoid having to copy the 827 y vector into the bjac->x vector. The reason for 828 the bjac->x vector is that we need a sequential vector 829 for the sequential solve. 830 */ 831 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 832 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 833 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 834 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 835 ierr = KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr); 836 ierr = VecResetArray(bjac->x);CHKERRQ(ierr); 837 ierr = VecResetArray(bjac->y);CHKERRQ(ierr); 838 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 839 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842 843 #undef __FUNCT__ 844 #define __FUNCT__ "PCSetUp_BJacobi_Singleblock" 845 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat) 846 { 847 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 848 PetscErrorCode ierr; 849 PetscInt m; 850 KSP ksp; 851 PC_BJacobi_Singleblock *bjac; 852 PetscBool wasSetup; 853 854 PetscFunctionBegin; 855 856 if (!pc->setupcalled) { 857 const char *prefix; 858 859 if (!jac->ksp) { 860 wasSetup = PETSC_FALSE; 861 ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 862 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 863 ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 864 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 865 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 866 ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 867 ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 868 869 pc->ops->reset = PCReset_BJacobi_Singleblock; 870 pc->ops->destroy = PCDestroy_BJacobi_Singleblock; 871 pc->ops->apply = PCApply_BJacobi_Singleblock; 872 pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock; 873 pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock; 874 pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock; 875 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock; 876 877 ierr = PetscMalloc(sizeof(KSP),&jac->ksp);CHKERRQ(ierr); 878 jac->ksp[0] = ksp; 879 880 ierr = PetscNewLog(pc,PC_BJacobi_Singleblock,&bjac);CHKERRQ(ierr); 881 jac->data = (void*)bjac; 882 } else { 883 ksp = jac->ksp[0]; 884 bjac = (PC_BJacobi_Singleblock *)jac->data; 885 } 886 887 /* 888 The reason we need to generate these vectors is to serve 889 as the right-hand side and solution vector for the solve on the 890 block. We do not need to allocate space for the vectors since 891 that is provided via VecPlaceArray() just before the call to 892 KSPSolve() on the block. 893 */ 894 ierr = MatGetSize(pmat,&m,&m);CHKERRQ(ierr); 895 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&bjac->x);CHKERRQ(ierr); 896 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&bjac->y);CHKERRQ(ierr); 897 ierr = PetscLogObjectParent(pc,bjac->x);CHKERRQ(ierr); 898 ierr = PetscLogObjectParent(pc,bjac->y);CHKERRQ(ierr); 899 } else { 900 wasSetup = PETSC_TRUE; 901 ksp = jac->ksp[0]; 902 bjac = (PC_BJacobi_Singleblock *)jac->data; 903 } 904 if (jac->use_true_local) { 905 ierr = KSPSetOperators(ksp,mat,pmat,pc->flag);CHKERRQ(ierr); 906 } else { 907 ierr = KSPSetOperators(ksp,pmat,pmat,pc->flag);CHKERRQ(ierr); 908 } 909 if (!wasSetup && pc->setfromoptionscalled) { 910 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 911 } 912 PetscFunctionReturn(0); 913 } 914 915 /* ---------------------------------------------------------------------------------------------*/ 916 #undef __FUNCT__ 917 #define __FUNCT__ "PCReset_BJacobi_Multiblock" 918 PetscErrorCode PCReset_BJacobi_Multiblock(PC pc) 919 { 920 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 921 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 922 PetscErrorCode ierr; 923 PetscInt i; 924 925 PetscFunctionBegin; 926 if (bjac && bjac->pmat) { 927 ierr = MatDestroyMatrices(jac->n_local,&bjac->pmat);CHKERRQ(ierr); 928 if (jac->use_true_local) { 929 ierr = MatDestroyMatrices(jac->n_local,&bjac->mat);CHKERRQ(ierr); 930 } 931 } 932 933 for (i=0; i<jac->n_local; i++) { 934 ierr = KSPReset(jac->ksp[i]);CHKERRQ(ierr); 935 if (bjac && bjac->x) { 936 ierr = VecDestroy(&bjac->x[i]);CHKERRQ(ierr); 937 ierr = VecDestroy(&bjac->y[i]);CHKERRQ(ierr); 938 ierr = ISDestroy(&bjac->is[i]);CHKERRQ(ierr); 939 } 940 } 941 ierr = PetscFree(jac->l_lens);CHKERRQ(ierr); 942 ierr = PetscFree(jac->g_lens);CHKERRQ(ierr); 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNCT__ 947 #define __FUNCT__ "PCDestroy_BJacobi_Multiblock" 948 PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc) 949 { 950 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 951 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 952 PetscErrorCode ierr; 953 PetscInt i; 954 955 PetscFunctionBegin; 956 ierr = PCReset_BJacobi_Multiblock(pc);CHKERRQ(ierr); 957 if (bjac) { 958 ierr = PetscFree2(bjac->x,bjac->y);CHKERRQ(ierr); 959 ierr = PetscFree(bjac->starts);CHKERRQ(ierr); 960 ierr = PetscFree(bjac->is);CHKERRQ(ierr); 961 } 962 ierr = PetscFree(jac->data);CHKERRQ(ierr); 963 for (i=0; i<jac->n_local; i++) { 964 ierr = KSPDestroy(&jac->ksp[i]);CHKERRQ(ierr); 965 } 966 ierr = PetscFree(jac->ksp);CHKERRQ(ierr); 967 ierr = PetscFree(pc->data);CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 971 #undef __FUNCT__ 972 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Multiblock" 973 PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc) 974 { 975 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 976 PetscErrorCode ierr; 977 PetscInt i,n_local = jac->n_local; 978 979 PetscFunctionBegin; 980 for (i=0; i<n_local; i++) { 981 ierr = KSPSetUp(jac->ksp[i]);CHKERRQ(ierr); 982 } 983 PetscFunctionReturn(0); 984 } 985 986 /* 987 Preconditioner for block Jacobi 988 */ 989 #undef __FUNCT__ 990 #define __FUNCT__ "PCApply_BJacobi_Multiblock" 991 PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y) 992 { 993 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 994 PetscErrorCode ierr; 995 PetscInt i,n_local = jac->n_local; 996 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 997 PetscScalar *xin,*yin; 998 999 PetscFunctionBegin; 1000 ierr = VecGetArray(x,&xin);CHKERRQ(ierr); 1001 ierr = VecGetArray(y,&yin);CHKERRQ(ierr); 1002 for (i=0; i<n_local; i++) { 1003 /* 1004 To avoid copying the subvector from x into a workspace we instead 1005 make the workspace vector array point to the subpart of the array of 1006 the global vector. 1007 */ 1008 ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr); 1009 ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr); 1010 1011 ierr = PetscLogEventBegin(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1012 ierr = KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr); 1013 ierr = PetscLogEventEnd(PC_SetUpOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1014 1015 ierr = VecResetArray(bjac->x[i]);CHKERRQ(ierr); 1016 ierr = VecResetArray(bjac->y[i]);CHKERRQ(ierr); 1017 } 1018 ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr); 1019 ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr); 1020 PetscFunctionReturn(0); 1021 } 1022 1023 /* 1024 Preconditioner for block Jacobi 1025 */ 1026 #undef __FUNCT__ 1027 #define __FUNCT__ "PCApplyTranspose_BJacobi_Multiblock" 1028 PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y) 1029 { 1030 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1031 PetscErrorCode ierr; 1032 PetscInt i,n_local = jac->n_local; 1033 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 1034 PetscScalar *xin,*yin; 1035 1036 PetscFunctionBegin; 1037 ierr = VecGetArray(x,&xin);CHKERRQ(ierr); 1038 ierr = VecGetArray(y,&yin);CHKERRQ(ierr); 1039 for (i=0; i<n_local; i++) { 1040 /* 1041 To avoid copying the subvector from x into a workspace we instead 1042 make the workspace vector array point to the subpart of the array of 1043 the global vector. 1044 */ 1045 ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr); 1046 ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr); 1047 1048 ierr = PetscLogEventBegin(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1049 ierr = KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr); 1050 ierr = PetscLogEventEnd(PC_ApplyTransposeOnBlocks,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1051 1052 ierr = VecResetArray(bjac->x[i]);CHKERRQ(ierr); 1053 ierr = VecResetArray(bjac->y[i]);CHKERRQ(ierr); 1054 } 1055 ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr); 1056 ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr); 1057 PetscFunctionReturn(0); 1058 } 1059 1060 #undef __FUNCT__ 1061 #define __FUNCT__ "PCSetUp_BJacobi_Multiblock" 1062 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat) 1063 { 1064 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1065 PetscErrorCode ierr; 1066 PetscInt m,n_local,N,M,start,i; 1067 const char *prefix,*pprefix,*mprefix; 1068 KSP ksp; 1069 Vec x,y; 1070 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 1071 PC subpc; 1072 IS is; 1073 MatReuse scall; 1074 1075 PetscFunctionBegin; 1076 ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr); 1077 1078 n_local = jac->n_local; 1079 1080 if (jac->use_true_local) { 1081 PetscBool same; 1082 ierr = PetscTypeCompare((PetscObject)mat,((PetscObject)pmat)->type_name,&same);CHKERRQ(ierr); 1083 if (!same) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Matrices not of same type"); 1084 } 1085 1086 if (!pc->setupcalled) { 1087 scall = MAT_INITIAL_MATRIX; 1088 1089 if (!jac->ksp) { 1090 pc->ops->reset = PCReset_BJacobi_Multiblock; 1091 pc->ops->destroy = PCDestroy_BJacobi_Multiblock; 1092 pc->ops->apply = PCApply_BJacobi_Multiblock; 1093 pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock; 1094 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock; 1095 1096 ierr = PetscNewLog(pc,PC_BJacobi_Multiblock,&bjac);CHKERRQ(ierr); 1097 ierr = PetscMalloc(n_local*sizeof(KSP),&jac->ksp);CHKERRQ(ierr); 1098 ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP)));CHKERRQ(ierr); 1099 ierr = PetscMalloc2(n_local,Vec,&bjac->x,n_local,Vec,&bjac->y);CHKERRQ(ierr); 1100 ierr = PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);CHKERRQ(ierr); 1101 ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar)));CHKERRQ(ierr); 1102 1103 jac->data = (void*)bjac; 1104 ierr = PetscMalloc(n_local*sizeof(IS),&bjac->is);CHKERRQ(ierr); 1105 ierr = PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS)));CHKERRQ(ierr); 1106 1107 for (i=0; i<n_local; i++) { 1108 ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 1109 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1110 ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 1111 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 1112 ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 1113 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1114 ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 1115 ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 1116 jac->ksp[i] = ksp; 1117 } 1118 } else { 1119 bjac = (PC_BJacobi_Multiblock*)jac->data; 1120 } 1121 1122 start = 0; 1123 for (i=0; i<n_local; i++) { 1124 m = jac->l_lens[i]; 1125 /* 1126 The reason we need to generate these vectors is to serve 1127 as the right-hand side and solution vector for the solve on the 1128 block. We do not need to allocate space for the vectors since 1129 that is provided via VecPlaceArray() just before the call to 1130 KSPSolve() on the block. 1131 1132 */ 1133 ierr = VecCreateSeq(PETSC_COMM_SELF,m,&x);CHKERRQ(ierr); 1134 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);CHKERRQ(ierr); 1135 ierr = PetscLogObjectParent(pc,x);CHKERRQ(ierr); 1136 ierr = PetscLogObjectParent(pc,y);CHKERRQ(ierr); 1137 bjac->x[i] = x; 1138 bjac->y[i] = y; 1139 bjac->starts[i] = start; 1140 1141 ierr = ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);CHKERRQ(ierr); 1142 bjac->is[i] = is; 1143 ierr = PetscLogObjectParent(pc,is);CHKERRQ(ierr); 1144 1145 start += m; 1146 } 1147 } else { 1148 bjac = (PC_BJacobi_Multiblock*)jac->data; 1149 /* 1150 Destroy the blocks from the previous iteration 1151 */ 1152 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1153 ierr = MatDestroyMatrices(n_local,&bjac->pmat);CHKERRQ(ierr); 1154 if (jac->use_true_local) { 1155 ierr = MatDestroyMatrices(n_local,&bjac->mat);CHKERRQ(ierr); 1156 } 1157 scall = MAT_INITIAL_MATRIX; 1158 } else { 1159 scall = MAT_REUSE_MATRIX; 1160 } 1161 } 1162 1163 ierr = MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);CHKERRQ(ierr); 1164 if (jac->use_true_local) { 1165 ierr = PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);CHKERRQ(ierr); 1166 ierr = MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);CHKERRQ(ierr); 1167 } 1168 /* Return control to the user so that the submatrices can be modified (e.g., to apply 1169 different boundary conditions for the submatrices than for the global problem) */ 1170 ierr = PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 1171 1172 ierr = PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);CHKERRQ(ierr); 1173 for (i=0; i<n_local; i++) { 1174 ierr = PetscLogObjectParent(pc,bjac->pmat[i]);CHKERRQ(ierr); 1175 ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);CHKERRQ(ierr); 1176 if (jac->use_true_local) { 1177 ierr = PetscLogObjectParent(pc,bjac->mat[i]);CHKERRQ(ierr); 1178 ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);CHKERRQ(ierr); 1179 ierr = KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr); 1180 } else { 1181 ierr = KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr); 1182 } 1183 if(pc->setfromoptionscalled){ 1184 ierr = KSPSetFromOptions(jac->ksp[i]);CHKERRQ(ierr); 1185 } 1186 } 1187 PetscFunctionReturn(0); 1188 } 1189 1190 /* ---------------------------------------------------------------------------------------------*/ 1191 /* 1192 These are for a single block with multiple processes; 1193 */ 1194 #undef __FUNCT__ 1195 #define __FUNCT__ "PCReset_BJacobi_Multiproc" 1196 static PetscErrorCode PCReset_BJacobi_Multiproc(PC pc) 1197 { 1198 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1199 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data; 1200 PetscErrorCode ierr; 1201 1202 PetscFunctionBegin; 1203 ierr = VecDestroy(&mpjac->ysub);CHKERRQ(ierr); 1204 ierr = VecDestroy(&mpjac->xsub);CHKERRQ(ierr); 1205 ierr = MatDestroy(&mpjac->submats);CHKERRQ(ierr); 1206 if (mpjac->ksp){ierr = KSPReset(mpjac->ksp);CHKERRQ(ierr);} 1207 PetscFunctionReturn(0); 1208 } 1209 1210 #undef __FUNCT__ 1211 #define __FUNCT__ "PCDestroy_BJacobi_Multiproc" 1212 static PetscErrorCode PCDestroy_BJacobi_Multiproc(PC pc) 1213 { 1214 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1215 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data; 1216 PetscErrorCode ierr; 1217 1218 PetscFunctionBegin; 1219 ierr = PCReset_BJacobi_Multiproc(pc);CHKERRQ(ierr); 1220 ierr = KSPDestroy(&mpjac->ksp);CHKERRQ(ierr); 1221 ierr = PetscSubcommDestroy(&mpjac->psubcomm);CHKERRQ(ierr); 1222 1223 ierr = PetscFree(mpjac);CHKERRQ(ierr); 1224 ierr = PetscFree(pc->data);CHKERRQ(ierr); 1225 PetscFunctionReturn(0); 1226 } 1227 1228 #undef __FUNCT__ 1229 #define __FUNCT__ "PCApply_BJacobi_Multiproc" 1230 static PetscErrorCode PCApply_BJacobi_Multiproc(PC pc,Vec x,Vec y) 1231 { 1232 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1233 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data; 1234 PetscErrorCode ierr; 1235 PetscScalar *xarray,*yarray; 1236 1237 PetscFunctionBegin; 1238 /* place x's and y's local arrays into xsub and ysub */ 1239 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 1240 ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 1241 ierr = VecPlaceArray(mpjac->xsub,xarray);CHKERRQ(ierr); 1242 ierr = VecPlaceArray(mpjac->ysub,yarray);CHKERRQ(ierr); 1243 1244 /* apply preconditioner on each matrix block */ 1245 ierr = KSPSolve(mpjac->ksp,mpjac->xsub,mpjac->ysub);CHKERRQ(ierr); 1246 1247 ierr = VecResetArray(mpjac->xsub);CHKERRQ(ierr); 1248 ierr = VecResetArray(mpjac->ysub);CHKERRQ(ierr); 1249 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 1250 ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 1251 PetscFunctionReturn(0); 1252 } 1253 1254 extern PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*); 1255 #undef __FUNCT__ 1256 #define __FUNCT__ "PCSetUp_BJacobi_Multiproc" 1257 static PetscErrorCode PCSetUp_BJacobi_Multiproc(PC pc) 1258 { 1259 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1260 PC_BJacobi_Multiproc *mpjac = (PC_BJacobi_Multiproc*)jac->data; 1261 PetscErrorCode ierr; 1262 PetscInt m,n; 1263 MPI_Comm comm = ((PetscObject)pc)->comm,subcomm=0; 1264 const char *prefix; 1265 1266 PetscFunctionBegin; 1267 if (jac->n_local > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Only a single block in a subcommunicator is supported"); 1268 jac->n_local = 1; /* currently only a single block is supported for a subcommunicator */ 1269 if (!pc->setupcalled) { 1270 ierr = PetscNewLog(pc,PC_BJacobi_Multiproc,&mpjac);CHKERRQ(ierr); 1271 jac->data = (void*)mpjac; 1272 1273 /* initialize datastructure mpjac */ 1274 if (!jac->psubcomm) { 1275 /* Create default contiguous subcommunicatiors if user does not provide them */ 1276 ierr = PetscSubcommCreate(comm,&jac->psubcomm);CHKERRQ(ierr); 1277 ierr = PetscSubcommSetNumber(jac->psubcomm,jac->n);CHKERRQ(ierr); 1278 ierr = PetscSubcommSetType(jac->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr); 1279 ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 1280 } 1281 mpjac->psubcomm = jac->psubcomm; 1282 subcomm = mpjac->psubcomm->comm; 1283 1284 /* Get matrix blocks of pmat */ 1285 ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);CHKERRQ(ierr); 1286 1287 /* create a new PC that processors in each subcomm have copy of */ 1288 ierr = KSPCreate(subcomm,&mpjac->ksp);CHKERRQ(ierr); 1289 ierr = PetscObjectIncrementTabLevel((PetscObject)mpjac->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 1290 ierr = PetscLogObjectParent(pc,mpjac->ksp);CHKERRQ(ierr); 1291 ierr = KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,pc->flag);CHKERRQ(ierr); 1292 ierr = KSPGetPC(mpjac->ksp,&mpjac->pc);CHKERRQ(ierr); 1293 1294 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1295 ierr = KSPSetOptionsPrefix(mpjac->ksp,prefix);CHKERRQ(ierr); 1296 ierr = KSPAppendOptionsPrefix(mpjac->ksp,"sub_");CHKERRQ(ierr); 1297 /* 1298 PetscMPIInt rank,subsize,subrank; 1299 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1300 ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr); 1301 ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr); 1302 1303 ierr = MatGetLocalSize(mpjac->submats,&m,PETSC_NULL);CHKERRQ(ierr); 1304 ierr = MatGetSize(mpjac->submats,&n,PETSC_NULL);CHKERRQ(ierr); 1305 ierr = PetscSynchronizedPrintf(comm,"[%d], sub-size %d,sub-rank %d\n",rank,subsize,subrank); 1306 ierr = PetscSynchronizedFlush(comm);CHKERRQ(ierr); 1307 */ 1308 1309 /* create dummy vectors xsub and ysub */ 1310 ierr = MatGetLocalSize(mpjac->submats,&m,&n);CHKERRQ(ierr); 1311 ierr = VecCreateMPIWithArray(subcomm,n,PETSC_DECIDE,PETSC_NULL,&mpjac->xsub);CHKERRQ(ierr); 1312 ierr = VecCreateMPIWithArray(subcomm,m,PETSC_DECIDE,PETSC_NULL,&mpjac->ysub);CHKERRQ(ierr); 1313 ierr = PetscLogObjectParent(pc,mpjac->xsub);CHKERRQ(ierr); 1314 ierr = PetscLogObjectParent(pc,mpjac->ysub);CHKERRQ(ierr); 1315 1316 pc->ops->reset = PCReset_BJacobi_Multiproc; 1317 pc->ops->destroy = PCDestroy_BJacobi_Multiproc; 1318 pc->ops->apply = PCApply_BJacobi_Multiproc; 1319 } else { /* pc->setupcalled */ 1320 subcomm = mpjac->psubcomm->comm; 1321 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1322 /* destroy old matrix blocks, then get new matrix blocks */ 1323 if (mpjac->submats){ierr = MatDestroy(&mpjac->submats);CHKERRQ(ierr);} 1324 ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,MAT_INITIAL_MATRIX,&mpjac->submats);CHKERRQ(ierr); 1325 } else { 1326 ierr = MatGetMultiProcBlock_MPIAIJ(pc->pmat,subcomm,MAT_REUSE_MATRIX,&mpjac->submats);CHKERRQ(ierr); 1327 } 1328 ierr = KSPSetOperators(mpjac->ksp,mpjac->submats,mpjac->submats,pc->flag);CHKERRQ(ierr); 1329 } 1330 1331 if (pc->setfromoptionscalled){ 1332 ierr = KSPSetFromOptions(mpjac->ksp);CHKERRQ(ierr); 1333 } 1334 PetscFunctionReturn(0); 1335 } 1336