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