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