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