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