1 #define PETSCKSP_DLL 2 3 /* 4 This file defines an additive Schwarz preconditioner for any Mat implementation. 5 6 Note that each processor may have any number of subdomains. But in order to 7 deal easily with the VecScatter(), we treat each processor as if it has the 8 same number of subdomains. 9 10 n - total number of true subdomains on all processors 11 n_local_true - actual number of subdomains on this processor 12 n_local = maximum over all processors of n_local_true 13 */ 14 #include "private/pcimpl.h" /*I "petscpc.h" I*/ 15 16 typedef struct { 17 PetscInt n,n_local,n_local_true; 18 PetscInt overlap; /* overlap requested by user */ 19 KSP *ksp; /* linear solvers for each block */ 20 VecScatter *scat; /* mapping to subregion */ 21 Vec *x,*y; /* work vectors */ 22 IS *is; /* index set that defines each subdomain */ 23 Mat *mat,*pmat; /* mat is not currently used */ 24 PCASMType type; /* use reduced interpolation, restriction or both */ 25 PetscTruth type_set; /* if user set this value (so won't change it for symmetric problems) */ 26 PetscTruth same_local_solves; /* flag indicating whether all local solvers are same */ 27 PetscTruth inplace; /* indicates that the sub-matrices are deleted after 28 PCSetUpOnBlocks() is done. Similar to inplace 29 factorization in the case of LU and ILU */ 30 } PC_ASM; 31 32 #undef __FUNCT__ 33 #define __FUNCT__ "PCView_ASM" 34 static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 35 { 36 PC_ASM *osm = (PC_ASM*)pc->data; 37 PetscErrorCode ierr; 38 PetscMPIInt rank; 39 PetscInt i; 40 PetscTruth iascii,isstring; 41 PetscViewer sviewer; 42 43 44 PetscFunctionBegin; 45 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 46 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 47 if (iascii) { 48 if (osm->n > 0) { 49 ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: total subdomain blocks = %D, amount of overlap = %D\n",osm->n,osm->overlap);CHKERRQ(ierr); 50 } else { 51 ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: total subdomain blocks not yet set, amount of overlap = %D\n",osm->overlap);CHKERRQ(ierr); 52 } 53 ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr); 54 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 55 if (osm->same_local_solves) { 56 ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr); 57 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 58 if (!rank && osm->ksp) { 59 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 60 ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 61 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 62 } 63 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 64 } else { 65 ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 66 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] number of local blocks = %D\n",rank,osm->n_local_true);CHKERRQ(ierr); 67 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 68 for (i=0; i<osm->n_local; i++) { 69 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 70 if (i < osm->n_local_true) { 71 ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D\n",rank,i);CHKERRQ(ierr); 72 ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 73 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 74 } 75 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 76 } 77 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 78 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 79 } 80 } else if (isstring) { 81 ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%D",osm->n,osm->overlap,osm->type);CHKERRQ(ierr); 82 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 83 if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 84 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 85 } else { 86 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCASM",((PetscObject)viewer)->type_name); 87 } 88 PetscFunctionReturn(0); 89 } 90 91 #undef __FUNCT__ 92 #define __FUNCT__ "PCASMPrintSubdomains" 93 static PetscErrorCode PCASMPrintSubdomains(PC pc) 94 { 95 PC_ASM *osm = (PC_ASM*)pc->data; 96 const char *prefix; 97 char fname[PETSC_MAX_PATH_LEN+1]; 98 PetscViewer viewer; 99 PetscInt i,j,nidx; 100 const PetscInt *idx; 101 PetscErrorCode ierr; 102 PetscFunctionBegin; 103 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 104 ierr = PetscOptionsGetString(prefix,"-pc_asm_print_subdomains", 105 fname,PETSC_MAX_PATH_LEN,PETSC_NULL);CHKERRQ(ierr); 106 if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 107 ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);CHKERRQ(ierr); 108 for (i=0;i<osm->n_local_true;i++) { 109 ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr); 110 ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr); 111 for (j=0; j<nidx; j++) { 112 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"%D ",idx[j]);CHKERRQ(ierr); 113 } 114 ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 115 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 116 } 117 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 118 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 119 PetscFunctionReturn(0); 120 } 121 122 #undef __FUNCT__ 123 #define __FUNCT__ "PCSetUp_ASM" 124 static PetscErrorCode PCSetUp_ASM(PC pc) 125 { 126 PC_ASM *osm = (PC_ASM*)pc->data; 127 PetscErrorCode ierr; 128 PetscTruth symset,flg; 129 PetscInt i,m; 130 PetscMPIInt size; 131 MatReuse scall = MAT_REUSE_MATRIX; 132 IS isl; 133 KSP ksp; 134 PC subpc; 135 const char *prefix,*pprefix; 136 Vec vec; 137 138 PetscFunctionBegin; 139 if (!pc->setupcalled) { 140 141 if (!osm->type_set) { 142 ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 143 if (symset && flg) { osm->type = PC_ASM_BASIC; } 144 } 145 146 if (osm->n == PETSC_DECIDE && osm->n_local_true < 1) { 147 /* no subdomains given, use one per processor */ 148 osm->n_local = osm->n_local_true = 1; 149 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 150 osm->n = size; 151 } else if (osm->n == PETSC_DECIDE) { 152 /* determine global number of subdomains */ 153 PetscInt inwork[2],outwork[2]; 154 inwork[0] = inwork[1] = osm->n_local_true; 155 ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr); 156 osm->n_local = outwork[0]; 157 osm->n = outwork[1]; 158 } 159 160 if (!osm->is){ /* create the index sets */ 161 ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr); 162 } 163 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 164 ierr = PetscOptionsHasName(prefix,"-pc_asm_print_subdomains",&flg);CHKERRQ(ierr); 165 if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); } 166 167 /* Extend the "overlapping" regions by a number of steps */ 168 ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr); 169 for (i=0; i<osm->n_local_true; i++) { 170 ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 171 } 172 173 /* Create the local work vectors and scatter contexts */ 174 ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 175 ierr = PetscMalloc(osm->n_local*sizeof(VecScatter *),&osm->scat);CHKERRQ(ierr); 176 ierr = PetscMalloc(osm->n_local*sizeof(Vec *),&osm->x);CHKERRQ(ierr); 177 ierr = PetscMalloc(osm->n_local*sizeof(Vec *),&osm->y);CHKERRQ(ierr); 178 for (i=0; i<osm->n_local_true; i++) { 179 ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 180 ierr = VecCreateSeq(PETSC_COMM_SELF,m,&osm->x[i]);CHKERRQ(ierr); 181 ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 182 ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 183 ierr = VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->scat[i]);CHKERRQ(ierr); 184 ierr = ISDestroy(isl);CHKERRQ(ierr); 185 } 186 for (i=osm->n_local_true; i<osm->n_local; i++) { 187 ierr = VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);CHKERRQ(ierr); 188 ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 189 ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);CHKERRQ(ierr); 190 ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->scat[i]);CHKERRQ(ierr); 191 ierr = ISDestroy(isl);CHKERRQ(ierr); 192 } 193 ierr = VecDestroy(vec);CHKERRQ(ierr); 194 195 /* Create the local solvers */ 196 ierr = PetscMalloc(osm->n_local_true*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr); 197 for (i=0; i<osm->n_local_true; i++) { 198 ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 199 ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 200 ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 201 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 202 ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 203 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 204 ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 205 ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 206 osm->ksp[i] = ksp; 207 } 208 scall = MAT_INITIAL_MATRIX; 209 210 } else { 211 /* 212 Destroy the blocks from the previous iteration 213 */ 214 if (pc->flag == DIFFERENT_NONZERO_PATTERN || osm->inplace) { 215 if (!osm->inplace) { 216 ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 217 } 218 scall = MAT_INITIAL_MATRIX; 219 } 220 } 221 222 /* 223 Extract out the submatrices 224 */ 225 ierr = MatGetSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 226 if (scall == MAT_INITIAL_MATRIX) { 227 ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 228 for (i=0; i<osm->n_local_true; i++) { 229 ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr); 230 ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 231 } 232 } 233 234 /* Return control to the user so that the submatrices can be modified (e.g., to apply 235 different boundary conditions for the submatrices than for the global problem) */ 236 ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 237 238 /* 239 Loop over subdomains putting them into local ksp 240 */ 241 for (i=0; i<osm->n_local_true; i++) { 242 ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr); 243 if (!pc->setupcalled) { 244 ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 245 } 246 } 247 248 PetscFunctionReturn(0); 249 } 250 251 #undef __FUNCT__ 252 #define __FUNCT__ "PCSetUpOnBlocks_ASM" 253 static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 254 { 255 PC_ASM *osm = (PC_ASM*)pc->data; 256 PetscErrorCode ierr; 257 PetscInt i; 258 259 PetscFunctionBegin; 260 for (i=0; i<osm->n_local_true; i++) { 261 ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 262 } 263 /* 264 If inplace flag is set, then destroy the matrix after the setup 265 on blocks is done. 266 */ 267 if (osm->inplace && osm->n_local_true > 0) { 268 ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 269 } 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "PCApply_ASM" 275 static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 276 { 277 PC_ASM *osm = (PC_ASM*)pc->data; 278 PetscErrorCode ierr; 279 PetscInt i,n_local = osm->n_local,n_local_true = osm->n_local_true; 280 ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 281 282 PetscFunctionBegin; 283 /* 284 Support for limiting the restriction or interpolation to only local 285 subdomain values (leaving the other values 0). 286 */ 287 if (!(osm->type & PC_ASM_RESTRICT)) { 288 forward = SCATTER_FORWARD_LOCAL; 289 /* have to zero the work RHS since scatter may leave some slots empty */ 290 for (i=0; i<n_local_true; i++) { 291 ierr = VecSet(osm->x[i],0.0);CHKERRQ(ierr); 292 } 293 } 294 if (!(osm->type & PC_ASM_INTERPOLATE)) { 295 reverse = SCATTER_REVERSE_LOCAL; 296 } 297 298 for (i=0; i<n_local; i++) { 299 ierr = VecScatterBegin(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 300 } 301 ierr = VecSet(y,0.0);CHKERRQ(ierr); 302 /* do the local solves */ 303 for (i=0; i<n_local_true; i++) { 304 ierr = VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 305 ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 306 ierr = VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 307 } 308 /* handle the rest of the scatters that do not have local solves */ 309 for (i=n_local_true; i<n_local; i++) { 310 ierr = VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 311 ierr = VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 312 } 313 for (i=0; i<n_local; i++) { 314 ierr = VecScatterEnd(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 315 } 316 PetscFunctionReturn(0); 317 } 318 319 #undef __FUNCT__ 320 #define __FUNCT__ "PCApplyTranspose_ASM" 321 static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 322 { 323 PC_ASM *osm = (PC_ASM*)pc->data; 324 PetscErrorCode ierr; 325 PetscInt i,n_local = osm->n_local,n_local_true = osm->n_local_true; 326 ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 327 328 PetscFunctionBegin; 329 /* 330 Support for limiting the restriction or interpolation to only local 331 subdomain values (leaving the other values 0). 332 333 Note: these are reversed from the PCApply_ASM() because we are applying the 334 transpose of the three terms 335 */ 336 if (!(osm->type & PC_ASM_INTERPOLATE)) { 337 forward = SCATTER_FORWARD_LOCAL; 338 /* have to zero the work RHS since scatter may leave some slots empty */ 339 for (i=0; i<n_local_true; i++) { 340 ierr = VecSet(osm->x[i],0.0);CHKERRQ(ierr); 341 } 342 } 343 if (!(osm->type & PC_ASM_RESTRICT)) { 344 reverse = SCATTER_REVERSE_LOCAL; 345 } 346 347 for (i=0; i<n_local; i++) { 348 ierr = VecScatterBegin(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 349 } 350 ierr = VecSet(y,0.0);CHKERRQ(ierr); 351 /* do the local solves */ 352 for (i=0; i<n_local_true; i++) { 353 ierr = VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 354 ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 355 ierr = VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 356 } 357 /* handle the rest of the scatters that do not have local solves */ 358 for (i=n_local_true; i<n_local; i++) { 359 ierr = VecScatterEnd(osm->scat[i],x,osm->x[i],INSERT_VALUES,forward);CHKERRQ(ierr); 360 ierr = VecScatterBegin(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 361 } 362 for (i=0; i<n_local; i++) { 363 ierr = VecScatterEnd(osm->scat[i],osm->y[i],y,ADD_VALUES,reverse);CHKERRQ(ierr); 364 } 365 PetscFunctionReturn(0); 366 } 367 368 #undef __FUNCT__ 369 #define __FUNCT__ "PCDestroy_ASM" 370 static PetscErrorCode PCDestroy_ASM(PC pc) 371 { 372 PC_ASM *osm = (PC_ASM*)pc->data; 373 PetscErrorCode ierr; 374 PetscInt i; 375 376 PetscFunctionBegin; 377 if (osm->ksp) { 378 for (i=0; i<osm->n_local_true; i++) { 379 ierr = KSPDestroy(osm->ksp[i]);CHKERRQ(ierr); 380 } 381 ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 382 } 383 if (osm->pmat && !osm->inplace) { 384 if (osm->n_local_true > 0) { 385 ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 386 } 387 } 388 if (osm->scat) { 389 for (i=0; i<osm->n_local; i++) { 390 ierr = VecScatterDestroy(osm->scat[i]);CHKERRQ(ierr); 391 ierr = VecDestroy(osm->x[i]);CHKERRQ(ierr); 392 ierr = VecDestroy(osm->y[i]);CHKERRQ(ierr); 393 } 394 ierr = PetscFree(osm->scat);CHKERRQ(ierr); 395 ierr = PetscFree(osm->x);CHKERRQ(ierr); 396 ierr = PetscFree(osm->y);CHKERRQ(ierr); 397 } 398 if (osm->is) { 399 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is);CHKERRQ(ierr); 400 } 401 ierr = PetscFree(osm);CHKERRQ(ierr); 402 PetscFunctionReturn(0); 403 } 404 405 #undef __FUNCT__ 406 #define __FUNCT__ "PCSetFromOptions_ASM" 407 static PetscErrorCode PCSetFromOptions_ASM(PC pc) 408 { 409 PC_ASM *osm = (PC_ASM*)pc->data; 410 PetscErrorCode ierr; 411 PetscInt blocks,ovl; 412 PetscTruth symset,flg; 413 PCASMType asmtype; 414 415 PetscFunctionBegin; 416 /* set the type to symmetric if matrix is symmetric */ 417 if (!osm->type_set && pc->pmat) { 418 ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 419 if (symset && flg) { osm->type = PC_ASM_BASIC; } 420 } 421 ierr = PetscOptionsHead("Additive Schwarz options");CHKERRQ(ierr); 422 ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains", 423 osm->n,&blocks,&flg);CHKERRQ(ierr); 424 if (flg) {ierr = PCASMSetTotalSubdomains(pc,blocks,PETSC_NULL);CHKERRQ(ierr); } 425 ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap", 426 osm->overlap,&ovl,&flg);CHKERRQ(ierr); 427 if (flg) {ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); } 428 ierr = PetscOptionsName("-pc_asm_in_place","Perform matrix factorization inplace","PCASMSetUseInPlace",&flg);CHKERRQ(ierr); 429 if (flg) {ierr = PCASMSetUseInPlace(pc);CHKERRQ(ierr); } 430 ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType", 431 PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 432 if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 433 ierr = PetscOptionsTail();CHKERRQ(ierr); 434 PetscFunctionReturn(0); 435 } 436 437 /*------------------------------------------------------------------------------------*/ 438 439 EXTERN_C_BEGIN 440 #undef __FUNCT__ 441 #define __FUNCT__ "PCASMSetLocalSubdomains_ASM" 442 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[]) 443 { 444 PC_ASM *osm = (PC_ASM*)pc->data; 445 PetscErrorCode ierr; 446 PetscInt i; 447 448 PetscFunctionBegin; 449 if (n < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 450 if (pc->setupcalled && (n != osm->n_local_true || is)) { 451 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetup()."); 452 } 453 if (!pc->setupcalled) { 454 if (is) { 455 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 456 } 457 if (osm->is) { 458 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is);CHKERRQ(ierr); 459 } 460 osm->n_local_true = n; 461 osm->is = 0; 462 if (is) { 463 ierr = PetscMalloc(n*sizeof(IS *),&osm->is);CHKERRQ(ierr); 464 for (i=0; i<n; i++) { osm->is[i] = is[i]; } 465 } 466 } 467 PetscFunctionReturn(0); 468 } 469 EXTERN_C_END 470 471 EXTERN_C_BEGIN 472 #undef __FUNCT__ 473 #define __FUNCT__ "PCASMSetTotalSubdomains_ASM" 474 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is) 475 { 476 PC_ASM *osm = (PC_ASM*)pc->data; 477 PetscErrorCode ierr; 478 PetscMPIInt rank,size; 479 PetscInt n; 480 481 PetscFunctionBegin; 482 if (N < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 483 if (is) SETERRQ(PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet."); 484 485 /* 486 Split the subdomains equally amoung all processors 487 */ 488 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 489 ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 490 n = N/size + ((N % size) > rank); 491 if (!n) SETERRQ3(PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N); 492 if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetup()."); 493 if (!pc->setupcalled) { 494 if (osm->is) { 495 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is);CHKERRQ(ierr); 496 } 497 osm->n_local_true = n; 498 osm->is = 0; 499 } 500 PetscFunctionReturn(0); 501 } 502 EXTERN_C_END 503 504 EXTERN_C_BEGIN 505 #undef __FUNCT__ 506 #define __FUNCT__ "PCASMSetOverlap_ASM" 507 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 508 { 509 PC_ASM *osm = (PC_ASM*)pc->data; 510 511 PetscFunctionBegin; 512 if (ovl < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 513 if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetup()."); 514 if (!pc->setupcalled) { 515 osm->overlap = ovl; 516 } 517 PetscFunctionReturn(0); 518 } 519 EXTERN_C_END 520 521 EXTERN_C_BEGIN 522 #undef __FUNCT__ 523 #define __FUNCT__ "PCASMSetType_ASM" 524 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType_ASM(PC pc,PCASMType type) 525 { 526 PC_ASM *osm = (PC_ASM*)pc->data; 527 528 PetscFunctionBegin; 529 osm->type = type; 530 osm->type_set = PETSC_TRUE; 531 PetscFunctionReturn(0); 532 } 533 EXTERN_C_END 534 535 EXTERN_C_BEGIN 536 #undef __FUNCT__ 537 #define __FUNCT__ "PCASMGetSubKSP_ASM" 538 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 539 { 540 PC_ASM *osm = (PC_ASM*)pc->data; 541 PetscErrorCode ierr; 542 543 PetscFunctionBegin; 544 if (osm->n_local_true < 1) { 545 SETERRQ(PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 546 } 547 548 if (n_local) { 549 *n_local = osm->n_local_true; 550 } 551 if (first_local) { 552 ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 553 *first_local -= osm->n_local_true; 554 } 555 if (ksp) { 556 /* Assume that local solves are now different; not necessarily 557 true though! This flag is used only for PCView_ASM() */ 558 *ksp = osm->ksp; 559 osm->same_local_solves = PETSC_FALSE; 560 } 561 PetscFunctionReturn(0); 562 } 563 EXTERN_C_END 564 565 EXTERN_C_BEGIN 566 #undef __FUNCT__ 567 #define __FUNCT__ "PCASMSetUseInPlace_ASM" 568 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace_ASM(PC pc) 569 { 570 PC_ASM *osm = (PC_ASM*)pc->data; 571 572 PetscFunctionBegin; 573 osm->inplace = PETSC_TRUE; 574 PetscFunctionReturn(0); 575 } 576 EXTERN_C_END 577 578 /*----------------------------------------------------------------------------*/ 579 #undef __FUNCT__ 580 #define __FUNCT__ "PCASMSetUseInPlace" 581 /*@ 582 PCASMSetUseInPlace - Tells the system to destroy the matrix after setup is done. 583 584 Collective on PC 585 586 Input Parameters: 587 . pc - the preconditioner context 588 589 Options Database Key: 590 . -pc_asm_in_place - Activates in-place factorization 591 592 Note: 593 PCASMSetUseInplace() can only be used with the KSP method KSPPREONLY, and 594 when the original matrix is not required during the Solve process. 595 This destroys the matrix, early thus, saving on memory usage. 596 597 Level: intermediate 598 599 .keywords: PC, set, factorization, direct, inplace, in-place, ASM 600 601 .seealso: PCFactorSetUseInPlace() 602 @*/ 603 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetUseInPlace(PC pc) 604 { 605 PetscErrorCode ierr,(*f)(PC); 606 607 PetscFunctionBegin; 608 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 609 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 610 if (f) { 611 ierr = (*f)(pc);CHKERRQ(ierr); 612 } 613 PetscFunctionReturn(0); 614 } 615 /*----------------------------------------------------------------------------*/ 616 617 #undef __FUNCT__ 618 #define __FUNCT__ "PCASMSetLocalSubdomains" 619 /*@C 620 PCASMSetLocalSubdomains - Sets the local subdomains (for this processor 621 only) for the additive Schwarz preconditioner. 622 623 Collective on PC 624 625 Input Parameters: 626 + pc - the preconditioner context 627 . n - the number of subdomains for this processor (default value = 1) 628 - is - the index sets that define the subdomains for this processor 629 (or PETSC_NULL for PETSc to determine subdomains) 630 631 Notes: 632 The IS numbering is in the parallel, global numbering of the vector. 633 634 By default the ASM preconditioner uses 1 block per processor. 635 636 These index sets cannot be destroyed until after completion of the 637 linear solves for which the ASM preconditioner is being used. 638 639 Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 640 641 Level: advanced 642 643 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 644 645 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 646 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 647 @*/ 648 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[]) 649 { 650 PetscErrorCode ierr,(*f)(PC,PetscInt,IS[]); 651 652 PetscFunctionBegin; 653 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 654 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 655 if (f) { 656 ierr = (*f)(pc,n,is);CHKERRQ(ierr); 657 } 658 PetscFunctionReturn(0); 659 } 660 661 #undef __FUNCT__ 662 #define __FUNCT__ "PCASMSetTotalSubdomains" 663 /*@C 664 PCASMSetTotalSubdomains - Sets the subdomains for all processor for the 665 additive Schwarz preconditioner. Either all or no processors in the 666 PC communicator must call this routine, with the same index sets. 667 668 Collective on PC 669 670 Input Parameters: 671 + pc - the preconditioner context 672 . n - the number of subdomains for all processors 673 - is - the index sets that define the subdomains for all processor 674 (or PETSC_NULL for PETSc to determine subdomains) 675 676 Options Database Key: 677 To set the total number of subdomain blocks rather than specify the 678 index sets, use the option 679 . -pc_asm_blocks <blks> - Sets total blocks 680 681 Notes: 682 Currently you cannot use this to set the actual subdomains with the argument is. 683 684 By default the ASM preconditioner uses 1 block per processor. 685 686 These index sets cannot be destroyed until after completion of the 687 linear solves for which the ASM preconditioner is being used. 688 689 Use PCASMSetLocalSubdomains() to set local subdomains. 690 691 Level: advanced 692 693 .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 694 695 .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 696 PCASMCreateSubdomains2D() 697 @*/ 698 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetTotalSubdomains(PC pc,PetscInt N,IS *is) 699 { 700 PetscErrorCode ierr,(*f)(PC,PetscInt,IS *); 701 702 PetscFunctionBegin; 703 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 704 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 705 if (f) { 706 ierr = (*f)(pc,N,is);CHKERRQ(ierr); 707 } 708 PetscFunctionReturn(0); 709 } 710 711 #undef __FUNCT__ 712 #define __FUNCT__ "PCASMSetOverlap" 713 /*@ 714 PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 715 additive Schwarz preconditioner. Either all or no processors in the 716 PC communicator must call this routine. 717 718 Collective on PC 719 720 Input Parameters: 721 + pc - the preconditioner context 722 - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 723 724 Options Database Key: 725 . -pc_asm_overlap <ovl> - Sets overlap 726 727 Notes: 728 By default the ASM preconditioner uses 1 block per processor. To use 729 multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 730 PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 731 732 The overlap defaults to 1, so if one desires that no additional 733 overlap be computed beyond what may have been set with a call to 734 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 735 must be set to be 0. In particular, if one does not explicitly set 736 the subdomains an application code, then all overlap would be computed 737 internally by PETSc, and using an overlap of 0 would result in an ASM 738 variant that is equivalent to the block Jacobi preconditioner. 739 740 Note that one can define initial index sets with any overlap via 741 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine 742 PCASMSetOverlap() merely allows PETSc to extend that overlap further 743 if desired. 744 745 Level: intermediate 746 747 .keywords: PC, ASM, set, overlap 748 749 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 750 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 751 @*/ 752 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetOverlap(PC pc,PetscInt ovl) 753 { 754 PetscErrorCode ierr,(*f)(PC,PetscInt); 755 756 PetscFunctionBegin; 757 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 758 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetOverlap_C",(void (**)(void))&f);CHKERRQ(ierr); 759 if (f) { 760 ierr = (*f)(pc,ovl);CHKERRQ(ierr); 761 } 762 PetscFunctionReturn(0); 763 } 764 765 #undef __FUNCT__ 766 #define __FUNCT__ "PCASMSetType" 767 /*@ 768 PCASMSetType - Sets the type of restriction and interpolation used 769 for local problems in the additive Schwarz method. 770 771 Collective on PC 772 773 Input Parameters: 774 + pc - the preconditioner context 775 - type - variant of ASM, one of 776 .vb 777 PC_ASM_BASIC - full interpolation and restriction 778 PC_ASM_RESTRICT - full restriction, local processor interpolation 779 PC_ASM_INTERPOLATE - full interpolation, local processor restriction 780 PC_ASM_NONE - local processor restriction and interpolation 781 .ve 782 783 Options Database Key: 784 . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 785 786 Level: intermediate 787 788 .keywords: PC, ASM, set, type 789 790 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 791 PCASMCreateSubdomains2D() 792 @*/ 793 PetscErrorCode PETSCKSP_DLLEXPORT PCASMSetType(PC pc,PCASMType type) 794 { 795 PetscErrorCode ierr,(*f)(PC,PCASMType); 796 797 PetscFunctionBegin; 798 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 799 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 800 if (f) { 801 ierr = (*f)(pc,type);CHKERRQ(ierr); 802 } 803 PetscFunctionReturn(0); 804 } 805 806 #undef __FUNCT__ 807 #define __FUNCT__ "PCASMGetSubKSP" 808 /*@C 809 PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 810 this processor. 811 812 Collective on PC iff first_local is requested 813 814 Input Parameter: 815 . pc - the preconditioner context 816 817 Output Parameters: 818 + n_local - the number of blocks on this processor or PETSC_NULL 819 . first_local - the global number of the first block on this processor or PETSC_NULL, 820 all processors must request or all must pass PETSC_NULL 821 - ksp - the array of KSP contexts 822 823 Note: 824 After PCASMGetSubKSP() the array of KSPes is not to be freed 825 826 Currently for some matrix implementations only 1 block per processor 827 is supported. 828 829 You must call KSPSetUp() before calling PCASMGetSubKSP(). 830 831 Level: advanced 832 833 .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 834 835 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 836 PCASMCreateSubdomains2D(), 837 @*/ 838 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 839 { 840 PetscErrorCode ierr,(*f)(PC,PetscInt*,PetscInt*,KSP **); 841 842 PetscFunctionBegin; 843 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 844 ierr = PetscObjectQueryFunction((PetscObject)pc,"PCASMGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 845 if (f) { 846 ierr = (*f)(pc,n_local,first_local,ksp);CHKERRQ(ierr); 847 } else { 848 SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 849 } 850 851 PetscFunctionReturn(0); 852 } 853 854 /* -------------------------------------------------------------------------------------*/ 855 /*MC 856 PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 857 its own KSP object. 858 859 Options Database Keys: 860 + -pc_asm_truelocal - Activates PCASMSetUseTrueLocal() 861 . -pc_asm_in_place - Activates in-place factorization 862 . -pc_asm_blocks <blks> - Sets total blocks 863 . -pc_asm_overlap <ovl> - Sets overlap 864 - -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 865 866 IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 867 will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 868 -pc_asm_type basic to use the standard ASM. 869 870 Notes: Each processor can have one or more blocks, but a block cannot be shared by more 871 than one processor. Defaults to one block per processor. 872 873 To set options on the solvers for each block append -sub_ to all the KSP, and PC 874 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 875 876 To set the options on the solvers separate for each block call PCASMGetSubKSP() 877 and set the options directly on the resulting KSP object (you can access its PC 878 with KSPGetPC()) 879 880 881 Level: beginner 882 883 Concepts: additive Schwarz method 884 885 References: 886 An additive variant of the Schwarz alternating method for the case of many subregions 887 M Dryja, OB Widlund - Courant Institute, New York University Technical report 888 889 Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 890 Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 891 892 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 893 PCBJACOBI, PCASMSetUseTrueLocal(), PCASMGetSubKSP(), PCASMSetLocalSubdomains(), 894 PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), 895 PCASMSetUseInPlace() 896 M*/ 897 898 EXTERN_C_BEGIN 899 #undef __FUNCT__ 900 #define __FUNCT__ "PCCreate_ASM" 901 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ASM(PC pc) 902 { 903 PetscErrorCode ierr; 904 PC_ASM *osm; 905 906 PetscFunctionBegin; 907 ierr = PetscNewLog(pc,PC_ASM,&osm);CHKERRQ(ierr); 908 osm->n = PETSC_DECIDE; 909 osm->n_local = 0; 910 osm->n_local_true = 0; 911 osm->overlap = 1; 912 osm->ksp = 0; 913 osm->scat = 0; 914 osm->x = 0; 915 osm->y = 0; 916 osm->is = 0; 917 osm->mat = 0; 918 osm->pmat = 0; 919 osm->type = PC_ASM_RESTRICT; 920 osm->same_local_solves = PETSC_TRUE; 921 osm->inplace = PETSC_FALSE; 922 923 pc->data = (void*)osm; 924 pc->ops->apply = PCApply_ASM; 925 pc->ops->applytranspose = PCApplyTranspose_ASM; 926 pc->ops->setup = PCSetUp_ASM; 927 pc->ops->destroy = PCDestroy_ASM; 928 pc->ops->setfromoptions = PCSetFromOptions_ASM; 929 pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 930 pc->ops->view = PCView_ASM; 931 pc->ops->applyrichardson = 0; 932 933 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM", 934 PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 935 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM", 936 PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 937 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM", 938 PCASMSetOverlap_ASM);CHKERRQ(ierr); 939 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM", 940 PCASMSetType_ASM);CHKERRQ(ierr); 941 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubKSP_C","PCASMGetSubKSP_ASM", 942 PCASMGetSubKSP_ASM);CHKERRQ(ierr); 943 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetUseInPlace_C","PCASMSetUseInPlace_ASM", 944 PCASMSetUseInPlace_ASM);CHKERRQ(ierr); 945 PetscFunctionReturn(0); 946 } 947 EXTERN_C_END 948 949 950 #undef __FUNCT__ 951 #define __FUNCT__ "PCASMCreateSubdomains" 952 /*@C 953 PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 954 preconditioner for a any problem on a general grid. 955 956 Collective 957 958 Input Parameters: 959 + A - The global matrix operator 960 - n - the number of local blocks 961 962 Output Parameters: 963 . outis - the array of index sets defining the subdomains 964 965 Level: advanced 966 967 Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 968 from these if you use PCASMSetLocalSubdomains() 969 970 In the Fortran version you must provide the array outis[] already allocated of length n. 971 972 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 973 974 .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 975 @*/ 976 PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 977 { 978 MatPartitioning mpart; 979 const char *prefix; 980 PetscErrorCode (*f)(Mat,PetscTruth*,MatReuse,Mat*); 981 PetscMPIInt size; 982 PetscInt i,j,rstart,rend,bs; 983 PetscTruth iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 984 Mat Ad = PETSC_NULL, adj; 985 IS ispart,isnumb,*is; 986 PetscErrorCode ierr; 987 988 PetscFunctionBegin; 989 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 990 PetscValidPointer(outis,3); 991 if (n < 1) {SETERRQ1(PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n);} 992 993 /* Get prefix, row distribution, and block size */ 994 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 995 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 996 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 997 if (rstart/bs*bs != rstart || rend/bs*bs != rend) { 998 SETERRQ3(PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs); 999 } 1000 /* Get diagonal block from matrix if possible */ 1001 ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1002 ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 1003 if (f) { 1004 ierr = (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);CHKERRQ(ierr); 1005 } else if (size == 1) { 1006 iscopy = PETSC_FALSE; Ad = A; 1007 } else { 1008 iscopy = PETSC_FALSE; Ad = PETSC_NULL; 1009 } 1010 if (Ad) { 1011 ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1012 if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1013 } 1014 if (Ad && n > 1) { 1015 PetscTruth match,done; 1016 /* Try to setup a good matrix partitioning if available */ 1017 ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1018 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1019 ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1020 ierr = PetscTypeCompare((PetscObject)mpart,MAT_PARTITIONING_CURRENT,&match);CHKERRQ(ierr); 1021 if (!match) { 1022 ierr = PetscTypeCompare((PetscObject)mpart,MAT_PARTITIONING_SQUARE,&match);CHKERRQ(ierr); 1023 } 1024 if (!match) { /* assume a "good" partitioner is available */ 1025 PetscInt na,*ia,*ja; 1026 ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1027 if (done) { 1028 /* Build adjacency matrix by hand. Unfortunately a call to 1029 MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1030 remove the block-aij structure and we cannot expect 1031 MatPartitioning to split vertices as we need */ 1032 PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0; 1033 nnz = 0; 1034 for (i=0; i<na; i++) { /* count number of nonzeros */ 1035 len = ia[i+1] - ia[i]; 1036 row = ja + ia[i]; 1037 for (j=0; j<len; j++) { 1038 if (row[j] == i) { /* don't count diagonal */ 1039 len--; break; 1040 } 1041 } 1042 nnz += len; 1043 } 1044 ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1045 ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1046 nnz = 0; 1047 iia[0] = 0; 1048 for (i=0; i<na; i++) { /* fill adjacency */ 1049 cnt = 0; 1050 len = ia[i+1] - ia[i]; 1051 row = ja + ia[i]; 1052 for (j=0; j<len; j++) { 1053 if (row[j] != i) { /* if not diagonal */ 1054 jja[nnz+cnt++] = row[j]; 1055 } 1056 } 1057 nnz += cnt; 1058 iia[i+1] = nnz; 1059 } 1060 /* Partitioning of the adjacency matrix */ 1061 ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr); 1062 ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1063 ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1064 ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1065 ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1066 ierr = MatDestroy(adj);CHKERRQ(ierr); 1067 foundpart = PETSC_TRUE; 1068 } 1069 ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1070 } 1071 ierr = MatPartitioningDestroy(mpart);CHKERRQ(ierr); 1072 } 1073 if (iscopy) {ierr = MatDestroy(Ad);CHKERRQ(ierr);} 1074 1075 ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1076 *outis = is; 1077 1078 if (!foundpart) { 1079 1080 /* Partitioning by contiguous chunks of rows */ 1081 1082 PetscInt mbs = (rend-rstart)/bs; 1083 PetscInt start = rstart; 1084 for (i=0; i<n; i++) { 1085 PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1086 ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1087 start += count; 1088 } 1089 1090 } else { 1091 1092 /* Partitioning by adjacency of diagonal block */ 1093 1094 const PetscInt *numbering; 1095 PetscInt *count,nidx,*indices,*newidx,start=0; 1096 /* Get node count in each partition */ 1097 ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1098 ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1099 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1100 for (i=0; i<n; i++) count[i] *= bs; 1101 } 1102 /* Build indices from node numbering */ 1103 ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1104 ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1105 for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1106 ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1107 ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1108 ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1109 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1110 ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1111 for (i=0; i<nidx; i++) 1112 for (j=0; j<bs; j++) 1113 newidx[i*bs+j] = indices[i]*bs + j; 1114 ierr = PetscFree(indices);CHKERRQ(ierr); 1115 nidx *= bs; 1116 indices = newidx; 1117 } 1118 /* Shift to get global indices */ 1119 for (i=0; i<nidx; i++) indices[i] += rstart; 1120 1121 /* Build the index sets for each block */ 1122 for (i=0; i<n; i++) { 1123 ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],&is[i]);CHKERRQ(ierr); 1124 ierr = ISSort(is[i]);CHKERRQ(ierr); 1125 start += count[i]; 1126 } 1127 1128 ierr = PetscFree(count); 1129 ierr = PetscFree(indices); 1130 ierr = ISDestroy(isnumb);CHKERRQ(ierr); 1131 ierr = ISDestroy(ispart);CHKERRQ(ierr); 1132 1133 } 1134 1135 PetscFunctionReturn(0); 1136 } 1137 1138 #undef __FUNCT__ 1139 #define __FUNCT__ "PCASMDestroySubdomains" 1140 /*@C 1141 PCASMDestroySubdomains - Destroys the index sets created with 1142 PCASMCreateSubdomains(). Should be called after setting subdomains 1143 with PCASMSetLocalSubdomains(). 1144 1145 Collective 1146 1147 Input Parameters: 1148 + n - the number of index sets 1149 - is - the array of index sets 1150 1151 Level: advanced 1152 1153 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 1154 1155 .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 1156 @*/ 1157 PetscErrorCode PETSCKSP_DLLEXPORT PCASMDestroySubdomains(PetscInt n, IS is[]) 1158 { 1159 PetscInt i; 1160 PetscErrorCode ierr; 1161 PetscFunctionBegin; 1162 if (n <= 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n); 1163 PetscValidPointer(is,2); 1164 for (i=0; i<n; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); } 1165 ierr = PetscFree(is);CHKERRQ(ierr); 1166 PetscFunctionReturn(0); 1167 } 1168 1169 #undef __FUNCT__ 1170 #define __FUNCT__ "PCASMCreateSubdomains2D" 1171 /*@ 1172 PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1173 preconditioner for a two-dimensional problem on a regular grid. 1174 1175 Not Collective 1176 1177 Input Parameters: 1178 + m, n - the number of mesh points in the x and y directions 1179 . M, N - the number of subdomains in the x and y directions 1180 . dof - degrees of freedom per node 1181 - overlap - overlap in mesh lines 1182 1183 Output Parameters: 1184 + Nsub - the number of subdomains created 1185 - is - the array of index sets defining the subdomains 1186 1187 Note: 1188 Presently PCAMSCreateSubdomains2d() is valid only for sequential 1189 preconditioners. More general related routines are 1190 PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 1191 1192 Level: advanced 1193 1194 .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 1195 1196 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1197 PCASMSetOverlap() 1198 @*/ 1199 PetscErrorCode PETSCKSP_DLLEXPORT PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is) 1200 { 1201 PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outter; 1202 PetscErrorCode ierr; 1203 PetscInt nidx,*idx,loc,ii,jj,count; 1204 1205 PetscFunctionBegin; 1206 if (dof != 1) SETERRQ(PETSC_ERR_SUP," "); 1207 1208 *Nsub = N*M; 1209 ierr = PetscMalloc((*Nsub)*sizeof(IS *),is);CHKERRQ(ierr); 1210 ystart = 0; 1211 loc_outter = 0; 1212 for (i=0; i<N; i++) { 1213 height = n/N + ((n % N) > i); /* height of subdomain */ 1214 if (height < 2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 1215 yleft = ystart - overlap; if (yleft < 0) yleft = 0; 1216 yright = ystart + height + overlap; if (yright > n) yright = n; 1217 xstart = 0; 1218 for (j=0; j<M; j++) { 1219 width = m/M + ((m % M) > j); /* width of subdomain */ 1220 if (width < 2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 1221 xleft = xstart - overlap; if (xleft < 0) xleft = 0; 1222 xright = xstart + width + overlap; if (xright > m) xright = m; 1223 nidx = (xright - xleft)*(yright - yleft); 1224 ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1225 loc = 0; 1226 for (ii=yleft; ii<yright; ii++) { 1227 count = m*ii + xleft; 1228 for (jj=xleft; jj<xright; jj++) { 1229 idx[loc++] = count++; 1230 } 1231 } 1232 ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,(*is)+loc_outter++);CHKERRQ(ierr); 1233 ierr = PetscFree(idx);CHKERRQ(ierr); 1234 xstart += width; 1235 } 1236 ystart += height; 1237 } 1238 for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 1239 PetscFunctionReturn(0); 1240 } 1241 1242 #undef __FUNCT__ 1243 #define __FUNCT__ "PCASMGetLocalSubdomains" 1244 /*@C 1245 PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1246 only) for the additive Schwarz preconditioner. 1247 1248 Collective on PC 1249 1250 Input Parameter: 1251 . pc - the preconditioner context 1252 1253 Output Parameters: 1254 + n - the number of subdomains for this processor (default value = 1) 1255 - is - the index sets that define the subdomains for this processor 1256 1257 1258 Notes: 1259 The IS numbering is in the parallel, global numbering of the vector. 1260 1261 Level: advanced 1262 1263 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 1264 1265 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1266 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 1267 @*/ 1268 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[]) 1269 { 1270 PC_ASM *osm; 1271 PetscErrorCode ierr; 1272 PetscTruth match; 1273 1274 PetscFunctionBegin; 1275 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1276 PetscValidIntPointer(n,2); 1277 if (is) PetscValidPointer(is,3); 1278 ierr = PetscTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1279 if (!match) { 1280 if (n) *n = 0; 1281 if (is) *is = PETSC_NULL; 1282 } else { 1283 osm = (PC_ASM*)pc->data; 1284 if (n) *n = osm->n_local_true; 1285 if (is) *is = osm->is; 1286 } 1287 PetscFunctionReturn(0); 1288 } 1289 1290 #undef __FUNCT__ 1291 #define __FUNCT__ "PCASMGetLocalSubmatrices" 1292 /*@C 1293 PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1294 only) for the additive Schwarz preconditioner. 1295 1296 Collective on PC 1297 1298 Input Parameter: 1299 . pc - the preconditioner context 1300 1301 Output Parameters: 1302 + n - the number of matrices for this processor (default value = 1) 1303 - mat - the matrices 1304 1305 1306 Level: advanced 1307 1308 .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 1309 1310 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1311 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1312 @*/ 1313 PetscErrorCode PETSCKSP_DLLEXPORT PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1314 { 1315 PC_ASM *osm; 1316 PetscErrorCode ierr; 1317 PetscTruth match; 1318 1319 PetscFunctionBegin; 1320 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1321 PetscValidIntPointer(n,2); 1322 if (mat) PetscValidPointer(mat,3); 1323 if (!pc->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1324 ierr = PetscTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1325 if (!match) { 1326 if (n) *n = 0; 1327 if (mat) *mat = PETSC_NULL; 1328 } else { 1329 osm = (PC_ASM*)pc->data; 1330 if (n) *n = osm->n_local_true; 1331 if (mat) *mat = osm->pmat; 1332 } 1333 PetscFunctionReturn(0); 1334 } 1335