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