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