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