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