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