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