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