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