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 dm_subdomains; /* 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(PetscObjectComm((PetscObject)pc),&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(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 111 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &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(PetscObjectComm((PetscObject)pc),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->dm_subdomains && 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,PetscObjectComm((PetscObject)pc));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 = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 589 ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 590 if (flg) { 591 ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 592 osm->dm_subdomains = PETSC_FALSE; 593 } 594 ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 595 if (flg) { 596 ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); 597 osm->dm_subdomains = PETSC_FALSE; 598 } 599 flg = PETSC_FALSE; 600 ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 601 if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 602 ierr = PetscOptionsTail();CHKERRQ(ierr); 603 PetscFunctionReturn(0); 604 } 605 606 /*------------------------------------------------------------------------------------*/ 607 608 EXTERN_C_BEGIN 609 #undef __FUNCT__ 610 #define __FUNCT__ "PCASMSetLocalSubdomains_ASM" 611 PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 612 { 613 PC_ASM *osm = (PC_ASM*)pc->data; 614 PetscErrorCode ierr; 615 PetscInt i; 616 617 PetscFunctionBegin; 618 if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 619 if (pc->setupcalled && (n != osm->n_local_true || is)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 620 621 if (!pc->setupcalled) { 622 if (is) { 623 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 624 } 625 if (is_local) { 626 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 627 } 628 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 629 630 osm->n_local_true = n; 631 osm->is = 0; 632 osm->is_local = 0; 633 if (is) { 634 ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr); 635 for (i=0; i<n; i++) osm->is[i] = is[i]; 636 /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 637 osm->overlap = -1; 638 } 639 if (is_local) { 640 ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr); 641 for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 642 if (!is) { 643 ierr = PetscMalloc(osm->n_local_true*sizeof(IS),&osm->is);CHKERRQ(ierr); 644 for (i=0; i<osm->n_local_true; i++) { 645 if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 646 ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr); 647 ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr); 648 } else { 649 ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr); 650 osm->is[i] = osm->is_local[i]; 651 } 652 } 653 } 654 } 655 } 656 PetscFunctionReturn(0); 657 } 658 EXTERN_C_END 659 660 EXTERN_C_BEGIN 661 #undef __FUNCT__ 662 #define __FUNCT__ "PCASMSetTotalSubdomains_ASM" 663 PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 664 { 665 PC_ASM *osm = (PC_ASM*)pc->data; 666 PetscErrorCode ierr; 667 PetscMPIInt rank,size; 668 PetscInt n; 669 670 PetscFunctionBegin; 671 if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 672 if (is || is_local) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet."); 673 674 /* 675 Split the subdomains equally among all processors 676 */ 677 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 678 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 679 n = N/size + ((N % size) > rank); 680 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); 681 if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 682 if (!pc->setupcalled) { 683 ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 684 685 osm->n_local_true = n; 686 osm->is = 0; 687 osm->is_local = 0; 688 } 689 PetscFunctionReturn(0); 690 } 691 EXTERN_C_END 692 693 EXTERN_C_BEGIN 694 #undef __FUNCT__ 695 #define __FUNCT__ "PCASMSetOverlap_ASM" 696 PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 697 { 698 PC_ASM *osm = (PC_ASM*)pc->data; 699 700 PetscFunctionBegin; 701 if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 702 if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 703 if (!pc->setupcalled) osm->overlap = ovl; 704 PetscFunctionReturn(0); 705 } 706 EXTERN_C_END 707 708 EXTERN_C_BEGIN 709 #undef __FUNCT__ 710 #define __FUNCT__ "PCASMSetType_ASM" 711 PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 712 { 713 PC_ASM *osm = (PC_ASM*)pc->data; 714 715 PetscFunctionBegin; 716 osm->type = type; 717 osm->type_set = PETSC_TRUE; 718 PetscFunctionReturn(0); 719 } 720 EXTERN_C_END 721 722 EXTERN_C_BEGIN 723 #undef __FUNCT__ 724 #define __FUNCT__ "PCASMSetSortIndices_ASM" 725 PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 726 { 727 PC_ASM *osm = (PC_ASM*)pc->data; 728 729 PetscFunctionBegin; 730 osm->sort_indices = doSort; 731 PetscFunctionReturn(0); 732 } 733 EXTERN_C_END 734 735 EXTERN_C_BEGIN 736 #undef __FUNCT__ 737 #define __FUNCT__ "PCASMGetSubKSP_ASM" 738 PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 739 { 740 PC_ASM *osm = (PC_ASM*)pc->data; 741 PetscErrorCode ierr; 742 743 PetscFunctionBegin; 744 if (osm->n_local_true < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 745 746 if (n_local) *n_local = osm->n_local_true; 747 if (first_local) { 748 ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 749 *first_local -= osm->n_local_true; 750 } 751 if (ksp) { 752 /* Assume that local solves are now different; not necessarily 753 true though! This flag is used only for PCView_ASM() */ 754 *ksp = osm->ksp; 755 osm->same_local_solves = PETSC_FALSE; 756 } 757 PetscFunctionReturn(0); 758 } 759 EXTERN_C_END 760 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "PCASMSetLocalSubdomains" 764 /*@C 765 PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 766 767 Collective on PC 768 769 Input Parameters: 770 + pc - the preconditioner context 771 . n - the number of subdomains for this processor (default value = 1) 772 . is - the index set that defines the subdomains for this processor 773 (or NULL for PETSc to determine subdomains) 774 - is_local - the index sets that define the local part of the subdomains for this processor 775 (or NULL to use the default of 1 subdomain per process) 776 777 Notes: 778 The IS numbering is in the parallel, global numbering of the vector for both is and is_local 779 780 By default the ASM preconditioner uses 1 block per processor. 781 782 Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 783 784 Level: advanced 785 786 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 787 788 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 789 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 790 @*/ 791 PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 792 { 793 PetscErrorCode ierr; 794 795 PetscFunctionBegin; 796 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 797 ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 798 PetscFunctionReturn(0); 799 } 800 801 #undef __FUNCT__ 802 #define __FUNCT__ "PCASMSetTotalSubdomains" 803 /*@C 804 PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 805 additive Schwarz preconditioner. Either all or no processors in the 806 PC communicator must call this routine, with the same index sets. 807 808 Collective on PC 809 810 Input Parameters: 811 + pc - the preconditioner context 812 . N - the number of subdomains for all processors 813 . is - the index sets that define the subdomains for all processors 814 (or NULL to ask PETSc to compe up with subdomains) 815 - is_local - the index sets that define the local part of the subdomains for this processor 816 (or NULL to use the default of 1 subdomain per process) 817 818 Options Database Key: 819 To set the total number of subdomain blocks rather than specify the 820 index sets, use the option 821 . -pc_asm_blocks <blks> - Sets total blocks 822 823 Notes: 824 Currently you cannot use this to set the actual subdomains with the argument is. 825 826 By default the ASM preconditioner uses 1 block per processor. 827 828 These index sets cannot be destroyed until after completion of the 829 linear solves for which the ASM preconditioner is being used. 830 831 Use PCASMSetLocalSubdomains() to set local subdomains. 832 833 The IS numbering is in the parallel, global numbering of the vector for both is and is_local 834 835 Level: advanced 836 837 .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 838 839 .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 840 PCASMCreateSubdomains2D() 841 @*/ 842 PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 843 { 844 PetscErrorCode ierr; 845 846 PetscFunctionBegin; 847 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 848 ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 849 PetscFunctionReturn(0); 850 } 851 852 #undef __FUNCT__ 853 #define __FUNCT__ "PCASMSetOverlap" 854 /*@ 855 PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 856 additive Schwarz preconditioner. Either all or no processors in the 857 PC communicator must call this routine. 858 859 Logically Collective on PC 860 861 Input Parameters: 862 + pc - the preconditioner context 863 - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 864 865 Options Database Key: 866 . -pc_asm_overlap <ovl> - Sets overlap 867 868 Notes: 869 By default the ASM preconditioner uses 1 block per processor. To use 870 multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 871 PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 872 873 The overlap defaults to 1, so if one desires that no additional 874 overlap be computed beyond what may have been set with a call to 875 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 876 must be set to be 0. In particular, if one does not explicitly set 877 the subdomains an application code, then all overlap would be computed 878 internally by PETSc, and using an overlap of 0 would result in an ASM 879 variant that is equivalent to the block Jacobi preconditioner. 880 881 Note that one can define initial index sets with any overlap via 882 PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(); the routine 883 PCASMSetOverlap() merely allows PETSc to extend that overlap further 884 if desired. 885 886 Level: intermediate 887 888 .keywords: PC, ASM, set, overlap 889 890 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 891 PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains() 892 @*/ 893 PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 894 { 895 PetscErrorCode ierr; 896 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 899 PetscValidLogicalCollectiveInt(pc,ovl,2); 900 ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 901 PetscFunctionReturn(0); 902 } 903 904 #undef __FUNCT__ 905 #define __FUNCT__ "PCASMSetType" 906 /*@ 907 PCASMSetType - Sets the type of restriction and interpolation used 908 for local problems in the additive Schwarz method. 909 910 Logically Collective on PC 911 912 Input Parameters: 913 + pc - the preconditioner context 914 - type - variant of ASM, one of 915 .vb 916 PC_ASM_BASIC - full interpolation and restriction 917 PC_ASM_RESTRICT - full restriction, local processor interpolation 918 PC_ASM_INTERPOLATE - full interpolation, local processor restriction 919 PC_ASM_NONE - local processor restriction and interpolation 920 .ve 921 922 Options Database Key: 923 . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 924 925 Level: intermediate 926 927 .keywords: PC, ASM, set, type 928 929 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 930 PCASMCreateSubdomains2D() 931 @*/ 932 PetscErrorCode PCASMSetType(PC pc,PCASMType type) 933 { 934 PetscErrorCode ierr; 935 936 PetscFunctionBegin; 937 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 938 PetscValidLogicalCollectiveEnum(pc,type,2); 939 ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 940 PetscFunctionReturn(0); 941 } 942 943 #undef __FUNCT__ 944 #define __FUNCT__ "PCASMSetSortIndices" 945 /*@ 946 PCASMSetSortIndices - Determines whether subdomain indices are sorted. 947 948 Logically Collective on PC 949 950 Input Parameters: 951 + pc - the preconditioner context 952 - doSort - sort the subdomain indices 953 954 Level: intermediate 955 956 .keywords: PC, ASM, set, type 957 958 .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 959 PCASMCreateSubdomains2D() 960 @*/ 961 PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 962 { 963 PetscErrorCode ierr; 964 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 967 PetscValidLogicalCollectiveBool(pc,doSort,2); 968 ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 969 PetscFunctionReturn(0); 970 } 971 972 #undef __FUNCT__ 973 #define __FUNCT__ "PCASMGetSubKSP" 974 /*@C 975 PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 976 this processor. 977 978 Collective on PC iff first_local is requested 979 980 Input Parameter: 981 . pc - the preconditioner context 982 983 Output Parameters: 984 + n_local - the number of blocks on this processor or NULL 985 . first_local - the global number of the first block on this processor or NULL, 986 all processors must request or all must pass NULL 987 - ksp - the array of KSP contexts 988 989 Note: 990 After PCASMGetSubKSP() the array of KSPes is not to be freed. 991 992 Currently for some matrix implementations only 1 block per processor 993 is supported. 994 995 You must call KSPSetUp() before calling PCASMGetSubKSP(). 996 997 Fortran note: 998 The output argument 'ksp' must be an array of sufficient length or NULL_OBJECT. The latter can be used to learn the necessary length. 999 1000 Level: advanced 1001 1002 .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 1003 1004 .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 1005 PCASMCreateSubdomains2D(), 1006 @*/ 1007 PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1008 { 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1013 ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1014 PetscFunctionReturn(0); 1015 } 1016 1017 /* -------------------------------------------------------------------------------------*/ 1018 /*MC 1019 PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1020 its own KSP object. 1021 1022 Options Database Keys: 1023 + -pc_asm_truelocal - Activates PCASMSetUseTrueLocal() 1024 . -pc_asm_blocks <blks> - Sets total blocks 1025 . -pc_asm_overlap <ovl> - Sets overlap 1026 - -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1027 1028 IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1029 will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 1030 -pc_asm_type basic to use the standard ASM. 1031 1032 Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1033 than one processor. Defaults to one block per processor. 1034 1035 To set options on the solvers for each block append -sub_ to all the KSP, and PC 1036 options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1037 1038 To set the options on the solvers separate for each block call PCASMGetSubKSP() 1039 and set the options directly on the resulting KSP object (you can access its PC 1040 with KSPGetPC()) 1041 1042 1043 Level: beginner 1044 1045 Concepts: additive Schwarz method 1046 1047 References: 1048 An additive variant of the Schwarz alternating method for the case of many subregions 1049 M Dryja, OB Widlund - Courant Institute, New York University Technical report 1050 1051 Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1052 Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1053 1054 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1055 PCBJACOBI, PCASMSetUseTrueLocal(), PCASMGetSubKSP(), PCASMSetLocalSubdomains(), 1056 PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType() 1057 1058 M*/ 1059 1060 EXTERN_C_BEGIN 1061 #undef __FUNCT__ 1062 #define __FUNCT__ "PCCreate_ASM" 1063 PetscErrorCode PCCreate_ASM(PC pc) 1064 { 1065 PetscErrorCode ierr; 1066 PC_ASM *osm; 1067 1068 PetscFunctionBegin; 1069 ierr = PetscNewLog(pc,PC_ASM,&osm);CHKERRQ(ierr); 1070 1071 osm->n = PETSC_DECIDE; 1072 osm->n_local = 0; 1073 osm->n_local_true = PETSC_DECIDE; 1074 osm->overlap = 1; 1075 osm->ksp = 0; 1076 osm->restriction = 0; 1077 osm->localization = 0; 1078 osm->prolongation = 0; 1079 osm->x = 0; 1080 osm->y = 0; 1081 osm->y_local = 0; 1082 osm->is = 0; 1083 osm->is_local = 0; 1084 osm->mat = 0; 1085 osm->pmat = 0; 1086 osm->type = PC_ASM_RESTRICT; 1087 osm->same_local_solves = PETSC_TRUE; 1088 osm->sort_indices = PETSC_TRUE; 1089 osm->dm_subdomains = PETSC_FALSE; 1090 1091 pc->data = (void*)osm; 1092 pc->ops->apply = PCApply_ASM; 1093 pc->ops->applytranspose = PCApplyTranspose_ASM; 1094 pc->ops->setup = PCSetUp_ASM; 1095 pc->ops->reset = PCReset_ASM; 1096 pc->ops->destroy = PCDestroy_ASM; 1097 pc->ops->setfromoptions = PCSetFromOptions_ASM; 1098 pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 1099 pc->ops->view = PCView_ASM; 1100 pc->ops->applyrichardson = 0; 1101 1102 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetLocalSubdomains_C","PCASMSetLocalSubdomains_ASM",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1103 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetTotalSubdomains_C","PCASMSetTotalSubdomains_ASM",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1104 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetOverlap_C","PCASMSetOverlap_ASM",PCASMSetOverlap_ASM);CHKERRQ(ierr); 1105 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetType_C","PCASMSetType_ASM",PCASMSetType_ASM);CHKERRQ(ierr); 1106 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMSetSortIndices_C","PCASMSetSortIndices_ASM",PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1107 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCASMGetSubKSP_C","PCASMGetSubKSP_ASM",PCASMGetSubKSP_ASM);CHKERRQ(ierr); 1108 PetscFunctionReturn(0); 1109 } 1110 EXTERN_C_END 1111 1112 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "PCASMCreateSubdomains" 1115 /*@C 1116 PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1117 preconditioner for a any problem on a general grid. 1118 1119 Collective 1120 1121 Input Parameters: 1122 + A - The global matrix operator 1123 - n - the number of local blocks 1124 1125 Output Parameters: 1126 . outis - the array of index sets defining the subdomains 1127 1128 Level: advanced 1129 1130 Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 1131 from these if you use PCASMSetLocalSubdomains() 1132 1133 In the Fortran version you must provide the array outis[] already allocated of length n. 1134 1135 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 1136 1137 .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 1138 @*/ 1139 PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 1140 { 1141 MatPartitioning mpart; 1142 const char *prefix; 1143 PetscErrorCode (*f)(Mat,Mat*); 1144 PetscMPIInt size; 1145 PetscInt i,j,rstart,rend,bs; 1146 PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 1147 Mat Ad = NULL, adj; 1148 IS ispart,isnumb,*is; 1149 PetscErrorCode ierr; 1150 1151 PetscFunctionBegin; 1152 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1153 PetscValidPointer(outis,3); 1154 if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 1155 1156 /* Get prefix, row distribution, and block size */ 1157 ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1158 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1159 ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1160 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); 1161 1162 /* Get diagonal block from matrix if possible */ 1163 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 1164 ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**) (void))&f);CHKERRQ(ierr); 1165 if (f) { 1166 ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1167 } else if (size == 1) { 1168 Ad = A; 1169 } 1170 if (Ad) { 1171 ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1172 if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1173 } 1174 if (Ad && n > 1) { 1175 PetscBool match,done; 1176 /* Try to setup a good matrix partitioning if available */ 1177 ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1178 ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1179 ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1180 ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1181 if (!match) { 1182 ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1183 } 1184 if (!match) { /* assume a "good" partitioner is available */ 1185 PetscInt na; 1186 const PetscInt *ia,*ja; 1187 ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1188 if (done) { 1189 /* Build adjacency matrix by hand. Unfortunately a call to 1190 MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1191 remove the block-aij structure and we cannot expect 1192 MatPartitioning to split vertices as we need */ 1193 PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 1194 const PetscInt *row; 1195 nnz = 0; 1196 for (i=0; i<na; i++) { /* count number of nonzeros */ 1197 len = ia[i+1] - ia[i]; 1198 row = ja + ia[i]; 1199 for (j=0; j<len; j++) { 1200 if (row[j] == i) { /* don't count diagonal */ 1201 len--; break; 1202 } 1203 } 1204 nnz += len; 1205 } 1206 ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1207 ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1208 nnz = 0; 1209 iia[0] = 0; 1210 for (i=0; i<na; i++) { /* fill adjacency */ 1211 cnt = 0; 1212 len = ia[i+1] - ia[i]; 1213 row = ja + ia[i]; 1214 for (j=0; j<len; j++) { 1215 if (row[j] != i) { /* if not diagonal */ 1216 jja[nnz+cnt++] = row[j]; 1217 } 1218 } 1219 nnz += cnt; 1220 iia[i+1] = nnz; 1221 } 1222 /* Partitioning of the adjacency matrix */ 1223 ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1224 ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1225 ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1226 ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1227 ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1228 ierr = MatDestroy(&adj);CHKERRQ(ierr); 1229 foundpart = PETSC_TRUE; 1230 } 1231 ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1232 } 1233 ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1234 } 1235 1236 ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1237 *outis = is; 1238 1239 if (!foundpart) { 1240 1241 /* Partitioning by contiguous chunks of rows */ 1242 1243 PetscInt mbs = (rend-rstart)/bs; 1244 PetscInt start = rstart; 1245 for (i=0; i<n; i++) { 1246 PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1247 ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1248 start += count; 1249 } 1250 1251 } else { 1252 1253 /* Partitioning by adjacency of diagonal block */ 1254 1255 const PetscInt *numbering; 1256 PetscInt *count,nidx,*indices,*newidx,start=0; 1257 /* Get node count in each partition */ 1258 ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1259 ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1260 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1261 for (i=0; i<n; i++) count[i] *= bs; 1262 } 1263 /* Build indices from node numbering */ 1264 ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1265 ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1266 for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1267 ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1268 ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1269 ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1270 if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1271 ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1272 for (i=0; i<nidx; i++) { 1273 for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 1274 } 1275 ierr = PetscFree(indices);CHKERRQ(ierr); 1276 nidx *= bs; 1277 indices = newidx; 1278 } 1279 /* Shift to get global indices */ 1280 for (i=0; i<nidx; i++) indices[i] += rstart; 1281 1282 /* Build the index sets for each block */ 1283 for (i=0; i<n; i++) { 1284 ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1285 ierr = ISSort(is[i]);CHKERRQ(ierr); 1286 start += count[i]; 1287 } 1288 1289 ierr = PetscFree(count);CHKERRQ(ierr); 1290 ierr = PetscFree(indices);CHKERRQ(ierr); 1291 ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1292 ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1293 1294 } 1295 PetscFunctionReturn(0); 1296 } 1297 1298 #undef __FUNCT__ 1299 #define __FUNCT__ "PCASMDestroySubdomains" 1300 /*@C 1301 PCASMDestroySubdomains - Destroys the index sets created with 1302 PCASMCreateSubdomains(). Should be called after setting subdomains 1303 with PCASMSetLocalSubdomains(). 1304 1305 Collective 1306 1307 Input Parameters: 1308 + n - the number of index sets 1309 . is - the array of index sets 1310 - is_local - the array of local index sets, can be NULL 1311 1312 Level: advanced 1313 1314 .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 1315 1316 .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 1317 @*/ 1318 PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 1319 { 1320 PetscInt i; 1321 PetscErrorCode ierr; 1322 1323 PetscFunctionBegin; 1324 if (n <= 0) PetscFunctionReturn(0); 1325 if (is) { 1326 PetscValidPointer(is,2); 1327 for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 1328 ierr = PetscFree(is);CHKERRQ(ierr); 1329 } 1330 if (is_local) { 1331 PetscValidPointer(is_local,3); 1332 for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 1333 ierr = PetscFree(is_local);CHKERRQ(ierr); 1334 } 1335 PetscFunctionReturn(0); 1336 } 1337 1338 #undef __FUNCT__ 1339 #define __FUNCT__ "PCASMCreateSubdomains2D" 1340 /*@ 1341 PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1342 preconditioner for a two-dimensional problem on a regular grid. 1343 1344 Not Collective 1345 1346 Input Parameters: 1347 + m, n - the number of mesh points in the x and y directions 1348 . M, N - the number of subdomains in the x and y directions 1349 . dof - degrees of freedom per node 1350 - overlap - overlap in mesh lines 1351 1352 Output Parameters: 1353 + Nsub - the number of subdomains created 1354 . is - array of index sets defining overlapping (if overlap > 0) subdomains 1355 - is_local - array of index sets defining non-overlapping subdomains 1356 1357 Note: 1358 Presently PCAMSCreateSubdomains2d() is valid only for sequential 1359 preconditioners. More general related routines are 1360 PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 1361 1362 Level: advanced 1363 1364 .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 1365 1366 .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1367 PCASMSetOverlap() 1368 @*/ 1369 PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 1370 { 1371 PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 1372 PetscErrorCode ierr; 1373 PetscInt nidx,*idx,loc,ii,jj,count; 1374 1375 PetscFunctionBegin; 1376 if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 1377 1378 *Nsub = N*M; 1379 ierr = PetscMalloc((*Nsub)*sizeof(IS*),is);CHKERRQ(ierr); 1380 ierr = PetscMalloc((*Nsub)*sizeof(IS*),is_local);CHKERRQ(ierr); 1381 ystart = 0; 1382 loc_outer = 0; 1383 for (i=0; i<N; i++) { 1384 height = n/N + ((n % N) > i); /* height of subdomain */ 1385 if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 1386 yleft = ystart - overlap; if (yleft < 0) yleft = 0; 1387 yright = ystart + height + overlap; if (yright > n) yright = n; 1388 xstart = 0; 1389 for (j=0; j<M; j++) { 1390 width = m/M + ((m % M) > j); /* width of subdomain */ 1391 if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 1392 xleft = xstart - overlap; if (xleft < 0) xleft = 0; 1393 xright = xstart + width + overlap; if (xright > m) xright = m; 1394 nidx = (xright - xleft)*(yright - yleft); 1395 ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1396 loc = 0; 1397 for (ii=yleft; ii<yright; ii++) { 1398 count = m*ii + xleft; 1399 for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 1400 } 1401 ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 1402 if (overlap == 0) { 1403 ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 1404 1405 (*is_local)[loc_outer] = (*is)[loc_outer]; 1406 } else { 1407 for (loc=0,ii=ystart; ii<ystart+height; ii++) { 1408 for (jj=xstart; jj<xstart+width; jj++) { 1409 idx[loc++] = m*ii + jj; 1410 } 1411 } 1412 ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 1413 } 1414 ierr = PetscFree(idx);CHKERRQ(ierr); 1415 xstart += width; 1416 loc_outer++; 1417 } 1418 ystart += height; 1419 } 1420 for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 1421 PetscFunctionReturn(0); 1422 } 1423 1424 #undef __FUNCT__ 1425 #define __FUNCT__ "PCASMGetLocalSubdomains" 1426 /*@C 1427 PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1428 only) for the additive Schwarz preconditioner. 1429 1430 Not Collective 1431 1432 Input Parameter: 1433 . pc - the preconditioner context 1434 1435 Output Parameters: 1436 + n - the number of subdomains for this processor (default value = 1) 1437 . is - the index sets that define the subdomains for this processor 1438 - is_local - the index sets that define the local part of the subdomains for this processor (can be NULL) 1439 1440 1441 Notes: 1442 The IS numbering is in the parallel, global numbering of the vector. 1443 1444 Level: advanced 1445 1446 .keywords: PC, ASM, set, local, subdomains, additive Schwarz 1447 1448 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1449 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 1450 @*/ 1451 PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 1452 { 1453 PC_ASM *osm; 1454 PetscErrorCode ierr; 1455 PetscBool match; 1456 1457 PetscFunctionBegin; 1458 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1459 PetscValidIntPointer(n,2); 1460 if (is) PetscValidPointer(is,3); 1461 ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1462 if (!match) { 1463 if (n) *n = 0; 1464 if (is) *is = NULL; 1465 } else { 1466 osm = (PC_ASM*)pc->data; 1467 if (n) *n = osm->n_local_true; 1468 if (is) *is = osm->is; 1469 if (is_local) *is_local = osm->is_local; 1470 } 1471 PetscFunctionReturn(0); 1472 } 1473 1474 #undef __FUNCT__ 1475 #define __FUNCT__ "PCASMGetLocalSubmatrices" 1476 /*@C 1477 PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1478 only) for the additive Schwarz preconditioner. 1479 1480 Not Collective 1481 1482 Input Parameter: 1483 . pc - the preconditioner context 1484 1485 Output Parameters: 1486 + n - the number of matrices for this processor (default value = 1) 1487 - mat - the matrices 1488 1489 1490 Level: advanced 1491 1492 Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks()) 1493 1494 Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner. 1495 1496 .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 1497 1498 .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1499 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices() 1500 @*/ 1501 PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1502 { 1503 PC_ASM *osm; 1504 PetscErrorCode ierr; 1505 PetscBool match; 1506 1507 PetscFunctionBegin; 1508 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1509 PetscValidIntPointer(n,2); 1510 if (mat) PetscValidPointer(mat,3); 1511 if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1512 ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1513 if (!match) { 1514 if (n) *n = 0; 1515 if (mat) *mat = NULL; 1516 } else { 1517 osm = (PC_ASM*)pc->data; 1518 if (n) *n = osm->n_local_true; 1519 if (mat) *mat = osm->pmat; 1520 } 1521 PetscFunctionReturn(0); 1522 } 1523 1524 #undef __FUNCT__ 1525 #define __FUNCT__ "PCASMSetDMSubdomains" 1526 /*@ 1527 PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1528 Logically Collective 1529 1530 Input Parameter: 1531 + pc - the preconditioner 1532 - flg - boolean indicating whether to use subdomains defined by the DM 1533 1534 Options Database Key: 1535 . -pc_asm_dm_subdomains 1536 1537 Level: intermediate 1538 1539 Notes: 1540 PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1541 so setting either of the first two effectively turns the latter off. 1542 1543 .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1544 1545 .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1546 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1547 @*/ 1548 PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1549 { 1550 PC_ASM *osm = (PC_ASM*)pc->data; 1551 PetscErrorCode ierr; 1552 PetscBool match; 1553 1554 PetscFunctionBegin; 1555 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1556 PetscValidLogicalCollectiveBool(pc,flg,2); 1557 if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1558 ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1559 if (match) { 1560 osm->dm_subdomains = flg; 1561 } 1562 PetscFunctionReturn(0); 1563 } 1564 1565 #undef __FUNCT__ 1566 #define __FUNCT__ "PCASMGetDMSubdomains" 1567 /*@ 1568 PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1569 Not Collective 1570 1571 Input Parameter: 1572 . pc - the preconditioner 1573 1574 Output Parameter: 1575 . flg - boolean indicating whether to use subdomains defined by the DM 1576 1577 Level: intermediate 1578 1579 .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1580 1581 .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1582 PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1583 @*/ 1584 PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1585 { 1586 PC_ASM *osm = (PC_ASM*)pc->data; 1587 PetscErrorCode ierr; 1588 PetscBool match; 1589 1590 PetscFunctionBegin; 1591 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1592 PetscValidPointer(flg,2); 1593 ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1594 if (match) { 1595 if (flg) *flg = osm->dm_subdomains; 1596 } 1597 PetscFunctionReturn(0); 1598 } 1599