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