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