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