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