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