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