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