1 /* 2 Defines a block Jacobi preconditioner. 3 */ 4 #include "src/mat/matimpl.h" 5 #include "src/ksp/pc/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 11 #undef __FUNCT__ 12 #define __FUNCT__ "PCSetUp_BJacobi" 13 static PetscErrorCode PCSetUp_BJacobi(PC pc) 14 { 15 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 16 Mat mat = pc->mat,pmat = pc->pmat; 17 PetscErrorCode ierr,(*f)(Mat,PetscTruth*,MatReuse,Mat*); 18 PetscInt N,M,start,i,sum,end; 19 PetscInt bs,i_start=-1,i_end=-1; 20 PetscMPIInt rank,size; 21 char *pprefix,*mprefix; 22 23 PetscFunctionBegin; 24 ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr); 25 ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 26 ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr); 27 ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 28 29 /* ---------- 30 Determines the number of blocks assigned to each processor 31 */ 32 33 /* local block count given */ 34 if (jac->n_local > 0 && jac->n < 0) { 35 ierr = MPI_Allreduce(&jac->n_local,&jac->n,1,MPIU_INT,MPI_SUM,pc->comm);CHKERRQ(ierr); 36 if (jac->l_lens) { /* check that user set these correctly */ 37 sum = 0; 38 for (i=0; i<jac->n_local; i++) { 39 if (jac->l_lens[i]/bs*bs !=jac->l_lens[i]) { 40 SETERRQ(PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout"); 41 } 42 sum += jac->l_lens[i]; 43 } 44 if (sum != M) SETERRQ(PETSC_ERR_ARG_SIZ,"Local lens sent incorrectly"); 45 } else { 46 ierr = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr); 47 for (i=0; i<jac->n_local; i++) { 48 jac->l_lens[i] = bs*((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i)); 49 } 50 } 51 } else if (jac->n > 0 && jac->n_local < 0) { /* global block count given */ 52 /* global blocks given: determine which ones are local */ 53 if (jac->g_lens) { 54 /* check if the g_lens is has valid entries */ 55 for (i=0; i<jac->n; i++) { 56 if (!jac->g_lens[i]) SETERRQ(PETSC_ERR_ARG_SIZ,"Zero block not allowed"); 57 if (jac->g_lens[i]/bs*bs != jac->g_lens[i]) { 58 SETERRQ(PETSC_ERR_ARG_SIZ,"Mat blocksize doesn't match block Jacobi layout"); 59 } 60 } 61 if (size == 1) { 62 jac->n_local = jac->n; 63 ierr = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr); 64 ierr = PetscMemcpy(jac->l_lens,jac->g_lens,jac->n_local*sizeof(PetscInt));CHKERRQ(ierr); 65 /* check that user set these correctly */ 66 sum = 0; 67 for (i=0; i<jac->n_local; i++) sum += jac->l_lens[i]; 68 if (sum != M) SETERRQ(PETSC_ERR_ARG_SIZ,"Global lens sent incorrectly"); 69 } else { 70 ierr = MatGetOwnershipRange(pc->pmat,&start,&end);CHKERRQ(ierr); 71 /* loop over blocks determing first one owned by me */ 72 sum = 0; 73 for (i=0; i<jac->n+1; i++) { 74 if (sum == start) { i_start = i; goto start_1;} 75 if (i < jac->n) sum += jac->g_lens[i]; 76 } 77 SETERRQ(PETSC_ERR_ARG_SIZ,"Block sizes\n\ 78 used in PCBJacobiSetTotalBlocks()\n\ 79 are 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_ERR_ARG_SIZ,"Block sizes\n\ 86 used in PCBJacobiSetTotalBlocks()\n\ 87 are not compatible with parallel matrix layout"); 88 end_1: 89 jac->n_local = i_end - i_start; 90 ierr = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr); 91 ierr = PetscMemcpy(jac->l_lens,jac->g_lens+i_start,jac->n_local*sizeof(PetscInt));CHKERRQ(ierr); 92 } 93 } else { /* no global blocks given, determine then using default layout */ 94 jac->n_local = jac->n/size + ((jac->n % size) > rank); 95 ierr = PetscMalloc(jac->n_local*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr); 96 for (i=0; i<jac->n_local; i++) { 97 jac->l_lens[i] = ((M/bs)/jac->n_local + (((M/bs) % jac->n_local) > i))*bs; 98 if (!jac->l_lens[i]) SETERRQ(PETSC_ERR_ARG_SIZ,"Too many blocks given"); 99 } 100 } 101 } else if (jac->n < 0 && jac->n_local < 0) { /* no blocks given */ 102 jac->n = size; 103 jac->n_local = 1; 104 ierr = PetscMalloc(sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr); 105 jac->l_lens[0] = M; 106 } 107 108 ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr); 109 ierr = PetscObjectQueryFunction((PetscObject)pc->mat,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 110 if (size == 1 && !f) { 111 mat = pc->mat; 112 pmat = pc->pmat; 113 } else { 114 PetscTruth iscopy; 115 MatReuse scall; 116 117 if (jac->use_true_local) { 118 scall = MAT_INITIAL_MATRIX; 119 if (pc->setupcalled) { 120 if (pc->flag == SAME_NONZERO_PATTERN) { 121 if (jac->tp_mat) { 122 scall = MAT_REUSE_MATRIX; 123 mat = jac->tp_mat; 124 } 125 } else { 126 if (jac->tp_mat) { 127 ierr = MatDestroy(jac->tp_mat);CHKERRQ(ierr); 128 } 129 } 130 } 131 if (!f) { 132 SETERRQ(PETSC_ERR_SUP,"This matrix does not support getting diagonal block"); 133 } 134 ierr = (*f)(pc->mat,&iscopy,scall,&mat);CHKERRQ(ierr); 135 /* make submatrix have same prefix as entire matrix */ 136 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->mat,&mprefix);CHKERRQ(ierr); 137 ierr = PetscObjectSetOptionsPrefix((PetscObject)mat,mprefix);CHKERRQ(ierr); 138 if (iscopy) { 139 jac->tp_mat = mat; 140 } 141 } 142 if (pc->pmat != pc->mat || !jac->use_true_local) { 143 scall = MAT_INITIAL_MATRIX; 144 if (pc->setupcalled) { 145 if (pc->flag == SAME_NONZERO_PATTERN) { 146 if (jac->tp_pmat) { 147 scall = MAT_REUSE_MATRIX; 148 pmat = jac->tp_pmat; 149 } 150 } else { 151 if (jac->tp_pmat) { 152 ierr = MatDestroy(jac->tp_pmat);CHKERRQ(ierr); 153 } 154 } 155 } 156 ierr = PetscObjectQueryFunction((PetscObject)pc->pmat,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 157 if (!f) { 158 SETERRQ(PETSC_ERR_SUP,"This matrix does not support getting diagonal block"); 159 } 160 ierr = (*f)(pc->pmat,&iscopy,scall,&pmat);CHKERRQ(ierr); 161 /* make submatrix have same prefix as entire matrix */ 162 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 163 ierr = PetscObjectSetOptionsPrefix((PetscObject)pmat,pprefix);CHKERRQ(ierr); 164 if (iscopy) { 165 jac->tp_pmat = pmat; 166 } 167 } else { 168 pmat = mat; 169 } 170 } 171 172 /* ------ 173 Setup code depends on the number of blocks 174 */ 175 if (jac->n_local == 1) { 176 ierr = PCSetUp_BJacobi_Singleblock(pc,mat,pmat);CHKERRQ(ierr); 177 } else { 178 ierr = PCSetUp_BJacobi_Multiblock(pc,mat,pmat);CHKERRQ(ierr); 179 } 180 PetscFunctionReturn(0); 181 } 182 183 /* Default destroy, if it has never been setup */ 184 #undef __FUNCT__ 185 #define __FUNCT__ "PCDestroy_BJacobi" 186 static PetscErrorCode PCDestroy_BJacobi(PC pc) 187 { 188 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 189 PetscErrorCode ierr; 190 191 PetscFunctionBegin; 192 if (jac->g_lens) {ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);} 193 if (jac->l_lens) {ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);} 194 ierr = PetscFree(jac);CHKERRQ(ierr); 195 PetscFunctionReturn(0); 196 } 197 198 #undef __FUNCT__ 199 #define __FUNCT__ "PCSetFromOptions_BJacobi" 200 static PetscErrorCode PCSetFromOptions_BJacobi(PC pc) 201 { 202 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 203 PetscErrorCode ierr; 204 PetscInt blocks; 205 PetscTruth flg; 206 207 PetscFunctionBegin; 208 ierr = PetscOptionsHead("Block Jacobi options");CHKERRQ(ierr); 209 ierr = PetscOptionsInt("-pc_bjacobi_blocks","Total number of blocks","PCBJacobiSetTotalBlocks",jac->n,&blocks,&flg);CHKERRQ(ierr); 210 if (flg) { 211 ierr = PCBJacobiSetTotalBlocks(pc,blocks,PETSC_NULL);CHKERRQ(ierr); 212 } 213 ierr = PetscOptionsName("-pc_bjacobi_truelocal","Use the true matrix, not preconditioner matrix to define matrix vector product in sub-problems","PCBJacobiSetUseTrueLocal",&flg);CHKERRQ(ierr); 214 if (flg) { 215 ierr = PCBJacobiSetUseTrueLocal(pc);CHKERRQ(ierr); 216 } 217 ierr = PetscOptionsTail();CHKERRQ(ierr); 218 PetscFunctionReturn(0); 219 } 220 221 #undef __FUNCT__ 222 #define __FUNCT__ "PCView_BJacobi" 223 static PetscErrorCode PCView_BJacobi(PC pc,PetscViewer viewer) 224 { 225 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 226 PetscErrorCode ierr; 227 PetscMPIInt rank; 228 PetscInt i; 229 PetscTruth iascii,isstring; 230 PetscViewer sviewer; 231 232 PetscFunctionBegin; 233 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 234 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 235 if (iascii) { 236 if (jac->use_true_local) { 237 ierr = PetscViewerASCIIPrintf(viewer," block Jacobi: using true local matrix, number of blocks = %D\n",jac->n);CHKERRQ(ierr); 238 } 239 ierr = PetscViewerASCIIPrintf(viewer," block Jacobi: number of blocks = %D\n",jac->n);CHKERRQ(ierr); 240 ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr); 241 if (jac->same_local_solves) { 242 ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr); 243 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 244 if (!rank && jac->ksp) { 245 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 246 ierr = KSPView(jac->ksp[0],sviewer);CHKERRQ(ierr); 247 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 248 } 249 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 250 } else { 251 252 ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 253 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D, first local block number = %D\n", 254 rank,jac->n_local,jac->first_local);CHKERRQ(ierr); 255 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 256 for (i=0; i<jac->n_local; i++) { 257 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] local block number %D\n",rank,i);CHKERRQ(ierr); 258 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 259 ierr = KSPView(jac->ksp[i],sviewer);CHKERRQ(ierr); 260 if (i != jac->n_local-1) { 261 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 262 } 263 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 264 } 265 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 266 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 267 } 268 } else if (isstring) { 269 ierr = PetscViewerStringSPrintf(viewer," blks=%D",jac->n);CHKERRQ(ierr); 270 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 271 if (jac->ksp) {ierr = KSPView(jac->ksp[0],sviewer);CHKERRQ(ierr);} 272 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 273 } else { 274 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for block Jacobi",((PetscObject)viewer)->type_name); 275 } 276 PetscFunctionReturn(0); 277 } 278 279 /* -------------------------------------------------------------------------------------*/ 280 281 EXTERN_C_BEGIN 282 #undef __FUNCT__ 283 #define __FUNCT__ "PCBJacobiSetUseTrueLocal_BJacobi" 284 PetscErrorCode PCBJacobiSetUseTrueLocal_BJacobi(PC pc) 285 { 286 PC_BJacobi *jac; 287 288 PetscFunctionBegin; 289 jac = (PC_BJacobi*)pc->data; 290 jac->use_true_local = PETSC_TRUE; 291 PetscFunctionReturn(0); 292 } 293 EXTERN_C_END 294 295 EXTERN_C_BEGIN 296 #undef __FUNCT__ 297 #define __FUNCT__ "PCBJacobiGetSubKSP_BJacobi" 298 PetscErrorCode PCBJacobiGetSubKSP_BJacobi(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 299 { 300 PC_BJacobi *jac = (PC_BJacobi*)pc->data;; 301 302 PetscFunctionBegin; 303 if (!pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call KSPSetUp() or PCSetUp() first"); 304 305 if (n_local) *n_local = jac->n_local; 306 if (first_local) *first_local = jac->first_local; 307 *ksp = jac->ksp; 308 jac->same_local_solves = PETSC_FALSE; /* Assume that local solves are now different; 309 not necessarily true though! This flag is 310 used only for PCView_BJacobi() */ 311 PetscFunctionReturn(0); 312 } 313 EXTERN_C_END 314 315 EXTERN_C_BEGIN 316 #undef __FUNCT__ 317 #define __FUNCT__ "PCBJacobiSetTotalBlocks_BJacobi" 318 PetscErrorCode PCBJacobiSetTotalBlocks_BJacobi(PC pc,PetscInt blocks,PetscInt *lens) 319 { 320 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 321 PetscErrorCode ierr; 322 323 PetscFunctionBegin; 324 325 if (pc->setupcalled > 0) SETERRQ(PETSC_ERR_ORDER,"Cannot set number of blocks after PCSetUp()/KSPSetUp() has been called"); 326 jac->n = blocks; 327 if (!lens) { 328 jac->g_lens = 0; 329 } else { 330 ierr = PetscMalloc(blocks*sizeof(PetscInt),&jac->g_lens);CHKERRQ(ierr); 331 PetscLogObjectMemory(pc,blocks*sizeof(PetscInt)); 332 ierr = PetscMemcpy(jac->g_lens,lens,blocks*sizeof(PetscInt));CHKERRQ(ierr); 333 } 334 PetscFunctionReturn(0); 335 } 336 EXTERN_C_END 337 338 EXTERN_C_BEGIN 339 #undef __FUNCT__ 340 #define __FUNCT__ "PCBJacobiGetTotalBlocks_BJacobi" 341 PetscErrorCode PCBJacobiGetTotalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) 342 { 343 PC_BJacobi *jac = (PC_BJacobi*) pc->data; 344 345 PetscFunctionBegin; 346 *blocks = jac->n; 347 if (lens) *lens = jac->g_lens; 348 PetscFunctionReturn(0); 349 } 350 EXTERN_C_END 351 352 EXTERN_C_BEGIN 353 #undef __FUNCT__ 354 #define __FUNCT__ "PCBJacobiSetLocalBlocks_BJacobi" 355 PetscErrorCode PCBJacobiSetLocalBlocks_BJacobi(PC pc,PetscInt blocks,const PetscInt lens[]) 356 { 357 PC_BJacobi *jac; 358 PetscErrorCode ierr; 359 360 PetscFunctionBegin; 361 jac = (PC_BJacobi*)pc->data; 362 363 jac->n_local = blocks; 364 if (!lens) { 365 jac->l_lens = 0; 366 } else { 367 ierr = PetscMalloc(blocks*sizeof(PetscInt),&jac->l_lens);CHKERRQ(ierr); 368 PetscLogObjectMemory(pc,blocks*sizeof(PetscInt)); 369 ierr = PetscMemcpy(jac->l_lens,lens,blocks*sizeof(PetscInt));CHKERRQ(ierr); 370 } 371 PetscFunctionReturn(0); 372 } 373 EXTERN_C_END 374 375 EXTERN_C_BEGIN 376 #undef __FUNCT__ 377 #define __FUNCT__ "PCBJacobiGetLocalBlocks_BJacobi" 378 PetscErrorCode PCBJacobiGetLocalBlocks_BJacobi(PC pc, PetscInt *blocks, const PetscInt *lens[]) 379 { 380 PC_BJacobi *jac = (PC_BJacobi*) pc->data; 381 382 PetscFunctionBegin; 383 *blocks = jac->n_local; 384 if (lens) *lens = jac->l_lens; 385 PetscFunctionReturn(0); 386 } 387 EXTERN_C_END 388 389 /* -------------------------------------------------------------------------------------*/ 390 391 #undef __FUNCT__ 392 #define __FUNCT__ "PCBJacobiSetUseTrueLocal" 393 /*@ 394 PCBJacobiSetUseTrueLocal - Sets a flag to indicate that the block 395 problem is associated with the linear system matrix instead of the 396 default (where it is associated with the preconditioning matrix). 397 That is, if the local system is solved iteratively then it iterates 398 on the block from the matrix using the block from the preconditioner 399 as the preconditioner for the local block. 400 401 Collective on PC 402 403 Input Parameters: 404 . pc - the preconditioner context 405 406 Options Database Key: 407 . -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal() 408 409 Notes: 410 For the common case in which the preconditioning and linear 411 system matrices are identical, this routine is unnecessary. 412 413 Level: intermediate 414 415 .keywords: block, Jacobi, set, true, local, flag 416 417 .seealso: PCSetOperators(), PCBJacobiSetLocalBlocks() 418 @*/ 419 PetscErrorCode PCBJacobiSetUseTrueLocal(PC pc) 420 { 421 PetscErrorCode ierr,(*f)(PC); 422 423 PetscFunctionBegin; 424 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 425 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 426 if (f) { 427 ierr = (*f)(pc);CHKERRQ(ierr); 428 } 429 430 PetscFunctionReturn(0); 431 } 432 433 #undef __FUNCT__ 434 #define __FUNCT__ "PCBJacobiGetSubKSP" 435 /*@C 436 PCBJacobiGetSubKSP - Gets the local KSP contexts for all blocks on 437 this processor. 438 439 Note Collective 440 441 Input Parameter: 442 . pc - the preconditioner context 443 444 Output Parameters: 445 + n_local - the number of blocks on this processor, or PETSC_NULL 446 . first_local - the global number of the first block on this processor, or PETSC_NULL 447 - ksp - the array of KSP contexts 448 449 Notes: 450 After PCBJacobiGetSubKSP() the array of KSP contexts is not to be freed. 451 452 Currently for some matrix implementations only 1 block per processor 453 is supported. 454 455 You must call KSPSetUp() or PCSetUp() before calling PCBJacobiGetSubKSP(). 456 457 Level: advanced 458 459 .keywords: block, Jacobi, get, sub, KSP, context 460 461 .seealso: PCBJacobiGetSubKSP() 462 @*/ 463 PetscErrorCode PCBJacobiGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 464 { 465 PetscErrorCode ierr,(*f)(PC,PetscInt *,PetscInt *,KSP **); 466 467 PetscFunctionBegin; 468 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 469 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 470 if (f) { 471 ierr = (*f)(pc,n_local,first_local,ksp);CHKERRQ(ierr); 472 } else { 473 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subsolvers for this preconditioner"); 474 } 475 PetscFunctionReturn(0); 476 } 477 478 #undef __FUNCT__ 479 #define __FUNCT__ "PCBJacobiSetTotalBlocks" 480 /*@ 481 PCBJacobiSetTotalBlocks - Sets the global number of blocks for the block 482 Jacobi preconditioner. 483 484 Collective on PC 485 486 Input Parameters: 487 + pc - the preconditioner context 488 . blocks - the number of blocks 489 - lens - [optional] integer array containing the size of each block 490 491 Options Database Key: 492 . -pc_bjacobi_blocks <blocks> - Sets the number of global blocks 493 494 Notes: 495 Currently only a limited number of blocking configurations are supported. 496 All processors sharing the PC must call this routine with the same data. 497 498 Level: intermediate 499 500 .keywords: set, number, Jacobi, global, total, blocks 501 502 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetLocalBlocks() 503 @*/ 504 PetscErrorCode PCBJacobiSetTotalBlocks(PC pc,PetscInt blocks,const PetscInt lens[]) 505 { 506 PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt[]); 507 508 PetscFunctionBegin; 509 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 510 if (blocks <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have positive blocks"); 511 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetTotalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 512 if (f) { 513 ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr); 514 } 515 PetscFunctionReturn(0); 516 } 517 518 #undef __FUNCT__ 519 #define __FUNCT__ "PCBJacobiGetTotalBlocks" 520 /*@C 521 PCBJacobiGetTotalBlocks - Gets the global number of blocks for the block 522 Jacobi preconditioner. 523 524 Collective on PC 525 526 Input Parameter: 527 . pc - the preconditioner context 528 529 Output parameters: 530 + blocks - the number of blocks 531 - lens - integer array containing the size of each block 532 533 Level: intermediate 534 535 .keywords: get, number, Jacobi, global, total, blocks 536 537 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetLocalBlocks() 538 @*/ 539 PetscErrorCode PCBJacobiGetTotalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) 540 { 541 PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]); 542 543 PetscFunctionBegin; 544 PetscValidHeaderSpecific(pc, PC_COOKIE,1); 545 PetscValidIntPointer(blocks,2); 546 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetTotalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 547 if (f) { 548 ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr); 549 } 550 PetscFunctionReturn(0); 551 } 552 553 #undef __FUNCT__ 554 #define __FUNCT__ "PCBJacobiSetLocalBlocks" 555 /*@ 556 PCBJacobiSetLocalBlocks - Sets the local number of blocks for the block 557 Jacobi preconditioner. 558 559 Not Collective 560 561 Input Parameters: 562 + pc - the preconditioner context 563 . blocks - the number of blocks 564 - lens - [optional] integer array containing size of each block 565 566 Note: 567 Currently only a limited number of blocking configurations are supported. 568 569 Level: intermediate 570 571 .keywords: PC, set, number, Jacobi, local, blocks 572 573 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiSetTotalBlocks() 574 @*/ 575 PetscErrorCode PCBJacobiSetLocalBlocks(PC pc,PetscInt blocks,const PetscInt lens[]) 576 { 577 PetscErrorCode ierr,(*f)(PC,PetscInt,const PetscInt []); 578 579 PetscFunctionBegin; 580 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 581 if (blocks < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Must have nonegative blocks"); 582 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiSetLocalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 583 if (f) { 584 ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr); 585 } 586 PetscFunctionReturn(0); 587 } 588 589 #undef __FUNCT__ 590 #define __FUNCT__ "PCBJacobiGetLocalBlocks" 591 /*@C 592 PCBJacobiGetLocalBlocks - Gets the local number of blocks for the block 593 Jacobi preconditioner. 594 595 Not Collective 596 597 Input Parameters: 598 + pc - the preconditioner context 599 . blocks - the number of blocks 600 - lens - [optional] integer array containing size of each block 601 602 Note: 603 Currently only a limited number of blocking configurations are supported. 604 605 Level: intermediate 606 607 .keywords: PC, get, number, Jacobi, local, blocks 608 609 .seealso: PCBJacobiSetUseTrueLocal(), PCBJacobiGetTotalBlocks() 610 @*/ 611 PetscErrorCode PCBJacobiGetLocalBlocks(PC pc, PetscInt *blocks, const PetscInt *lens[]) 612 { 613 PetscErrorCode ierr,(*f)(PC,PetscInt*, const PetscInt *[]); 614 615 PetscFunctionBegin; 616 PetscValidHeaderSpecific(pc, PC_COOKIE,1); 617 PetscValidIntPointer(blocks,2); 618 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCBJacobiGetLocalBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 619 if (f) { 620 ierr = (*f)(pc,blocks,lens);CHKERRQ(ierr); 621 } 622 PetscFunctionReturn(0); 623 } 624 625 /* -----------------------------------------------------------------------------------*/ 626 627 /*MC 628 PCBJACOBI - Use block Jacobi preconditioning, each block is (approximately) solved with 629 its own KSP object. 630 631 Options Database Keys: 632 . -pc_bjacobi_truelocal - Activates PCBJacobiSetUseTrueLocal() 633 634 Notes: Each processor can have one or more blocks, but a block cannot be shared by more 635 than one processor. Defaults to one block per processor. 636 637 To set options on the solvers for each block append -sub_ to all the KSP, KSP, and PC 638 options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1 -sub_ksp_type preonly 639 640 To set the options on the solvers seperate for each block call PCBJacobiGetSubKSP() 641 and set the options directly on the resulting KSP object (you can access its PC 642 KSPGetPC()) 643 644 Level: beginner 645 646 Concepts: block Jacobi 647 648 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 649 PCASM, PCBJacobiSetUseTrueLocal(), PCBJacobiGetSubKSP(), PCBJacobiSetTotalBlocks(), 650 PCBJacobiSetLocalBlocks(), PCSetModifySubmatrices() 651 M*/ 652 653 EXTERN_C_BEGIN 654 #undef __FUNCT__ 655 #define __FUNCT__ "PCCreate_BJacobi" 656 PetscErrorCode PCCreate_BJacobi(PC pc) 657 { 658 PetscErrorCode ierr; 659 PetscMPIInt rank; 660 PC_BJacobi *jac; 661 662 PetscFunctionBegin; 663 ierr = PetscNew(PC_BJacobi,&jac);CHKERRQ(ierr); 664 PetscLogObjectMemory(pc,sizeof(PC_BJacobi)); 665 ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr); 666 pc->ops->apply = 0; 667 pc->ops->applytranspose = 0; 668 pc->ops->setup = PCSetUp_BJacobi; 669 pc->ops->destroy = PCDestroy_BJacobi; 670 pc->ops->setfromoptions = PCSetFromOptions_BJacobi; 671 pc->ops->view = PCView_BJacobi; 672 pc->ops->applyrichardson = 0; 673 674 pc->data = (void*)jac; 675 jac->n = -1; 676 jac->n_local = -1; 677 jac->first_local = rank; 678 jac->ksp = 0; 679 jac->use_true_local = PETSC_FALSE; 680 jac->same_local_solves = PETSC_TRUE; 681 jac->g_lens = 0; 682 jac->l_lens = 0; 683 jac->tp_mat = 0; 684 jac->tp_pmat = 0; 685 686 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetUseTrueLocal_C", 687 "PCBJacobiSetUseTrueLocal_BJacobi", 688 PCBJacobiSetUseTrueLocal_BJacobi);CHKERRQ(ierr); 689 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetSubKSP_C","PCBJacobiGetSubKSP_BJacobi", 690 PCBJacobiGetSubKSP_BJacobi);CHKERRQ(ierr); 691 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetTotalBlocks_C","PCBJacobiSetTotalBlocks_BJacobi", 692 PCBJacobiSetTotalBlocks_BJacobi);CHKERRQ(ierr); 693 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetTotalBlocks_C","PCBJacobiGetTotalBlocks_BJacobi", 694 PCBJacobiGetTotalBlocks_BJacobi);CHKERRQ(ierr); 695 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiSetLocalBlocks_C","PCBJacobiSetLocalBlocks_BJacobi", 696 PCBJacobiSetLocalBlocks_BJacobi);CHKERRQ(ierr); 697 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCBJacobiGetLocalBlocks_C","PCBJacobiGetLocalBlocks_BJacobi", 698 PCBJacobiGetLocalBlocks_BJacobi);CHKERRQ(ierr); 699 700 PetscFunctionReturn(0); 701 } 702 EXTERN_C_END 703 704 /* --------------------------------------------------------------------------------------------*/ 705 /* 706 These are for a single block per processor; works for AIJ, BAIJ; Seq and MPI 707 */ 708 #undef __FUNCT__ 709 #define __FUNCT__ "PCDestroy_BJacobi_Singleblock" 710 PetscErrorCode PCDestroy_BJacobi_Singleblock(PC pc) 711 { 712 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 713 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 714 PetscErrorCode ierr; 715 716 PetscFunctionBegin; 717 /* 718 If the on processor block had to be generated via a MatGetDiagonalBlock() 719 that creates a copy (for example MPIBDiag matrices do), this frees the space 720 */ 721 if (jac->tp_mat) { 722 ierr = MatDestroy(jac->tp_mat);CHKERRQ(ierr); 723 } 724 if (jac->tp_pmat) { 725 ierr = MatDestroy(jac->tp_pmat);CHKERRQ(ierr); 726 } 727 728 ierr = KSPDestroy(jac->ksp[0]);CHKERRQ(ierr); 729 ierr = PetscFree(jac->ksp);CHKERRQ(ierr); 730 ierr = VecDestroy(bjac->x);CHKERRQ(ierr); 731 ierr = VecDestroy(bjac->y);CHKERRQ(ierr); 732 if (jac->l_lens) {ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);} 733 if (jac->g_lens) {ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);} 734 ierr = PetscFree(bjac);CHKERRQ(ierr); 735 ierr = PetscFree(jac);CHKERRQ(ierr); 736 PetscFunctionReturn(0); 737 } 738 739 #undef __FUNCT__ 740 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Singleblock" 741 PetscErrorCode PCSetUpOnBlocks_BJacobi_Singleblock(PC pc) 742 { 743 PetscErrorCode ierr; 744 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 745 746 PetscFunctionBegin; 747 ierr = KSPSetUp(jac->ksp[0]);CHKERRQ(ierr); 748 PetscFunctionReturn(0); 749 } 750 751 #undef __FUNCT__ 752 #define __FUNCT__ "PCApply_BJacobi_Singleblock" 753 PetscErrorCode PCApply_BJacobi_Singleblock(PC pc,Vec x,Vec y) 754 { 755 PetscErrorCode ierr; 756 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 757 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 758 PetscScalar *x_array,*y_array; 759 760 PetscFunctionBegin; 761 /* 762 The VecPlaceArray() is to avoid having to copy the 763 y vector into the bjac->x vector. The reason for 764 the bjac->x vector is that we need a sequential vector 765 for the sequential solve. 766 */ 767 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 768 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 769 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 770 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 771 ierr = KSPSolve(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr); 772 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 773 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 774 PetscFunctionReturn(0); 775 } 776 777 #undef __FUNCT__ 778 #define __FUNCT__ "PCApplySymmetricLeft_BJacobi_Singleblock" 779 PetscErrorCode PCApplySymmetricLeft_BJacobi_Singleblock(PC pc,Vec x,Vec y) 780 { 781 PetscErrorCode ierr; 782 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 783 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 784 PetscScalar *x_array,*y_array; 785 PC subpc; 786 787 PetscFunctionBegin; 788 /* 789 The VecPlaceArray() is to avoid having to copy the 790 y vector into the bjac->x vector. The reason for 791 the bjac->x vector is that we need a sequential vector 792 for the sequential solve. 793 */ 794 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 795 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 796 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 797 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 798 799 /* apply the symmetric left portion of the inner PC operator */ 800 /* note this by-passes the inner KSP and its options completely */ 801 802 ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr); 803 ierr = PCApplySymmetricLeft(subpc,bjac->x,bjac->y);CHKERRQ(ierr); 804 805 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 806 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 807 PetscFunctionReturn(0); 808 } 809 810 #undef __FUNCT__ 811 #define __FUNCT__ "PCApplySymmetricRight_BJacobi_Singleblock" 812 PetscErrorCode PCApplySymmetricRight_BJacobi_Singleblock(PC pc,Vec x,Vec y) 813 { 814 PetscErrorCode ierr; 815 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 816 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 817 PetscScalar *x_array,*y_array; 818 PC subpc; 819 820 PetscFunctionBegin; 821 /* 822 The VecPlaceArray() is to avoid having to copy the 823 y vector into the bjac->x vector. The reason for 824 the bjac->x vector is that we need a sequential vector 825 for the sequential solve. 826 */ 827 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 828 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 829 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 830 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 831 832 /* apply the symmetric right portion of the inner PC operator */ 833 /* note this by-passes the inner KSP and its options completely */ 834 835 ierr = KSPGetPC(jac->ksp[0],&subpc);CHKERRQ(ierr); 836 ierr = PCApplySymmetricRight(subpc,bjac->x,bjac->y);CHKERRQ(ierr); 837 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__ "PCApplyTranspose_BJacobi_Singleblock" 845 PetscErrorCode PCApplyTranspose_BJacobi_Singleblock(PC pc,Vec x,Vec y) 846 { 847 PetscErrorCode ierr; 848 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 849 PC_BJacobi_Singleblock *bjac = (PC_BJacobi_Singleblock*)jac->data; 850 PetscScalar *x_array,*y_array; 851 852 PetscFunctionBegin; 853 /* 854 The VecPlaceArray() is to avoid having to copy the 855 y vector into the bjac->x vector. The reason for 856 the bjac->x vector is that we need a sequential vector 857 for the sequential solve. 858 */ 859 ierr = VecGetArray(x,&x_array);CHKERRQ(ierr); 860 ierr = VecGetArray(y,&y_array);CHKERRQ(ierr); 861 ierr = VecPlaceArray(bjac->x,x_array);CHKERRQ(ierr); 862 ierr = VecPlaceArray(bjac->y,y_array);CHKERRQ(ierr); 863 ierr = KSPSolveTranspose(jac->ksp[0],bjac->x,bjac->y);CHKERRQ(ierr); 864 ierr = VecRestoreArray(x,&x_array);CHKERRQ(ierr); 865 ierr = VecRestoreArray(y,&y_array);CHKERRQ(ierr); 866 PetscFunctionReturn(0); 867 } 868 869 #undef __FUNCT__ 870 #define __FUNCT__ "PCSetUp_BJacobi_Singleblock" 871 static PetscErrorCode PCSetUp_BJacobi_Singleblock(PC pc,Mat mat,Mat pmat) 872 { 873 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 874 PetscErrorCode ierr; 875 PetscInt m; 876 KSP ksp; 877 Vec x,y; 878 PC_BJacobi_Singleblock *bjac; 879 PC subpc; 880 881 PetscFunctionBegin; 882 883 /* set default direct solver with no Krylov method */ 884 if (!pc->setupcalled) { 885 char *prefix; 886 ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 887 PetscLogObjectParent(pc,ksp); 888 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 889 ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 890 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 891 ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 892 ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 893 /* 894 The reason we need to generate these vectors is to serve 895 as the right-hand side and solution vector for the solve on the 896 block. We do not need to allocate space for the vectors since 897 that is provided via VecPlaceArray() just before the call to 898 KSPSolve() on the block. 899 */ 900 ierr = MatGetSize(pmat,&m,&m);CHKERRQ(ierr); 901 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&x);CHKERRQ(ierr); 902 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);CHKERRQ(ierr); 903 PetscLogObjectParent(pc,x); 904 PetscLogObjectParent(pc,y); 905 906 pc->ops->destroy = PCDestroy_BJacobi_Singleblock; 907 pc->ops->apply = PCApply_BJacobi_Singleblock; 908 pc->ops->applysymmetricleft = PCApplySymmetricLeft_BJacobi_Singleblock; 909 pc->ops->applysymmetricright = PCApplySymmetricRight_BJacobi_Singleblock; 910 pc->ops->applytranspose = PCApplyTranspose_BJacobi_Singleblock; 911 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Singleblock; 912 913 ierr = PetscMalloc(sizeof(PC_BJacobi_Singleblock),&bjac);CHKERRQ(ierr); 914 PetscLogObjectMemory(pc,sizeof(PC_BJacobi_Singleblock)); 915 bjac->x = x; 916 bjac->y = y; 917 918 ierr = PetscMalloc(sizeof(KSP),&jac->ksp);CHKERRQ(ierr); 919 jac->ksp[0] = ksp; 920 jac->data = (void*)bjac; 921 } else { 922 ksp = jac->ksp[0]; 923 bjac = (PC_BJacobi_Singleblock *)jac->data; 924 } 925 if (jac->use_true_local) { 926 ierr = KSPSetOperators(ksp,mat,pmat,pc->flag);CHKERRQ(ierr); 927 } else { 928 ierr = KSPSetOperators(ksp,pmat,pmat,pc->flag);CHKERRQ(ierr); 929 } 930 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 931 PetscFunctionReturn(0); 932 } 933 934 /* ---------------------------------------------------------------------------------------------*/ 935 936 #undef __FUNCT__ 937 #define __FUNCT__ "PCDestroy_BJacobi_Multiblock" 938 PetscErrorCode PCDestroy_BJacobi_Multiblock(PC pc) 939 { 940 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 941 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 942 PetscErrorCode ierr; 943 PetscInt i; 944 945 PetscFunctionBegin; 946 ierr = MatDestroyMatrices(jac->n_local,&bjac->pmat);CHKERRQ(ierr); 947 if (jac->use_true_local) { 948 ierr = MatDestroyMatrices(jac->n_local,&bjac->mat);CHKERRQ(ierr); 949 } 950 951 /* 952 If the on processor block had to be generated via a MatGetDiagonalBlock() 953 that creates a copy (for example MPIBDiag matrices do), this frees the space 954 */ 955 if (jac->tp_mat) { 956 ierr = MatDestroy(jac->tp_mat);CHKERRQ(ierr); 957 } 958 if (jac->tp_pmat) { 959 ierr = MatDestroy(jac->tp_pmat);CHKERRQ(ierr); 960 } 961 962 for (i=0; i<jac->n_local; i++) { 963 ierr = KSPDestroy(jac->ksp[i]);CHKERRQ(ierr); 964 ierr = VecDestroy(bjac->x[i]);CHKERRQ(ierr); 965 ierr = VecDestroy(bjac->y[i]);CHKERRQ(ierr); 966 ierr = ISDestroy(bjac->is[i]);CHKERRQ(ierr); 967 } 968 ierr = PetscFree(jac->ksp);CHKERRQ(ierr); 969 ierr = PetscFree(bjac->x);CHKERRQ(ierr); 970 ierr = PetscFree(bjac->starts);CHKERRQ(ierr); 971 ierr = PetscFree(bjac->is);CHKERRQ(ierr); 972 ierr = PetscFree(bjac);CHKERRQ(ierr); 973 if (jac->l_lens) {ierr = PetscFree(jac->l_lens);CHKERRQ(ierr);} 974 if (jac->g_lens) {ierr = PetscFree(jac->g_lens);CHKERRQ(ierr);} 975 ierr = PetscFree(jac);CHKERRQ(ierr); 976 PetscFunctionReturn(0); 977 } 978 979 #undef __FUNCT__ 980 #define __FUNCT__ "PCSetUpOnBlocks_BJacobi_Multiblock" 981 PetscErrorCode PCSetUpOnBlocks_BJacobi_Multiblock(PC pc) 982 { 983 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 984 PetscErrorCode ierr; 985 PetscInt i,n_local = jac->n_local; 986 987 PetscFunctionBegin; 988 for (i=0; i<n_local; i++) { 989 ierr = KSPSetUp(jac->ksp[i]);CHKERRQ(ierr); 990 } 991 PetscFunctionReturn(0); 992 } 993 994 /* 995 Preconditioner for block Jacobi 996 */ 997 #undef __FUNCT__ 998 #define __FUNCT__ "PCApply_BJacobi_Multiblock" 999 PetscErrorCode PCApply_BJacobi_Multiblock(PC pc,Vec x,Vec y) 1000 { 1001 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1002 PetscErrorCode ierr; 1003 PetscInt i,n_local = jac->n_local; 1004 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 1005 PetscScalar *xin,*yin; 1006 static PetscTruth flag = PETSC_TRUE; 1007 static PetscEvent SUBKspSolve; 1008 1009 PetscFunctionBegin; 1010 if (flag) { 1011 ierr = PetscLogEventRegister(&SUBKspSolve,"SubKspSolve",KSP_COOKIE);CHKERRQ(ierr); 1012 flag = PETSC_FALSE; 1013 } 1014 ierr = VecGetArray(x,&xin);CHKERRQ(ierr); 1015 ierr = VecGetArray(y,&yin);CHKERRQ(ierr); 1016 for (i=0; i<n_local; i++) { 1017 /* 1018 To avoid copying the subvector from x into a workspace we instead 1019 make the workspace vector array point to the subpart of the array of 1020 the global vector. 1021 */ 1022 ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr); 1023 ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr); 1024 1025 ierr = PetscLogEventBegin(SUBKspSolve,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1026 ierr = KSPSolve(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr); 1027 ierr = PetscLogEventEnd(SUBKspSolve,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1028 } 1029 ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr); 1030 ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 /* 1035 Preconditioner for block Jacobi 1036 */ 1037 #undef __FUNCT__ 1038 #define __FUNCT__ "PCApplyTranspose_BJacobi_Multiblock" 1039 PetscErrorCode PCApplyTranspose_BJacobi_Multiblock(PC pc,Vec x,Vec y) 1040 { 1041 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1042 PetscErrorCode ierr; 1043 PetscInt i,n_local = jac->n_local; 1044 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 1045 PetscScalar *xin,*yin; 1046 static PetscTruth flag = PETSC_TRUE; 1047 static PetscEvent SUBKspSolve; 1048 1049 PetscFunctionBegin; 1050 if (flag) { 1051 ierr = PetscLogEventRegister(&SUBKspSolve,"SubKspSolveTranspose",KSP_COOKIE);CHKERRQ(ierr); 1052 flag = PETSC_FALSE; 1053 } 1054 ierr = VecGetArray(x,&xin);CHKERRQ(ierr); 1055 ierr = VecGetArray(y,&yin);CHKERRQ(ierr); 1056 for (i=0; i<n_local; i++) { 1057 /* 1058 To avoid copying the subvector from x into a workspace we instead 1059 make the workspace vector array point to the subpart of the array of 1060 the global vector. 1061 */ 1062 ierr = VecPlaceArray(bjac->x[i],xin+bjac->starts[i]);CHKERRQ(ierr); 1063 ierr = VecPlaceArray(bjac->y[i],yin+bjac->starts[i]);CHKERRQ(ierr); 1064 1065 ierr = PetscLogEventBegin(SUBKspSolve,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1066 ierr = KSPSolveTranspose(jac->ksp[i],bjac->x[i],bjac->y[i]);CHKERRQ(ierr); 1067 ierr = PetscLogEventEnd(SUBKspSolve,jac->ksp[i],bjac->x[i],bjac->y[i],0);CHKERRQ(ierr); 1068 } 1069 ierr = VecRestoreArray(x,&xin);CHKERRQ(ierr); 1070 ierr = VecRestoreArray(y,&yin);CHKERRQ(ierr); 1071 PetscFunctionReturn(0); 1072 } 1073 1074 #undef __FUNCT__ 1075 #define __FUNCT__ "PCSetUp_BJacobi_Multiblock" 1076 static PetscErrorCode PCSetUp_BJacobi_Multiblock(PC pc,Mat mat,Mat pmat) 1077 { 1078 PC_BJacobi *jac = (PC_BJacobi*)pc->data; 1079 PetscErrorCode ierr; 1080 PetscInt m,n_local,N,M,start,i; 1081 char *prefix,*pprefix,*mprefix; 1082 KSP ksp; 1083 Vec x,y; 1084 PC_BJacobi_Multiblock *bjac = (PC_BJacobi_Multiblock*)jac->data; 1085 PC subpc; 1086 IS is; 1087 MatReuse scall = MAT_REUSE_MATRIX; 1088 1089 PetscFunctionBegin; 1090 ierr = MatGetLocalSize(pc->pmat,&M,&N);CHKERRQ(ierr); 1091 1092 n_local = jac->n_local; 1093 1094 if (jac->use_true_local) { 1095 if (mat->type != pmat->type) SETERRQ(PETSC_ERR_ARG_INCOMP,"Matrices not of same type"); 1096 } 1097 1098 /* set default direct solver with no Krylov method */ 1099 if (!pc->setupcalled) { 1100 scall = MAT_INITIAL_MATRIX; 1101 pc->ops->destroy = PCDestroy_BJacobi_Multiblock; 1102 pc->ops->apply = PCApply_BJacobi_Multiblock; 1103 pc->ops->applytranspose= PCApplyTranspose_BJacobi_Multiblock; 1104 pc->ops->setuponblocks = PCSetUpOnBlocks_BJacobi_Multiblock; 1105 1106 ierr = PetscMalloc(sizeof(PC_BJacobi_Multiblock),&bjac);CHKERRQ(ierr); 1107 PetscLogObjectMemory(pc,sizeof(PC_BJacobi_Multiblock)); 1108 ierr = PetscMalloc(n_local*sizeof(KSP),&jac->ksp);CHKERRQ(ierr); 1109 PetscLogObjectMemory(pc,sizeof(n_local*sizeof(KSP))); 1110 ierr = PetscMalloc(2*n_local*sizeof(Vec),&bjac->x);CHKERRQ(ierr); 1111 PetscLogObjectMemory(pc,sizeof(2*n_local*sizeof(Vec))); 1112 bjac->y = bjac->x + n_local; 1113 ierr = PetscMalloc(n_local*sizeof(PetscScalar),&bjac->starts);CHKERRQ(ierr); 1114 PetscLogObjectMemory(pc,sizeof(n_local*sizeof(PetscScalar))); 1115 1116 jac->data = (void*)bjac; 1117 ierr = PetscMalloc(n_local*sizeof(IS),&bjac->is);CHKERRQ(ierr); 1118 PetscLogObjectMemory(pc,sizeof(n_local*sizeof(IS))); 1119 1120 start = 0; 1121 for (i=0; i<n_local; i++) { 1122 ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 1123 PetscLogObjectParent(pc,ksp); 1124 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 1125 ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 1126 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1127 ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 1128 ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 1129 1130 m = jac->l_lens[i]; 1131 1132 /* 1133 The reason we need to generate these vectors is to serve 1134 as the right-hand side and solution vector for the solve on the 1135 block. We do not need to allocate space for the vectors since 1136 that is provided via VecPlaceArray() just before the call to 1137 KSPSolve() on the block. 1138 1139 */ 1140 ierr = VecCreateSeq(PETSC_COMM_SELF,m,&x);CHKERRQ(ierr); 1141 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,m,PETSC_NULL,&y);CHKERRQ(ierr); 1142 PetscLogObjectParent(pc,x); 1143 PetscLogObjectParent(pc,y); 1144 bjac->x[i] = x; 1145 bjac->y[i] = y; 1146 bjac->starts[i] = start; 1147 jac->ksp[i] = ksp; 1148 1149 ierr = ISCreateStride(PETSC_COMM_SELF,m,start,1,&is);CHKERRQ(ierr); 1150 bjac->is[i] = is; 1151 PetscLogObjectParent(pc,is); 1152 1153 start += m; 1154 } 1155 } else { 1156 bjac = (PC_BJacobi_Multiblock*)jac->data; 1157 /* 1158 Destroy the blocks from the previous iteration 1159 */ 1160 if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1161 ierr = MatDestroyMatrices(n_local,&bjac->pmat);CHKERRQ(ierr); 1162 if (jac->use_true_local) { 1163 ierr = MatDestroyMatrices(n_local,&bjac->mat);CHKERRQ(ierr); 1164 } 1165 scall = MAT_INITIAL_MATRIX; 1166 } 1167 } 1168 1169 ierr = MatGetSubMatrices(pmat,n_local,bjac->is,bjac->is,scall,&bjac->pmat);CHKERRQ(ierr); 1170 if (jac->use_true_local) { 1171 ierr = PetscObjectGetOptionsPrefix((PetscObject)mat,&mprefix);CHKERRQ(ierr); 1172 ierr = MatGetSubMatrices(mat,n_local,bjac->is,bjac->is,scall,&bjac->mat);CHKERRQ(ierr); 1173 } 1174 /* Return control to the user so that the submatrices can be modified (e.g., to apply 1175 different boundary conditions for the submatrices than for the global problem) */ 1176 ierr = PCModifySubMatrices(pc,n_local,bjac->is,bjac->is,bjac->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 1177 1178 ierr = PetscObjectGetOptionsPrefix((PetscObject)pmat,&pprefix);CHKERRQ(ierr); 1179 for (i=0; i<n_local; i++) { 1180 PetscLogObjectParent(pc,bjac->pmat[i]); 1181 ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->pmat[i],pprefix);CHKERRQ(ierr); 1182 if (jac->use_true_local) { 1183 PetscLogObjectParent(pc,bjac->mat[i]); 1184 ierr = PetscObjectSetOptionsPrefix((PetscObject)bjac->mat[i],mprefix);CHKERRQ(ierr); 1185 ierr = KSPSetOperators(jac->ksp[i],bjac->mat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr); 1186 } else { 1187 ierr = KSPSetOperators(jac->ksp[i],bjac->pmat[i],bjac->pmat[i],pc->flag);CHKERRQ(ierr); 1188 } 1189 ierr = KSPSetFromOptions(jac->ksp[i]);CHKERRQ(ierr); 1190 } 1191 1192 PetscFunctionReturn(0); 1193 } 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204