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