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